source: trunk/LMDZ.VENUS/libf/phyvenus/phyetat0.F90 @ 3892

Last change on this file since 3892 was 3892, checked in by emillour, 4 months ago

Venus PCM:
Some code cleanup:

  • remove clmain.old and clmain.classic
  • remove all the unnecessary "EXTERNAL" statements in physiq
  • turn hgardfou into a module (and remove obsolete cpp directives for CRAY computers)

EM

File size: 11.5 KB
Line 
1module phyetat0_mod
2
3implicit none
4
5contains
6
7subroutine phyetat0(fichnom)
8! Load initial state for the physics
9! and do some resulting initializations
10
11  use dimphy, only: klon
12  use iophy, only: init_iophy_new
13  use phys_state_var_mod, only: ftsol, ftsoil, age, dlw, falbe, fder, q2, &
14                                radsol, sollw, sollwdown, solsw, &
15                                ancien_ok, t_ancien, &
16                                zmea, zgam, zpic, zsig, zstd, zthe, zval
17  use iostart, only: get_var, get_field, open_startphy, close_startphy
18  use geometry_mod, only: longitude_deg, latitude_deg
19  use time_phylmdz_mod, only: itau_phy, raz_date, pdtphys
20  use ioipsl_getin_p_mod, only: getin_p
21  use nonoro_gwd_ran_mod, only: du_nonoro_gwd, dv_nonoro_gwd, &
22                                east_gwstress, west_gwstress
23  use tabcontrol_mod, only: dtime, radpas, length, tabcntr0
24  use soil_mod, only: nsoilmx
25  use clesphys_mod, only: startphy_file
26
27implicit none
28
29character(len=*),intent(in) :: fichnom ! input file name
30LOGICAL :: found
31REAL    :: tab_cntrl(length)
32integer :: i,isoil
33CHARACTER(len=2) :: str2
34REAL :: lon_startphy(klon), lat_startphy(klon)
35REAL :: surface_albedo
36character(len=8) :: modname="phyetat0"
37
38! global variables read in startphy.nc file
39
40! open physics initial state file:
41if (startphy_file) then
42  call open_startphy(fichnom)
43endif
44!
45! Load control parameters:
46!
47IF (startphy_file) THEN
48  CALL get_var("controle",tab_cntrl,found)
49  IF (.not.found) THEN
50    write(*,*) modname//': Field <controle> is missing'
51    call abort_physic(modname,"Missing <controle>",1)
52  ENDIF
53       
54  DO i = 1, length
55    tabcntr0( i ) = tab_cntrl( i )
56  ENDDO
57
58  dtime        = tab_cntrl(1)
59  radpas       = tab_cntrl(2)
60
61  itau_phy = tab_cntrl(15)
62
63! Warning, if raz_date is active :
64! itau_phy must be re-set to zero after phyetat0 !
65  IF (raz_date.eq.1) THEN
66    itau_phy=0
67  ENDIF
68
69ELSE
70  tabcntr0(:)=1 ! dummy initialization
71  ! Initialize parameter or get values from def files
72  dtime=pdtphys
73  radpas=1
74  itau_phy=0
75ENDIF ! of IF (startphy_file)
76
77IF (startphy_file) THEN
78  ! read latitudes and make a sanity check (because already known from dyn)
79  call get_field("latitude",lat_startphy,found)
80  IF (.not.found) THEN
81    write(*,*) modname//':: Field <latitude> is missing'
82    call abort_physic(modname,"Missing <latitude>",1)
83  ENDIF
84  DO i=1,klon
85    IF (ABS(lat_startphy(i)-latitude_deg(i))>=0.01) THEN
86      WRITE(*,*) modname//": Warning! Latitude discrepancy wrt startphy file:",&
87                 " i=",i," lat_startphy(i)=",lat_startphy(i),&
88                 " latitude_deg(i)=",latitude_deg(i)
89      call abort_physic(modname,"<latitude> values discrepency",1)
90    ENDIF
91  ENDDO
92
93  ! read longitudes and make a sanity check (because already known from dyn)
94  call get_field("longitude",lon_startphy,found)
95  IF (.not.found) THEN
96    write(*,*)'phyetat0: Field <longitude> is missing'
97    call abort_physic(modname,"Missing <longitude>",1)
98  ENDIF
99  DO i=1,klon
100    IF (ABS(lon_startphy(i)-longitude_deg(i))>=0.01) THEN
101      WRITE(*,*) modname//": Warning! Longitude discrepancy wrt startphy file:",&
102                 " i=",i," lon_startphy(i)=",lon_startphy(i),&
103                 " longitude_deg(i)=",longitude_deg(i)
104      call abort_physic(modname,"<longitude> values discrepency",1)
105    ENDIF
106  ENDDO
107ENDIF ! of IF (startphy_file)
108
109! read in other variables here ...
110
111IF (startphy_file) THEN
112  ! Load surface temperature:
113  CALL get_field("TS",ftsol(:),found)
114  IF (.not.found) THEN
115    WRITE(*,*) modname//": Field <TS> is missing"
116    call abort_physic(modname,"Missing <TS>",1)
117  ELSE
118    WRITE(*,*) modname//": Field <TS> is present"
119    WRITE(*,*) 'Surface temperature <TS>', minval(ftsol), maxval(ftsol)
120  ENDIF
121ELSE
122  ! Dummy initialization, but in fact this is later handled in physiq
123  ftsol(:)=0
124ENDIF ! of IF (startphy_file)
125
126IF (startphy_file) THEN
127  ! Load sub-surface temperatures:
128  DO isoil=1, nsoilmx
129    IF (isoil.GT.99) THEN
130       WRITE(*,*) "Too many layers!"
131       call abort_physic(modname,"Too many subsurface layers",1)
132    ENDIF
133    WRITE(str2,'(i2.2)') isoil
134    CALL get_field('Tsoil'//str2,ftsoil(:,isoil),found)
135    IF (.not.found) THEN
136      WRITE(*,*) modname//": Field <Tsoil"//str2//"> is missing"
137      WRITE(*,*) "           it will be nitialized with surface value"
138      DO i=1, klon
139             ftsoil(i,isoil)=ftsol(i)
140      ENDDO
141    ENDIF
142  ENDDO
143ELSE
144  ! Dummy initialization, but in fact this is later handled in physiq
145  ftsoil(:,:)=0
146ENDIF ! of IF (startphy_file)
147
148IF (startphy_file) THEN
149  ! Load surface albedo:
150  CALL get_field("ALBE", falbe,found)
151  IF (.not.found) THEN
152    WRITE(*,*) modname//": Field <ALBE> is missing"
153    call abort_physic(modname,"Missing <ALBE>",1)
154  ENDIF
155ELSE
156  ! Dummy initialization: read value from def file
157  surface_albedo=0.5 ! default
158  CALL getin_p("surface_albedo",surface_albedo)
159  falbe(:)=surface_albedo
160ENDIF ! of IF (startphy_file)
161WRITE(*,*) 'Surface albedo <ALBE>', minval(falbe), maxval(falbe)
162
163IF (startphy_file) THEN
164  ! Solar flux to the surface:
165  CALL get_field("solsw",solsw,found)
166  IF (.not.found) THEN
167    WRITE(*,*) modname//": Field <solsw> is missing"
168    WRITE(*,*) "set to zero"
169    solsw = 0.
170  ENDIF
171ELSE
172  ! Dummy initialization
173  solsw(:)=0
174ENDIF ! of IF (startphy_file)
175WRITE(*,*) 'Solar flux to the surface solsw:', minval(solsw), maxval(solsw)
176
177IF (startphy_file) THEN
178  ! IR flux to the surface:
179  CALL get_field("sollw",sollw,found)
180  IF (.not.found) THEN
181    WRITE(*,*) modname//": Field <sollw> is missing"
182    WRITE(*,*) "set to zero"
183    sollw = 0.
184  ENDIF
185ELSE
186  ! Dummy initialization
187  sollw(:)=0
188ENDIF ! of IF (startphy_file)
189WRITE(*,*) 'IR flux to the surface sollw:', minval(sollw), maxval(solsw)
190
191IF (startphy_file) THEN
192  ! Derivative of the sensible and latent fluxes:
193  CALL get_field("fder",fder,found)
194  IF (.not.found) THEN
195    WRITE(*,*) modname//": Field <fder> is missing"
196    WRITE(*,*) "set to zero"
197    fder = 0.
198  ENDIF
199ELSE
200  ! Dummy initialization
201  fder(:)=0
202ENDIF ! of IF (startphy_file)
203WRITE(*,*) 'Derivative of the flux fder:', minval(fder), maxval(fder)
204
205IF (startphy_file) THEN
206  ! Derivative of the IR flux:
207  CALL get_field("dlw",dlw,found)
208  IF (.not.found) THEN
209    WRITE(*,*) modname//": Field <dlw> is missing"
210    WRITE(*,*) "set to zero"
211    dlw = 0.
212  ENDIF
213ELSE
214  ! Dummy initialization
215  dlw(:)=0
216ENDIF ! of IF (startphy_file)
217WRITE(*,*) 'Derivative of the IR flux dlw:', minval(dlw), maxval(dlw)
218
219IF (startphy_file) THEN
220  ! Downward IR flux at the surface:
221  CALL get_field("sollwdown",sollwdown,found)
222  IF (.not.found) THEN
223    WRITE(*,*) modname//": Field <sollwdown> is missing"
224    WRITE(*,*) "set to zero"
225    sollwdown = 0.
226  ENDIF
227ELSE
228  ! Dummy initialization
229  sollwdown(:)=0
230ENDIF ! of IF (startphy_file)
231WRITE(*,*) 'Downward IR flux at the surface sollwdown:', minval(sollwdown), maxval(sollwdown)
232
233IF (startphy_file) THEN
234  ! Net flux at the surface:
235  CALL get_field("RADS",radsol,found)
236  IF (.not.found) THEN
237    WRITE(*,*) modname//": Field <RADS> is missing"
238    call abort_physic(modname,"Missing <RADS>",1)
239  ENDIF
240ELSE
241  ! Dummy initialization
242  radsol(:)=0
243ENDIF ! of IF (startphy_file)
244WRITE(*,*) 'Net flux at the surface radsol:', minval(radsol), maxval(radsol)
245
246IF (startphy_file) THEN
247  ! Load sub-grid scale orography parameters:
248  CALL get_field("ZMEA",zmea,found)
249  IF (.not.found) THEN
250    WRITE(*,*) modname//": Field <ZMEA> is missing"
251    WRITE(*,*) "set to zero"
252    zmea=0.
253  ENDIF
254ELSE
255  zmea(:)=0
256ENDIF ! of IF (startphy_file)
257WRITE(*,*) 'Subgrid scale orography zmea:', minval(zmea), maxval(zmea)
258
259IF (startphy_file) THEN
260  ! Load sub-grid scale orography parameters:
261  CALL get_field("ZSTD",zstd,found)
262  IF (.not.found) THEN
263    WRITE(*,*) modname//": Field <ZSTD> is missing"
264    WRITE(*,*) "set to zero"
265    zstd=0.
266  ENDIF
267ELSE
268  zstd(:)=0
269ENDIF ! of IF (startphy_file)
270WRITE(*,*) 'Subgrid scale orography zstd:', minval(zstd), maxval(zstd)
271
272IF (startphy_file) THEN
273  ! Load sub-grid scale orography parameters:
274  CALL get_field("ZSIG",zsig,found)
275  IF (.not.found) THEN
276    WRITE(*,*) modname//": Field <ZSIG> is missing"
277    WRITE(*,*) "set to zero"
278    zsig=0.
279  ENDIF
280ELSE
281  zsig(:)=0
282ENDIF ! of IF (startphy_file)
283WRITE(*,*) 'Subgrid scale orography zsig:', minval(zsig), maxval(zsig)
284
285IF (startphy_file) THEN
286  ! Load sub-grid scale orography parameters:
287  CALL get_field("ZGAM",zgam,found)
288  IF (.not.found) THEN
289    WRITE(*,*) modname//": Field <ZGAM> is missing"
290    WRITE(*,*) "set to zero"
291    zgam=0.
292  ENDIF
293ELSE
294  zgam(:)=0
295ENDIF ! of IF (startphy_file)
296WRITE(*,*) 'Subgrid scale orography zgam:', minval(zgam), maxval(zgam)
297
298IF (startphy_file) THEN
299  ! Load sub-grid scale orography parameters:
300  CALL get_field("ZTHE",zthe,found)
301  IF (.not.found) THEN
302    WRITE(*,*) modname//": Field <ZTHE> is missing"
303    WRITE(*,*) "set to zero"
304    zthe=0.
305  ENDIF
306ELSE
307  zthe(:)=0
308ENDIF ! of IF (startphy_file)
309WRITE(*,*) 'Subgrid scale orography zthe:', minval(zthe), maxval(zthe)
310
311IF (startphy_file) THEN
312  ! Load sub-grid scale orography parameters:
313  CALL get_field("ZPIC",zpic,found)
314  IF (.not.found) THEN
315    WRITE(*,*) modname//": Field <ZPIC> is missing"
316    WRITE(*,*) "set to zero"
317    zpic=0.
318  ENDIF
319ELSE
320  zpic(:)=0
321ENDIF ! of IF (startphy_file)
322WRITE(*,*) 'Subgrid scale orography zpic:', minval(zpic), maxval(zpic)
323
324IF (startphy_file) THEN
325  ! Load sub-grid scale orography parameters:
326  CALL get_field("ZVAL",zval,found)
327  IF (.not.found) THEN
328    WRITE(*,*) modname//": Field <ZVAL> is missing"
329    WRITE(*,*) "set to zero"
330    zval=0.
331  ENDIF
332ELSE
333  zval(:)=0
334ENDIF ! of IF (startphy_file)
335WRITE(*,*) 'Subgrid scale orography zval:', minval(zval), maxval(zval)
336
337IF (startphy_file) THEN
338  ! Lecture de TANCIEN:
339  ancien_ok = .TRUE.
340
341  CALL get_field("TANCIEN",t_ancien,found)
342  IF (.not.found) THEN
343    WRITE(*,*) modname//": Field <TANCIEN> is missing"
344    WRITE(*,*) "Slightly biaised start. But let's keep going."
345    ancien_ok = .FALSE.
346  ENDIF
347ELSE
348  ancien_ok=.false.
349ENDIF
350
351IF (startphy_file) THEN
352  CALL get_field("age",age,found)
353  IF (.not.found) THEN
354    PRINT*, "phyetat0: Age of air is missing"
355    PRINT*, "Reinitialising age of air to 0"
356    age(:,:) = 0.
357  ENDIF
358ELSE
359  age(:,:) = 0.
360ENDIF
361
362IF (startphy_file) THEN
363  ! Load Q2 the TKE at interlayer:
364  CALL get_field("Q2",q2,found)
365  IF (.not.found) THEN
366    WRITE(*,*) modname//": Field <Q2> is missing"
367    WRITE(*,*) "set to zero"
368    q2(:,:)=0.
369  ENDIF
370ELSE
371  ! Dummy initialization
372  q2(:,:)=0
373ENDIF ! of IF (startphy_file)
374WRITE(*,*) 'Turbulent Kinetic Energy', minval(q2), maxval(q2)
375
376! Non-orographic gravity waves
377if (startphy_file) then
378   call get_field("du_nonoro_gwd",du_nonoro_gwd,found)
379   if (.not.found) then
380      write(*,*) "phyetat0: <du_nonoro_gwd> not in file"
381      du_nonoro_gwd(:,:)=0.
382   endif
383endif ! of if (startphy_file)
384if (startphy_file) then
385   call get_field("dv_nonoro_gwd",dv_nonoro_gwd,found)
386   if (.not.found) then
387      write(*,*) "phyetat0: <dv_nonoro_gwd> not in file"
388      dv_nonoro_gwd(:,:)=0.
389   endif
390endif ! of if (startphy_file)
391if (startphy_file) then
392   call get_field("east_gwstress",east_gwstress,found)
393   if (.not.found) then
394      write(*,*) "phyetat0: <east_gwstress> not in file"
395      east_gwstress(:,:)=0.
396   endif
397endif ! of if (startphy_file)
398if (startphy_file) then
399   call get_field("west_gwstress",west_gwstress,found)
400   if (.not.found) then
401      write(*,*) "phyetat0: <west_gwstress> not in file"
402      west_gwstress(:,:)=0.
403   endif
404endif ! of if (startphy_file)
405
406! close file
407IF (startphy_file) call close_startphy
408
409! do some more initializations
410call init_iophy_new(latitude_deg,longitude_deg)
411
412end subroutine phyetat0
413
414end module phyetat0_mod
Note: See TracBrowser for help on using the repository browser.