1 | module phyetat0_mod |
---|
2 | |
---|
3 | implicit none |
---|
4 | |
---|
5 | contains |
---|
6 | |
---|
7 | subroutine 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 | |
---|
24 | implicit none |
---|
25 | !====================================================================== |
---|
26 | include "dimsoil.h" |
---|
27 | include "clesphys.h" |
---|
28 | include "tabcontrol.h" |
---|
29 | !====================================================================== |
---|
30 | |
---|
31 | character(len=*),intent(in) :: fichnom ! input file name |
---|
32 | LOGICAL :: found |
---|
33 | REAL :: tab_cntrl(length) |
---|
34 | integer :: i,isoil |
---|
35 | CHARACTER(len=2) :: str2 |
---|
36 | REAL :: lon_startphy(klon), lat_startphy(klon) |
---|
37 | REAL :: surface_albedo |
---|
38 | character(len=8) :: modname="phyetat0" |
---|
39 | |
---|
40 | ! global variables read in startphy.nc file |
---|
41 | |
---|
42 | ! open physics initial state file: |
---|
43 | if (startphy_file) then |
---|
44 | call open_startphy(fichnom) |
---|
45 | endif |
---|
46 | |
---|
47 | ! |
---|
48 | ! Load control parameters: |
---|
49 | ! |
---|
50 | IF (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 | |
---|
72 | ELSE |
---|
73 | tabcntr0(:)=1 ! dummy initialization |
---|
74 | ! Initialize parameter or get values from def files |
---|
75 | dtime=pdtphys |
---|
76 | radpas=1 |
---|
77 | itau_phy=0 |
---|
78 | ENDIF ! of IF (startphy_file) |
---|
79 | |
---|
80 | IF (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 |
---|
110 | ENDIF ! of IF (startphy_file) |
---|
111 | |
---|
112 | ! read in other variables here ... |
---|
113 | |
---|
114 | IF (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 |
---|
124 | ELSE |
---|
125 | ! Dummy initialization, but in fact this is later handled in physiq |
---|
126 | ftsol(:)=0 |
---|
127 | ENDIF ! of IF (startphy_file) |
---|
128 | |
---|
129 | IF (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 |
---|
146 | ELSE |
---|
147 | ! Dummy initialization, but in fact this is later handled in physiq |
---|
148 | ftsoil(:,:)=0 |
---|
149 | ENDIF ! of IF (startphy_file) |
---|
150 | |
---|
151 | IF (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 |
---|
158 | ELSE |
---|
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 |
---|
163 | ENDIF ! of IF (startphy_file) |
---|
164 | WRITE(*,*) 'Surface albedo <ALBE>', minval(falbe), maxval(falbe) |
---|
165 | |
---|
166 | IF (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 |
---|
174 | ELSE |
---|
175 | ! Dummy initialization |
---|
176 | solsw(:)=0 |
---|
177 | ENDIF ! of IF (startphy_file) |
---|
178 | WRITE(*,*) 'Solar flux to the surface solsw:', minval(solsw), maxval(solsw) |
---|
179 | |
---|
180 | IF (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 |
---|
188 | ELSE |
---|
189 | ! Dummy initialization |
---|
190 | sollw(:)=0 |
---|
191 | ENDIF ! of IF (startphy_file) |
---|
192 | WRITE(*,*) 'IR flux to the surface sollw:', minval(sollw), maxval(solsw) |
---|
193 | |
---|
194 | IF (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 |
---|
202 | ELSE |
---|
203 | ! Dummy initialization |
---|
204 | fder(:)=0 |
---|
205 | ENDIF ! of IF (startphy_file) |
---|
206 | WRITE(*,*) 'Derivative of the flux fder:', minval(fder), maxval(fder) |
---|
207 | |
---|
208 | IF (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 |
---|
216 | ELSE |
---|
217 | ! Dummy initialization |
---|
218 | dlw(:)=0 |
---|
219 | ENDIF ! of IF (startphy_file) |
---|
220 | WRITE(*,*) 'Derivative of the IR flux dlw:', minval(dlw), maxval(dlw) |
---|
221 | |
---|
222 | IF (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 |
---|
230 | ELSE |
---|
231 | ! Dummy initialization |
---|
232 | sollwdown(:)=0 |
---|
233 | ENDIF ! of IF (startphy_file) |
---|
234 | WRITE(*,*) 'Downward IR flux at the surface sollwdown:', minval(sollwdown), maxval(sollwdown) |
---|
235 | |
---|
236 | IF (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 |
---|
243 | ELSE |
---|
244 | ! Dummy initialization |
---|
245 | radsol(:)=0 |
---|
246 | ENDIF ! of IF (startphy_file) |
---|
247 | WRITE(*,*) 'Net flux at the surface radsol:', minval(radsol), maxval(radsol) |
---|
248 | |
---|
249 | IF (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 |
---|
257 | ELSE |
---|
258 | zmea(:)=0 |
---|
259 | ENDIF ! of IF (startphy_file) |
---|
260 | WRITE(*,*) 'Subgrid scale orography zmea:', minval(zmea), maxval(zmea) |
---|
261 | |
---|
262 | IF (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 |
---|
270 | ELSE |
---|
271 | zstd(:)=0 |
---|
272 | ENDIF ! of IF (startphy_file) |
---|
273 | WRITE(*,*) 'Subgrid scale orography zstd:', minval(zstd), maxval(zstd) |
---|
274 | |
---|
275 | IF (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 |
---|
283 | ELSE |
---|
284 | zsig(:)=0 |
---|
285 | ENDIF ! of IF (startphy_file) |
---|
286 | WRITE(*,*) 'Subgrid scale orography zsig:', minval(zsig), maxval(zsig) |
---|
287 | |
---|
288 | IF (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 |
---|
296 | ELSE |
---|
297 | zgam(:)=0 |
---|
298 | ENDIF ! of IF (startphy_file) |
---|
299 | WRITE(*,*) 'Subgrid scale orography zgam:', minval(zgam), maxval(zgam) |
---|
300 | |
---|
301 | IF (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 |
---|
309 | ELSE |
---|
310 | zthe(:)=0 |
---|
311 | ENDIF ! of IF (startphy_file) |
---|
312 | WRITE(*,*) 'Subgrid scale orography zthe:', minval(zthe), maxval(zthe) |
---|
313 | |
---|
314 | IF (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 |
---|
322 | ELSE |
---|
323 | zpic(:)=0 |
---|
324 | ENDIF ! of IF (startphy_file) |
---|
325 | WRITE(*,*) 'Subgrid scale orography zpic:', minval(zpic), maxval(zpic) |
---|
326 | |
---|
327 | IF (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 |
---|
335 | ELSE |
---|
336 | zval(:)=0 |
---|
337 | ENDIF ! of IF (startphy_file) |
---|
338 | WRITE(*,*) 'Subgrid scale orography zval:', minval(zval), maxval(zval) |
---|
339 | |
---|
340 | IF (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 |
---|
350 | ELSE |
---|
351 | ancien_ok=.false. |
---|
352 | ENDIF |
---|
353 | |
---|
354 | IF (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 |
---|
361 | ELSE |
---|
362 | age(:,:) = 0. |
---|
363 | ENDIF |
---|
364 | |
---|
365 | IF (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 |
---|
373 | ELSE |
---|
374 | ! Dummy initialization |
---|
375 | q2(:,:)=0 |
---|
376 | ENDIF ! of IF (startphy_file) |
---|
377 | WRITE(*,*) 'Turbulent Kinetic Energy', minval(q2), maxval(q2) |
---|
378 | |
---|
379 | ! Non-orographic gravity waves |
---|
380 | if (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 |
---|
386 | endif ! of if (startphy_file) |
---|
387 | if (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 |
---|
393 | endif ! of if (startphy_file) |
---|
394 | if (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 |
---|
400 | endif ! of if (startphy_file) |
---|
401 | if (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 |
---|
407 | endif ! of if (startphy_file) |
---|
408 | |
---|
409 | ! close file |
---|
410 | IF (startphy_file) call close_startphy |
---|
411 | |
---|
412 | ! do some more initializations |
---|
413 | call init_iophy_new(latitude_deg,longitude_deg) |
---|
414 | |
---|
415 | end subroutine phyetat0 |
---|
416 | |
---|
417 | end module phyetat0_mod |
---|