source: trunk/SOURCES/Old-sources/lect-Ant_gen2010_dat.f90 @ 334

Last change on this file since 334 was 236, checked in by aquiquet, 5 years ago

Deprecated modules moved to Old-sources

File size: 12.9 KB
Line 
1!> \file lect-Ant_CISM_gen_dat.f90
2!!Module pour la lecture de la topographie
3!<
4
5!> \namespace lect_topo_ant_CISM_gen
6!! Module pour la lecture de la topographie
7!! \author ...
8!! \date ...
9!! @note Used module
10!! @note   - use module3D_phy
11!! @note   - use param_phy_mod
12!<
13
14module lect_topo_ant_gen
15
16  use module3D_phy
17  use param_phy_mod
18  use interface_input
19  use io_netcdf_grisli
20
21  character(len=100) :: topo_surf      !< surface file name
22  character(len=100) :: correc_surf    !< ice-real correction file name
23  character(len=100) :: topo_thick     !< thickness file name
24  character(len=100) :: topo_bed       !< bedrock file name
25  character(len=100) :: mask_grounded  !< mask file name
26  character(len=100) :: longitude      !< longitude file name
27  character(len=100) :: latitude       !< latitude file name
28  character(len=100) :: heatflux       !< geothermal heat flux
29  logical            :: corr_firn      !< eventuellement correction firn (deja fait dans Le Brocq)
30
31!  for a start with a different topography than the reference one (without temperature inside)
32  integer            :: istart = 0     !< to read another topography than the present one > 0, (2 when  icptr=1)
33  character(len=100) :: surf_start     !< surface file name
34  character(len=100) :: thick_start    !< thickness file name
35  character(len=100) :: bed_start      !< bed file name
36  character(len=100) :: mask_start     !< masque file name
37
38  real,dimension(nx,ny) :: work_tab 
39  real*8,dimension(:),pointer        :: Tab1D           ! pointeur vers les tableaux des coordonnees
40
41  character(len=100) :: file_ncdf      !< fichier netcdf issue des fichiers .dat
42
43contains
44
45  subroutine input_topo
46
47    namelist/topo_ant_gen/topo_surf,correc_surf,topo_thick,topo_bed,mask_grounded,longitude,latitude,heatflux
48    namelist/topo_ant_startingpoint/surf_start,thick_start,bed_start,mask_start,sealevel
49
50    if (itracebug.eq.1)  call tracebug(' Entree dans routine input_topo')
51
52428 format(A)
53
54    rewind(num_param)                     ! pour revenir au debut du fichier param_list.dat
55    read(num_param,topo_ant_gen)
56
57    write(num_rep_42,428)'!___________________________________________________________' 
58    write(num_rep_42,428)'!  module  lect_topo_ant_gen                                '
59    write(num_rep_42,topo_ant_gen)
60    write(num_rep_42,428)'!___________________________________________________________' 
61
62
63    if ((trim(correc_surf).eq.'no').or.(trim(correc_surf).eq.'NO')) then
64       corr_firn=.false.
65    else
66       corr_firn=.true.
67    end if
68
69    topo_surf     = trim(dirnameinp)//trim(topo_surf)
70    correc_surf   = trim(dirnameinp)//trim(correc_surf)
71    topo_thick    = trim(dirnameinp)//trim(topo_thick)
72    topo_bed      = trim(dirnameinp)//trim(topo_bed)
73    mask_grounded = trim(dirnameinp)//trim(mask_grounded)
74    longitude     = trim(dirnameinp)//trim(longitude)
75    latitude      = trim(dirnameinp)//trim(latitude)
76    heatflux      = trim(dirnameinp)//trim(heatflux)
77
78
79    ! Transformation des fichiers .dat en .nc
80
81    file_ncdf= trim(dirnameinp)//trim(runname)//'.nc'
82
83    ! lecture adaptee aux fichiers ZBL.dat
84
85    call lect_input(1,'Bsoc',1,Bsoc,topo_bed,file_ncdf)
86    !call lect_datfile(nx,ny,Bsoc,1,topo_bed)               ! bedrock
87    Bsoc0(:,:) = Bsoc(:,:)
88
89    call lect_input(1,'S',1,S,topo_surf,file_ncdf)
90    !call lect_datfile(nx,ny,S,1,topo_surf)                 ! surface
91   
92    call lect_input(1,'H',1,H,topo_thick,file_ncdf)
93    !call lect_datfile(nx,ny,H,1,topo_thick)                ! thickness
94
95    if (itracebug.eq.1)  call tracebug(' avant corr firn')
96
97    if (corr_firn) then                                    ! density correction
98       !call lect_datfile(nx,ny,work_tab,1,correc_surf)     ! for CISM files but not
99       ! Cat treatment of Le Brocq
100       ! already done
101       
102       call lect_input(1,'work_tab1',1,work_tab,correc_surf,file_ncdf)
103
104       where (H(:,:).gt.0.)
105          work_tab(:,:)=min(H(:,:),work_tab(:,:))   
106          S(:,:) = S(:,:) - work_tab(:,:)                  ! corrected surface
107          H(:,:) = H(:,:) - work_tab(:,:)                  ! corrected thickness
108       end where
109    end if
110    if (itracebug.eq.1)  call tracebug(' apres corr firn')
111    H(:,:)=max(0.,H(:,:))
112    H0(:,:)= H(:,:)
113    S0(:,:)= S(:,:)
114
115
116    mk0(:,:)=0.
117    mk_init(:,:)=0.
118
119    call lect_input(1,'work_tab2',2,work_tab,mask_grounded,file_ncdf)
120    !call lect_datfile(nx,ny,work_tab,1,mask_grounded)      ! masque
121
122    ! definition du Mk_init
123    !__________________________________________
124    !  outcrops         -> -4
125    !  ocean libre      -> -3
126    !  iles a pb        -> -2
127    !  pose avec glace  ->  0
128    !  ice streams      ->  1 d'apres Le Brocq
129    !  ice rise         ->  2 pas attribue ici
130    !  ice shelves      ->  3
131
132    mk_init(:,:) = work_tab(:,:)
133
134    if (itracebug.eq.1)  call tracebug(' avant worktab')
135
136
137    ! to remove nasty islands transform them in ocean
138!!$    where (mk_init(:,:).eq.-2)    ! remove the test on H .and.(H(:,:).gt.1.))
139!!$       H(:,:) =1     
140!!$       S(:,:)= coef_Sflot * H(:,:)
141!!$       Bsoc(:,:)=min(Bsoc(:,:),-100.)
142!!$       B(:,:)= coef_Bflot * H(:,:)
143!!$    end where
144
145
146! test consistency of input data files (flotation, thickness ...)
147
148! Le masque (moyenne ou sampling)  n'est pas forcement exact
149! le calcul de flot est plus fiable et est base sur l'epaisseur
150
151    do j=1,ny
152       do i=1,nx
153          if ((Bsoc(i,j)+H(i,j)*ro/row -sealevel).lt.0.) then   ! flottaison
154             flot(i,j)=.true.
155
156             if (mk_init(i,j).eq.-2) then    ! ile a probleme : on fait une montagne de 100 m
157                Bsoc(i,j) = max(Bsoc(i,j),sealevel + 100.)
158                S(i,j)    = max(S(i,j),Bsoc(i,j))
159                H(i,j)    = S(i,j) - Bsoc(i,j)
160                B(i,j)    = Bsoc(i,j)
161                ! dans le dragging, drag_centre = beta_limgz
162
163             else
164
165                if (H(i,j).gt.0.) then                ! ice shelf
166                   mk_init(i,j) = 3       
167                   S(i,j)       = coef_Sflot * H(i,j) + sealevel
168                   B(i,j)       = S(i,j) - H(i,j)
169
170                else if (H(i,j).le.0.) then                   
171                   if (Bsoc(i,j).lt.sealevel) then
172                      mk_init(i,j) = -3                ! ocean libre
173                      H(i,j)       =  0.
174                      S(i,j)       =  0.
175                      B(i,j)       =  0.
176                   else
177                      mk_init(i,j) = -4                ! outcrop
178                      H(i,j)       =  0.
179                      S(i,j)       = Bsoc(i,j)
180                      B(i,j)       = Bsoc(i,j)
181                   end if
182
183                end if
184             endif
185
186          else                                       ! point pose
187             flot(i,j)=.false.
188             if ( (mk_init(i,j).eq.-4).or.(mk_init(i,j).eq.1).or.(mk_init(i,j).eq.2) ) then 
189                ! outcrop ou ice streams ou ice rise : on ne change pas
190
191             else if (mk_init(i,j).eq.-2) then    ! ile a probleme : on fait une montagne de 100 m
192                                                  ! dans le dragging, drag_centre = betalimgz
193
194                Bsoc(i,j) = max(Bsoc(i,j),sealevel + 100.)
195                S(i,j)    = min(S(i,j),Bsoc(i,j))
196                H(i,j)    = S(i,j) - Bsoc(i,j)
197                B(i,j)    = Bsoc(i,j)
198             else
199                if (H(i,j).gt.0) then
200                   mk_init(i,j) = 0               ! cas general du point pose
201                else
202                   mk_init(i,j) = -4              ! pas de glace -> outcrop
203                end if
204             end if
205
206                                                  ! dans tous les cas  poses           
207             S(i,j)    = Bsoc(i,j) + H(i,j)
208             B(i,j)    = Bsoc(i,j)
209
210          endif
211       enddo
212    enddo
213
214
215! on laisse ce passage mais il ne devrait plus y avoir de probleme.
216
217    do J=1,ny
218       do i=1,nx
219          if ((abs(mk_init(i,j)).ne.3)) then         ! points poses
220
221             if (abs((S(i,j)-Bsoc(i,j))-H(i,j)).ge.1.e-2) then
222                write(6,*) i,j,mk_init(i,j),'S-B =! H',(S(i,j)-Bsoc(i,j))-H(i,j)
223             endif
224
225             if (coef_Sflot*Bsoc(i,j)/coef_Bflot.gt.S(i,j)) then
226                write(6,*) i,j,mk_init(i,j),'flotte S:',coef_Sflot*Bsoc(i,j)/coef_Bflot, S(i,j)
227             end if
228
229          else if ((abs(mk_init(i,j)).eq.3)) then   ! points flottants
230             if (abs(S(i,j) - coef_Sflot * H(i,j)).gt.1.) then
231                write(6,*)
232                write(6,*) i,j,mk_init(i,j),'mauvaise flottaison',S(i,j) - coef_Sflot * H(i,j),'S H S2 B Bsoc B2'
233                write(6,'(7(f0.2,2x))') S(i,j), H(i,j), coef_Sflot * H(i,j), B(i,j), Bsoc(i,j),coef_Bflot* H(i,j)
234               
235             endif
236
237             if (coef_Bflot* H(i,j).lt.Bsoc(i,j)-1.) then
238                write(6,*) i,j,mk_init(i,j),'frotte B: ',coef_Bflot* H(i,j), Bsoc(i,j)
239             end if
240
241
242          end if
243       end do
244    end do
245    Bsoc0(:,:) = Bsoc(:,:)
246    H0(:,:)= H(:,:)
247    S0(:,:)= S(:,:)
248
249    if (itracebug.eq.1)  call tracebug(' apres worktab')
250
251! ecriture des variables apres consistency
252!!$do j=1,ny
253!!$   do i=1,nx
254!!$      write(122,*) S(i,j),H(i,j),Bsoc(i,j),B(i,j),Mk_init(i,j)
255!!$   end do
256!!$end do
257
258
259
260    ! pour l'Antarctique masque mko vrai partout (version 2006)
261    MK0(:,:)=1  ! mis apres la lecture cptr
262
263
264other_start:    if (istart.ge.1) then                 ! restart d'une nouvelle topo
265
266    rewind(num_param)                     ! pour revenir au debut du fichier param_list.dat
267    read(num_param,topo_ant_startingpoint)
268
269    write(num_rep_42,428)'!___________________________________________________________' 
270    write(num_rep_42,428)'!  module  lect_topo_ant_gen : partir differente topo       '
271    write(num_rep_42,topo_ant_gen)
272    write(num_rep_42,428)'!___________________________________________________________' 
273
274    surf_start   = trim(dirnameinp)//trim(surf_start)
275    thick_start  = trim(dirnameinp)//trim(thick_start)
276    bed_start    = trim(dirnameinp)//trim(bed_start)
277    mask_start   = trim(dirnameinp)//trim(mask_start)
278
279 endif other_start
280
281    if (istart.eq.1) then     
282 
283    call lect_input(1,'S',1,S,surf_start,file_ncdf)
284    call lect_input(1,'H',1,H,thick_start,file_ncdf)
285    call lect_input(1,'Bsoc',1,H,thick_start,file_ncdf)
286    call lect_input(1,'work_tab1',1,work_tab,mask_start,file_ncdf)
287
288    mk_init(:,:) = work_tab(:,:)
289
290    where ( mk_init(:,:).eq.1)
291       flot(:,:) = .False.
292    elsewhere
293       flot(:,:) = .True.
294    end where
295
296 else if (istart.eq.2) then                     ! juste lecture du masque
297
298   if (itracebug.eq.1)  call tracebug(' dans test istart = 2')
299
300    call lect_input(1,'work_tab1',1,work_tab,mask_start,file_ncdf)
301
302    mk_init(:,:) = work_tab(:,:)
303   if (itracebug.eq.1)  call tracebug(' dans test istart = 2')
304    where ( mk_init(:,:).eq.1)
305       flot(:,:) = .False.
306    elsewhere
307       flot(:,:) = .True.
308    end where
309
310
311 end if 
312
313
314    ! le calcul des courbures du socle devrait etre dans le dragging
315   if (itracebug.eq.1)  call tracebug(' avant lecture lat long')
316    call lect_input(1,'Xlong_1',1,Xlong,longitude,file_ncdf)
317    !call lect_datfile(nx,ny,Xlong,1,longitude)      ! lecture des coordonnées geographiques
318   
319    call lect_input(1,'Ylat_1',1,Ylat,latitude,file_ncdf)
320    !call lect_datfile(nx,ny,Ylat,1,latitude)        !
321
322   if (itracebug.eq.1)  call tracebug(' apres lecture lat long')
323
324! xmin ymin xmax ymax
325!------------------------------------------------
326! dans un des fichier grd, lire les xmin et xmax
327! Il faudra que ce fichier au moins aie les bonnes valeurs
328
329   write(6,*) 'topo_thick',topo_thick
330   call Read_Ncdf_var('x',trim(topo_thick),Tab1D)    ! lit le tableau x
331   xmin = minval(Tab1D(:))
332   xmax = maxval(Tab1d(:))
333
334   call Read_Ncdf_var('y',trim(topo_thick),Tab1D)    ! lit le tableau y
335   ymin = minval(Tab1D(:))
336   ymax = maxval(Tab1d(:))
337
338 
339
340!   longitude =trim(dirnameinp)//'long_ZBL_LBq_15km.dat'
341
342!    open(22,file=longitude)                         ! juste pour avoir l'en tete
343!    read(22,*) i,i,xmin,xmax,ymin,ymax
344    xmin=xmin/1000.
345    ymin=ymin/1000.
346    xmax=xmax/1000.
347    ymax=ymax/1000.
348!    close(22)
349!   write(6,*) 'bornes domaine a la lecture',xmin,xmax,ymin,ymax
350
351
352    !     lecture du fichier de reference pour le calcul du niveau des mers : etat actuel
353    !     pas de version correcte pour l'instant -> valeurs initiales
354    S_sealev(:,:) = S0(:,:)
355    H_sealev(:,:) = H0(:,:)
356    B_sealev(:,:) = Bsoc(:,:)
357    M_sealev(:,:) = Mk0(:,:)
358
359    if (itracebug.eq.1)  call tracebug(' avant lecture ghf')
360
361    call lect_input(1,'ghf',1,ghf,heatflux,file_ncdf)
362    !call lect_datfile(nx,ny,ghf,1,heatflux)
363
364    if (itracebug.eq.1)  call tracebug(' apres lecture ghf')
365
366    ! for the purposes of SeaRISE, GTHF should be capped at a value of -0.07 w/m^2,
367    ! and the field provided in the netCDF file does not include these updated values
368
369    ghf(:,:)=min(ghf(:,:),70.)           ! pour sea rise
370
371    ! pour passer les flux des W/m2 au J/m2/an     
372    ghf(:,:)=-SECYEAR*ghf(:,:)/1000.     ! /1000. parce que les flux sont donnes en mW/m2
373    ghf0(:,:)=ghf(:,:)
374
375
376  end subroutine input_topo
377
378
379end module lect_topo_ant_gen
380
Note: See TracBrowser for help on using the repository browser.