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

Last change on this file since 3719 was 3719, checked in by emillour, 2 months ago

Venus PCM:
Add reindexing of columns when reading/writing (re)startphy files. This is not
necessary with the lon-lat (LMDZ.COMMON) dynamical core, but required when
using DYNAMICO (where correspondance between dynamics and physics column
indexes changes with number of computing cores).
In addition cleaned up phyetat0: put comments in English, added "only" clauses
to used modules and turned it into a module.
EM

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