1 | ! |
---|
2 | ! $Header$ |
---|
3 | ! |
---|
4 | C |
---|
5 | C |
---|
6 | MODULE startvar |
---|
7 | ! |
---|
8 | ! |
---|
9 | ! There are three ways to access data from the database of atmospheric data which |
---|
10 | ! can be used to initialize the model. This depends on the type of field which needs |
---|
11 | ! to be extracted. In any case the call should come after a restget and should be of the type : |
---|
12 | ! CALL startget(...) |
---|
13 | ! |
---|
14 | ! We will details the possible arguments to startget here : |
---|
15 | ! |
---|
16 | ! - A 2D variable on the dynamical grid : |
---|
17 | ! CALL startget(varname, iml, jml, lon_in, lat_in, champ, val_ex, jml2, lon_in2, lat_in2, interbar ) |
---|
18 | ! |
---|
19 | ! - A 1D variable on the physical grid : |
---|
20 | ! CALL startget(varname, iml, jml, lon_in, lat_in, nbindex, champ, val_exp, jml2, lon_in2, lat_in2, interbar ) |
---|
21 | ! |
---|
22 | ! |
---|
23 | ! - A 3D variable on the dynamical grid : |
---|
24 | ! CALL startget(varname, iml, jml, lon_in, lat_in, lml, pls, workvar, champ, val_exp, jml2, lon_in2, lat_in2, interbar ) |
---|
25 | ! |
---|
26 | ! |
---|
27 | ! There is special constraint on the atmospheric data base except that the |
---|
28 | ! the data needs to be in netCDF and the variables should have the the following |
---|
29 | ! names in the file : |
---|
30 | ! |
---|
31 | ! 'RELIEF' : High resolution orography |
---|
32 | ! 'ST' : Surface temperature |
---|
33 | ! 'CDSW' : Soil moisture |
---|
34 | ! 'Z' : Surface geopotential |
---|
35 | ! 'SP' : Surface pressure |
---|
36 | ! 'U' : East ward wind |
---|
37 | ! 'V' : Northward wind |
---|
38 | ! 'TEMP' : Temperature |
---|
39 | ! 'R' : Relative humidity |
---|
40 | ! |
---|
41 | USE ioipsl |
---|
42 | ! |
---|
43 | ! |
---|
44 | IMPLICIT NONE |
---|
45 | ! |
---|
46 | ! |
---|
47 | PRIVATE |
---|
48 | PUBLIC startget |
---|
49 | ! |
---|
50 | ! |
---|
51 | INTERFACE startget |
---|
52 | MODULE PROCEDURE startget_phys2d, startget_phys1d, startget_dyn |
---|
53 | END INTERFACE |
---|
54 | ! |
---|
55 | INTEGER, SAVE :: fid_phys, fid_dyn |
---|
56 | INTEGER, SAVE :: iml_phys, iml_rel, iml_dyn |
---|
57 | INTEGER, SAVE :: jml_phys, jml_rel, jml_dyn |
---|
58 | INTEGER, SAVE :: llm_dyn, ttm_dyn |
---|
59 | REAL, ALLOCATABLE, SAVE, DIMENSION (:,:) :: lon_phys, lon_rug, |
---|
60 | . lon_alb, lon_rel, lon_dyn |
---|
61 | REAL, ALLOCATABLE, SAVE, DIMENSION (:,:) :: lat_phys, lat_rug, |
---|
62 | . lat_alb, lat_rel, lat_dyn |
---|
63 | REAL, ALLOCATABLE, SAVE, DIMENSION (:) :: levdyn_ini |
---|
64 | REAL, ALLOCATABLE, SAVE, DIMENSION (:,:) :: relief, zstd, zsig, |
---|
65 | . zgam, zthe, zpic, zval |
---|
66 | REAL, ALLOCATABLE, SAVE, DIMENSION (:,:) :: rugo, masque, phis |
---|
67 | ! |
---|
68 | REAL, ALLOCATABLE, SAVE, DIMENSION (:,:) :: tsol, qsol, psol_dyn |
---|
69 | REAL, ALLOCATABLE, SAVE, DIMENSION (:,:,:) :: var_ana3d |
---|
70 | ! |
---|
71 | CONTAINS |
---|
72 | ! |
---|
73 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
74 | ! |
---|
75 | SUBROUTINE startget_phys2d(varname, iml, jml, lon_in, lat_in, |
---|
76 | . champ, val_exp, jml2, lon_in2, lat_in2 , interbar, masque_lu ) |
---|
77 | ! |
---|
78 | ! There is a big mess with the size in logitude, should it be iml or iml+1. |
---|
79 | ! I have chosen to use the iml+1 as an argument to this routine and we declare |
---|
80 | ! internaly smaler fields when needed. This needs to be cleared once and for all in LMDZ. |
---|
81 | ! A convention is required. |
---|
82 | ! |
---|
83 | ! |
---|
84 | CHARACTER*(*), INTENT(in) :: varname |
---|
85 | INTEGER, INTENT(in) :: iml, jml ,jml2 |
---|
86 | REAL, INTENT(in) :: lon_in(iml), lat_in(jml) |
---|
87 | REAL, INTENT(in) :: lon_in2(iml), lat_in2(jml2) |
---|
88 | REAL, INTENT(inout) :: champ(iml,jml) |
---|
89 | REAL, INTENT(in) :: val_exp |
---|
90 | REAL, INTENT(in), optional :: masque_lu(iml,jml) |
---|
91 | LOGICAL interbar |
---|
92 | ! |
---|
93 | ! This routine only works if the variable does not exist or is constant |
---|
94 | ! |
---|
95 | IF ( MINVAL(champ(:,:)).EQ.MAXVAL(champ(:,:)) .AND. |
---|
96 | .MINVAL(champ(:,:)).EQ.val_exp ) THEN |
---|
97 | ! |
---|
98 | SELECTCASE(varname) |
---|
99 | ! |
---|
100 | CASE ('relief') |
---|
101 | ! |
---|
102 | ! If we do not have the orography we need to get it |
---|
103 | ! |
---|
104 | IF ( .NOT.ALLOCATED(relief)) THEN |
---|
105 | ! |
---|
106 | if (present(masque_lu)) then |
---|
107 | CALL start_init_orog( iml, jml, lon_in, lat_in, |
---|
108 | . jml2,lon_in2,lat_in2, interbar, masque_lu ) |
---|
109 | else |
---|
110 | CALL start_init_orog( iml, jml, lon_in, lat_in, |
---|
111 | . jml2,lon_in2,lat_in2, interbar) |
---|
112 | endif |
---|
113 | ! |
---|
114 | ENDIF |
---|
115 | ! |
---|
116 | IF (SIZE(relief) .NE. SIZE(champ)) THEN |
---|
117 | ! |
---|
118 | WRITE(*,*) 'STARTVAR module has been', |
---|
119 | .' initialized to the wrong size' |
---|
120 | STOP |
---|
121 | ! |
---|
122 | ENDIF |
---|
123 | ! |
---|
124 | champ(:,:) = relief(:,:) |
---|
125 | ! |
---|
126 | CASE ('rugosite') |
---|
127 | ! |
---|
128 | ! If we do not have the orography we need to get it |
---|
129 | ! |
---|
130 | IF ( .NOT.ALLOCATED(rugo)) THEN |
---|
131 | ! |
---|
132 | CALL start_init_orog( iml, jml, lon_in, lat_in, |
---|
133 | . jml2,lon_in2,lat_in2 , interbar ) |
---|
134 | ! |
---|
135 | ENDIF |
---|
136 | ! |
---|
137 | IF (SIZE(rugo) .NE. SIZE(champ)) THEN |
---|
138 | ! |
---|
139 | WRITE(*,*) |
---|
140 | . 'STARTVAR module has been initialized to the wrong size' |
---|
141 | STOP |
---|
142 | ! |
---|
143 | ENDIF |
---|
144 | ! |
---|
145 | champ(:,:) = rugo(:,:) |
---|
146 | ! |
---|
147 | CASE ('masque') |
---|
148 | ! |
---|
149 | ! If we do not have the orography we need to get it |
---|
150 | ! |
---|
151 | IF ( .NOT.ALLOCATED(masque)) THEN |
---|
152 | ! |
---|
153 | CALL start_init_orog( iml, jml, lon_in, lat_in, |
---|
154 | . jml2,lon_in2,lat_in2 , interbar ) |
---|
155 | ! |
---|
156 | ENDIF |
---|
157 | ! |
---|
158 | IF (SIZE(masque) .NE. SIZE(champ)) THEN |
---|
159 | ! |
---|
160 | WRITE(*,*) |
---|
161 | . 'STARTVAR module has been initialized to the wrong size' |
---|
162 | STOP |
---|
163 | ! |
---|
164 | ENDIF |
---|
165 | ! |
---|
166 | champ(:,:) = masque(:,:) |
---|
167 | ! |
---|
168 | CASE ('surfgeo') |
---|
169 | ! |
---|
170 | ! If we do not have the orography we need to get it |
---|
171 | ! |
---|
172 | IF ( .NOT.ALLOCATED(phis)) THEN |
---|
173 | ! |
---|
174 | CALL start_init_orog( iml, jml, lon_in, lat_in, |
---|
175 | . jml2,lon_in2, lat_in2 , interbar ) |
---|
176 | ! |
---|
177 | ENDIF |
---|
178 | ! |
---|
179 | IF (SIZE(phis) .NE. SIZE(champ)) THEN |
---|
180 | ! |
---|
181 | WRITE(*,*) |
---|
182 | . 'STARTVAR module has been initialized to the wrong size' |
---|
183 | STOP |
---|
184 | ! |
---|
185 | ENDIF |
---|
186 | ! |
---|
187 | champ(:,:) = phis(:,:) |
---|
188 | ! |
---|
189 | CASE ('psol') |
---|
190 | ! |
---|
191 | ! If we do not have the orography we need to get it |
---|
192 | ! |
---|
193 | IF ( .NOT.ALLOCATED(psol_dyn)) THEN |
---|
194 | ! |
---|
195 | CALL start_init_dyn( iml, jml, lon_in, lat_in, |
---|
196 | . jml2,lon_in2, lat_in2 , interbar ) |
---|
197 | ! |
---|
198 | ENDIF |
---|
199 | ! |
---|
200 | IF (SIZE(psol_dyn) .NE. SIZE(champ)) THEN |
---|
201 | ! |
---|
202 | WRITE(*,*) |
---|
203 | . 'STARTVAR module has been initialized to the wrong size' |
---|
204 | STOP |
---|
205 | ! |
---|
206 | ENDIF |
---|
207 | ! |
---|
208 | champ(:,:) = psol_dyn(:,:) |
---|
209 | ! |
---|
210 | CASE DEFAULT |
---|
211 | ! |
---|
212 | WRITE(*,*) 'startget_phys2d' |
---|
213 | WRITE(*,*) 'No rule is present to extract variable', |
---|
214 | . varname(:LEN_TRIM(varname)),' from any data set' |
---|
215 | STOP |
---|
216 | ! |
---|
217 | END SELECT |
---|
218 | ! |
---|
219 | ELSE |
---|
220 | ! |
---|
221 | ! There are a few fields we might need if we need to interpolate 3D filed. Thus if they come through here we |
---|
222 | ! will catch them |
---|
223 | ! |
---|
224 | SELECTCASE(varname) |
---|
225 | ! |
---|
226 | CASE ('surfgeo') |
---|
227 | ! |
---|
228 | IF ( .NOT.ALLOCATED(phis)) THEN |
---|
229 | ALLOCATE(phis(iml,jml)) |
---|
230 | ENDIF |
---|
231 | ! |
---|
232 | IF (SIZE(phis) .NE. SIZE(champ)) THEN |
---|
233 | ! |
---|
234 | WRITE(*,*) |
---|
235 | . 'STARTVAR module has been initialized to the wrong size' |
---|
236 | STOP |
---|
237 | ! |
---|
238 | ENDIF |
---|
239 | ! |
---|
240 | phis(:,:) = champ(:,:) |
---|
241 | ! |
---|
242 | END SELECT |
---|
243 | ! |
---|
244 | ENDIF |
---|
245 | ! |
---|
246 | END SUBROUTINE startget_phys2d |
---|
247 | ! |
---|
248 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
249 | ! |
---|
250 | SUBROUTINE start_init_orog ( iml,jml,lon_in, lat_in,jml2,lon_in2 , |
---|
251 | , lat_in2 , interbar, masque_lu ) |
---|
252 | ! |
---|
253 | INTEGER, INTENT(in) :: iml, jml, jml2 |
---|
254 | REAL, INTENT(in) :: lon_in(iml), lat_in(jml) |
---|
255 | REAL, INTENT(in) :: lon_in2(iml), lat_in2(jml2) |
---|
256 | REAL, intent(in), optional :: masque_lu(iml,jml) |
---|
257 | LOGICAL interbar |
---|
258 | ! |
---|
259 | ! LOCAL |
---|
260 | ! |
---|
261 | LOGICAL interbar2 |
---|
262 | REAL :: lev(1), date, dt,chmin,chmax |
---|
263 | INTEGER :: itau(1), fid |
---|
264 | INTEGER :: llm_tmp, ttm_tmp |
---|
265 | INTEGER :: i, j |
---|
266 | INTEGER :: iret |
---|
267 | CHARACTER*25 title |
---|
268 | REAL, ALLOCATABLE :: relief_hi(:,:) |
---|
269 | REAL, ALLOCATABLE :: lon_rad(:), lat_rad(:) |
---|
270 | REAL, ALLOCATABLE :: lon_ini(:), lat_ini(:) |
---|
271 | REAL, ALLOCATABLE :: tmp_var(:,:) |
---|
272 | INTEGER, ALLOCATABLE :: tmp_int(:,:) |
---|
273 | ! |
---|
274 | CHARACTER*120 :: orogfname |
---|
275 | LOGICAL :: check=.TRUE. |
---|
276 | ! |
---|
277 | ! |
---|
278 | orogfname = 'Relief.nc' |
---|
279 | ! |
---|
280 | IF ( check ) WRITE(*,*) 'Reading the high resolution orography' |
---|
281 | ! |
---|
282 | CALL flininfo(orogfname,iml_rel, jml_rel, llm_tmp, ttm_tmp, fid) |
---|
283 | ! |
---|
284 | ALLOCATE (lat_rel(iml_rel,jml_rel), stat=iret) |
---|
285 | ALLOCATE (lon_rel(iml_rel,jml_rel), stat=iret) |
---|
286 | ALLOCATE (relief_hi(iml_rel,jml_rel), stat=iret) |
---|
287 | ! |
---|
288 | CALL flinopen(orogfname, .FALSE., iml_rel, jml_rel, |
---|
289 | .llm_tmp, lon_rel, lat_rel, lev, ttm_tmp, |
---|
290 | . itau, date, dt, fid) |
---|
291 | ! |
---|
292 | CALL flinget(fid, 'RELIEF', iml_rel, jml_rel, llm_tmp, |
---|
293 | . ttm_tmp, 1, 1, relief_hi) |
---|
294 | ! |
---|
295 | CALL flinclo(fid) |
---|
296 | ! |
---|
297 | ! In case we have a file which is in degrees we do the transformation |
---|
298 | ! |
---|
299 | ALLOCATE(lon_rad(iml_rel)) |
---|
300 | ALLOCATE(lon_ini(iml_rel)) |
---|
301 | |
---|
302 | IF ( MAXVAL(lon_rel(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN |
---|
303 | lon_ini(:) = lon_rel(:,1) * 2.0 * ASIN(1.0) / 180.0 |
---|
304 | ELSE |
---|
305 | lon_ini(:) = lon_rel(:,1) |
---|
306 | ENDIF |
---|
307 | |
---|
308 | ALLOCATE(lat_rad(jml_rel)) |
---|
309 | ALLOCATE(lat_ini(jml_rel)) |
---|
310 | |
---|
311 | IF ( MAXVAL(lat_rel(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN |
---|
312 | lat_ini(:) = lat_rel(1,:) * 2.0 * ASIN(1.0) / 180.0 |
---|
313 | ELSE |
---|
314 | lat_ini(:) = lat_rel(1,:) |
---|
315 | ENDIF |
---|
316 | ! |
---|
317 | ! |
---|
318 | |
---|
319 | title='RELIEF' |
---|
320 | |
---|
321 | interbar2 = .FALSE. |
---|
322 | CALL conf_dat2d(title,iml_rel, jml_rel, lon_ini, lat_ini, |
---|
323 | . lon_rad, lat_rad, relief_hi , interbar2 ) |
---|
324 | |
---|
325 | IF ( check ) WRITE(*,*) 'Computes all the parameters needed', |
---|
326 | .' for the gravity wave drag code' |
---|
327 | ! |
---|
328 | ! Allocate the data we need to put in the interpolated fields |
---|
329 | ! |
---|
330 | ! RELIEF: orographie moyenne |
---|
331 | ALLOCATE(relief(iml,jml)) |
---|
332 | ! zphi : orographie moyenne |
---|
333 | ALLOCATE(phis(iml,jml)) |
---|
334 | ! zstd: deviation standard de l'orographie sous-maille |
---|
335 | ALLOCATE(zstd(iml,jml)) |
---|
336 | ! zsig: pente de l'orographie sous-maille |
---|
337 | ALLOCATE(zsig(iml,jml)) |
---|
338 | ! zgam: anisotropy de l'orographie sous maille |
---|
339 | ALLOCATE(zgam(iml,jml)) |
---|
340 | ! zthe: orientation de l'axe oriente dans la direction |
---|
341 | ! de plus grande pente de l'orographie sous maille |
---|
342 | ALLOCATE(zthe(iml,jml)) |
---|
343 | ! zpic: hauteur pics de la SSO |
---|
344 | ALLOCATE(zpic(iml,jml)) |
---|
345 | ! zval: hauteur vallees de la SSO |
---|
346 | ALLOCATE(zval(iml,jml)) |
---|
347 | ! masque : Masque terre ocean |
---|
348 | ALLOCATE(tmp_int(iml,jml)) |
---|
349 | ALLOCATE(masque(iml,jml)) |
---|
350 | |
---|
351 | masque = -99999. |
---|
352 | if (present(masque_lu)) then |
---|
353 | masque = masque_lu |
---|
354 | endif |
---|
355 | ! |
---|
356 | CALL grid_noro(iml_rel, jml_rel, lon_rad, lat_rad, relief_hi, |
---|
357 | . iml-1, jml, lon_in, lat_in, |
---|
358 | . phis, relief, zstd, zsig, zgam, zthe, zpic, zval, masque) |
---|
359 | phis = phis * 9.81 |
---|
360 | ! |
---|
361 | ! masque(:,:) = FLOAT(tmp_int(:,:)) |
---|
362 | ! |
---|
363 | ! Compute surface roughness |
---|
364 | ! |
---|
365 | IF ( check ) WRITE(*,*) |
---|
366 | .'Compute surface roughness induced by the orography' |
---|
367 | ! |
---|
368 | ALLOCATE(rugo(iml,jml)) |
---|
369 | ALLOCATE(tmp_var(iml-1,jml)) |
---|
370 | ! |
---|
371 | CALL rugsoro(iml_rel, jml_rel, lon_rad, lat_rad, relief_hi, |
---|
372 | . iml-1, jml, lon_in, lat_in, tmp_var) |
---|
373 | ! |
---|
374 | DO j = 1, jml |
---|
375 | DO i = 1, iml-1 |
---|
376 | rugo(i,j) = tmp_var(i,j) |
---|
377 | ENDDO |
---|
378 | rugo(iml,j) = tmp_var(1,j) |
---|
379 | ENDDO |
---|
380 | c |
---|
381 | cc *** rugo n'est pas utilise pour l'instant ****** |
---|
382 | ! |
---|
383 | ! Build land-sea mask |
---|
384 | ! |
---|
385 | ! |
---|
386 | RETURN |
---|
387 | ! |
---|
388 | END SUBROUTINE start_init_orog |
---|
389 | ! |
---|
390 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
391 | ! |
---|
392 | SUBROUTINE startget_phys1d(varname, iml, jml, lon_in, |
---|
393 | .lat_in, nbindex, champ, val_exp ,jml2, lon_in2, lat_in2,interbar) |
---|
394 | ! |
---|
395 | CHARACTER*(*), INTENT(in) :: varname |
---|
396 | INTEGER, INTENT(in) :: iml, jml, nbindex, jml2 |
---|
397 | REAL, INTENT(in) :: lon_in(iml), lat_in(jml) |
---|
398 | REAL, INTENT(in) :: lon_in2(iml), lat_in2(jml2) |
---|
399 | REAL, INTENT(inout) :: champ(nbindex) |
---|
400 | REAL, INTENT(in) :: val_exp |
---|
401 | LOGICAL interbar |
---|
402 | ! |
---|
403 | ! |
---|
404 | ! This routine only works if the variable does not exist or is constant |
---|
405 | ! |
---|
406 | IF ( MINVAL(champ(:)).EQ.MAXVAL(champ(:)) .AND. |
---|
407 | .MINVAL(champ(:)).EQ.val_exp ) THEN |
---|
408 | SELECTCASE(varname) |
---|
409 | CASE ('tsol') |
---|
410 | IF ( .NOT.ALLOCATED(tsol)) THEN |
---|
411 | CALL start_init_phys( iml, jml, lon_in, lat_in, |
---|
412 | . jml2, lon_in2, lat_in2, interbar ) |
---|
413 | ENDIF |
---|
414 | IF ( SIZE(tsol) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN |
---|
415 | WRITE(*,*) |
---|
416 | . 'STARTVAR module has been initialized to the wrong size' |
---|
417 | STOP |
---|
418 | ENDIF |
---|
419 | CALL gr_dyn_fi(1, iml, jml, nbindex, tsol, champ) |
---|
420 | CASE ('qsol') |
---|
421 | IF ( .NOT.ALLOCATED(qsol)) THEN |
---|
422 | CALL start_init_phys( iml, jml, lon_in, lat_in, |
---|
423 | . jml2, lon_in2,lat_in2 , interbar ) |
---|
424 | ENDIF |
---|
425 | IF ( SIZE(qsol) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN |
---|
426 | WRITE(*,*) |
---|
427 | . 'STARTVAR module has been initialized to the wrong size' |
---|
428 | STOP |
---|
429 | ENDIF |
---|
430 | CALL gr_dyn_fi(1, iml, jml, nbindex, qsol, champ) |
---|
431 | CASE ('psol') |
---|
432 | IF ( .NOT.ALLOCATED(psol_dyn)) THEN |
---|
433 | CALL start_init_dyn( iml, jml, lon_in, lat_in, |
---|
434 | . jml2, lon_in2,lat_in2 , interbar ) |
---|
435 | ENDIF |
---|
436 | IF (SIZE(psol_dyn) .NE. SIZE(lon_in)*SIZE(lat_in)) THEN |
---|
437 | WRITE(*,*) |
---|
438 | . 'STARTVAR module has been initialized to the wrong size' |
---|
439 | STOP |
---|
440 | ENDIF |
---|
441 | CALL gr_dyn_fi(1, iml, jml, nbindex, psol_dyn, champ) |
---|
442 | CASE ('zmea') |
---|
443 | IF ( .NOT.ALLOCATED(relief)) THEN |
---|
444 | CALL start_init_orog( iml, jml, lon_in, lat_in, |
---|
445 | . jml2, lon_in2,lat_in2 , interbar ) |
---|
446 | ENDIF |
---|
447 | IF ( SIZE(relief) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN |
---|
448 | WRITE(*,*) |
---|
449 | . 'STARTVAR module has been initialized to the wrong size' |
---|
450 | STOP |
---|
451 | ENDIF |
---|
452 | CALL gr_dyn_fi(1, iml, jml, nbindex, relief, champ) |
---|
453 | CASE ('zstd') |
---|
454 | IF ( .NOT.ALLOCATED(zstd)) THEN |
---|
455 | CALL start_init_orog( iml, jml, lon_in, lat_in, |
---|
456 | . jml2, lon_in2,lat_in2 , interbar ) |
---|
457 | ENDIF |
---|
458 | IF ( SIZE(zstd) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN |
---|
459 | WRITE(*,*) |
---|
460 | . 'STARTVAR module has been initialized to the wrong size' |
---|
461 | STOP |
---|
462 | ENDIF |
---|
463 | CALL gr_dyn_fi(1, iml, jml, nbindex,zstd, champ) |
---|
464 | CASE ('zsig') |
---|
465 | IF ( .NOT.ALLOCATED(zsig)) THEN |
---|
466 | CALL start_init_orog( iml, jml, lon_in, lat_in, |
---|
467 | . jml2, lon_in2,lat_in2 , interbar ) |
---|
468 | ENDIF |
---|
469 | IF ( SIZE(zsig) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN |
---|
470 | WRITE(*,*) |
---|
471 | . 'STARTVAR module has been initialized to the wrong size' |
---|
472 | STOP |
---|
473 | ENDIF |
---|
474 | CALL gr_dyn_fi(1, iml, jml, nbindex,zsig, champ) |
---|
475 | CASE ('zgam') |
---|
476 | IF ( .NOT.ALLOCATED(zgam)) THEN |
---|
477 | CALL start_init_orog( iml, jml, lon_in, lat_in, |
---|
478 | . jml2, lon_in2,lat_in2 , interbar ) |
---|
479 | ENDIF |
---|
480 | IF ( SIZE(zgam) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN |
---|
481 | WRITE(*,*) |
---|
482 | . 'STARTVAR module has been initialized to the wrong size' |
---|
483 | STOP |
---|
484 | ENDIF |
---|
485 | CALL gr_dyn_fi(1, iml, jml, nbindex,zgam, champ) |
---|
486 | CASE ('zthe') |
---|
487 | IF ( .NOT.ALLOCATED(zthe)) THEN |
---|
488 | CALL start_init_orog( iml, jml, lon_in, lat_in, |
---|
489 | . jml2, lon_in2,lat_in2 , interbar ) |
---|
490 | ENDIF |
---|
491 | IF ( SIZE(zthe) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN |
---|
492 | WRITE(*,*) |
---|
493 | . 'STARTVAR module has been initialized to the wrong size' |
---|
494 | STOP |
---|
495 | ENDIF |
---|
496 | CALL gr_dyn_fi(1, iml, jml, nbindex,zthe, champ) |
---|
497 | CASE ('zpic') |
---|
498 | IF ( .NOT.ALLOCATED(zpic)) THEN |
---|
499 | CALL start_init_orog( iml, jml, lon_in, lat_in, |
---|
500 | . jml2, lon_in2,lat_in2 , interbar ) |
---|
501 | ENDIF |
---|
502 | IF ( SIZE(zpic) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN |
---|
503 | WRITE(*,*) |
---|
504 | . 'STARTVAR module has been initialized to the wrong size' |
---|
505 | STOP |
---|
506 | ENDIF |
---|
507 | CALL gr_dyn_fi(1, iml, jml, nbindex,zpic, champ) |
---|
508 | CASE ('zval') |
---|
509 | IF ( .NOT.ALLOCATED(zval)) THEN |
---|
510 | CALL start_init_orog( iml, jml, lon_in, lat_in, |
---|
511 | . jml2, lon_in2,lat_in2 , interbar ) |
---|
512 | ENDIF |
---|
513 | IF ( SIZE(zval) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN |
---|
514 | WRITE(*,*) |
---|
515 | . 'STARTVAR module has been initialized to the wrong size' |
---|
516 | STOP |
---|
517 | ENDIF |
---|
518 | CALL gr_dyn_fi(1, iml, jml, nbindex,zval, champ) |
---|
519 | CASE ('rads') |
---|
520 | champ(:) = 0.0 |
---|
521 | CASE ('snow') |
---|
522 | champ(:) = 0.0 |
---|
523 | CASE ('deltat') |
---|
524 | champ(:) = 0.0 |
---|
525 | CASE ('rugmer') |
---|
526 | champ(:) = 0.001 |
---|
527 | CASE ('agsno') |
---|
528 | champ(:) = 50.0 |
---|
529 | CASE DEFAULT |
---|
530 | WRITE(*,*) 'startget_phys1d' |
---|
531 | WRITE(*,*) 'No rule is present to extract variable ', |
---|
532 | . varname(:LEN_TRIM(varname)),' from any data set' |
---|
533 | STOP |
---|
534 | END SELECT |
---|
535 | ELSE |
---|
536 | ! |
---|
537 | ! If we see tsol we catch it as we may need it for a 3D interpolation |
---|
538 | ! |
---|
539 | SELECTCASE(varname) |
---|
540 | CASE ('tsol') |
---|
541 | IF ( .NOT.ALLOCATED(tsol)) THEN |
---|
542 | ALLOCATE(tsol(SIZE(lon_in),SIZE(lat_in) )) |
---|
543 | ENDIF |
---|
544 | CALL gr_fi_dyn(1, iml, jml, nbindex, champ, tsol) |
---|
545 | END SELECT |
---|
546 | ENDIF |
---|
547 | END SUBROUTINE startget_phys1d |
---|
548 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
549 | ! |
---|
550 | SUBROUTINE start_init_phys( iml, jml, lon_in, lat_in, jml2, |
---|
551 | . lon_in2, lat_in2 , interbar ) |
---|
552 | ! |
---|
553 | INTEGER, INTENT(in) :: iml, jml ,jml2 |
---|
554 | REAL, INTENT(in) :: lon_in(iml), lat_in(jml) |
---|
555 | REAL, INTENT(in) :: lon_in2(iml), lat_in2(jml2) |
---|
556 | LOGICAL interbar |
---|
557 | ! |
---|
558 | ! LOCAL |
---|
559 | ! |
---|
560 | !ac REAL :: lev(1), date, dt |
---|
561 | REAL :: date, dt |
---|
562 | REAL, DIMENSION(:), ALLOCATABLE :: levphys_ini |
---|
563 | !ac |
---|
564 | INTEGER :: itau(1) |
---|
565 | INTEGER :: llm_tmp, ttm_tmp |
---|
566 | INTEGER :: i, j |
---|
567 | ! |
---|
568 | CHARACTER*25 title |
---|
569 | CHARACTER*120 :: physfname |
---|
570 | LOGICAL :: check=.TRUE. |
---|
571 | ! |
---|
572 | REAL, ALLOCATABLE :: lon_rad(:), lat_rad(:) |
---|
573 | REAL, ALLOCATABLE :: lon_ini(:), lat_ini(:) |
---|
574 | REAL, ALLOCATABLE :: var_ana(:,:), tmp_var(:,:) |
---|
575 | ! |
---|
576 | physfname = 'ECPHY.nc' |
---|
577 | ! |
---|
578 | IF ( check ) WRITE(*,*) 'Opening the surface analysis' |
---|
579 | ! |
---|
580 | CALL flininfo(physfname, iml_phys, jml_phys, llm_tmp, |
---|
581 | . ttm_tmp, fid_phys) |
---|
582 | ! |
---|
583 | ALLOCATE (lat_phys(iml_phys,jml_phys)) |
---|
584 | ALLOCATE (lon_phys(iml_phys,jml_phys)) |
---|
585 | !ac |
---|
586 | ALLOCATE (levphys_ini(llm_tmp)) |
---|
587 | ! |
---|
588 | ! CALL flinopen(physfname, .FALSE., iml_phys, jml_phys, |
---|
589 | ! . llm_tmp, lon_phys, lat_phys, lev, ttm_tmp, |
---|
590 | ! . itau, date, dt, fid_phys) |
---|
591 | ! |
---|
592 | CALL flinopen(physfname, .FALSE., iml_phys, jml_phys, |
---|
593 | . llm_tmp, lon_phys, lat_phys, levphys_ini, ttm_tmp, |
---|
594 | . itau, date, dt, fid_phys) |
---|
595 | ! |
---|
596 | DEALLOCATE (levphys_ini) |
---|
597 | !ac |
---|
598 | ! |
---|
599 | ! Allocate the space we will need to get the data out of this file |
---|
600 | ! |
---|
601 | ALLOCATE(var_ana(iml_phys, jml_phys)) |
---|
602 | ! |
---|
603 | ! In case we have a file which is in degrees we do the transformation |
---|
604 | ! |
---|
605 | ALLOCATE(lon_rad(iml_phys)) |
---|
606 | ALLOCATE(lon_ini(iml_phys)) |
---|
607 | |
---|
608 | IF ( MAXVAL(lon_phys(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN |
---|
609 | lon_ini(:) = lon_phys(:,1) * 2.0 * ASIN(1.0) / 180.0 |
---|
610 | ELSE |
---|
611 | lon_ini(:) = lon_phys(:,1) |
---|
612 | ENDIF |
---|
613 | |
---|
614 | ALLOCATE(lat_rad(jml_phys)) |
---|
615 | ALLOCATE(lat_ini(jml_phys)) |
---|
616 | |
---|
617 | IF ( MAXVAL(lat_phys(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN |
---|
618 | lat_ini(:) = lat_phys(1,:) * 2.0 * ASIN(1.0) / 180.0 |
---|
619 | ELSE |
---|
620 | lat_ini(:) = lat_phys(1,:) |
---|
621 | ENDIF |
---|
622 | |
---|
623 | |
---|
624 | ! |
---|
625 | ! We get the two standard varibales |
---|
626 | ! Surface temperature |
---|
627 | ! |
---|
628 | ALLOCATE(tsol(iml,jml)) |
---|
629 | ALLOCATE(tmp_var(iml-1,jml)) |
---|
630 | ! |
---|
631 | ! |
---|
632 | |
---|
633 | CALL flinget(fid_phys, 'ST', iml_phys, jml_phys, |
---|
634 | .llm_tmp, ttm_tmp, 1, 1, var_ana) |
---|
635 | |
---|
636 | title='ST' |
---|
637 | CALL conf_dat2d(title,iml_phys, jml_phys, lon_ini, lat_ini, |
---|
638 | . lon_rad, lat_rad, var_ana , interbar ) |
---|
639 | |
---|
640 | IF ( interbar ) THEN |
---|
641 | WRITE(6,*) '-------------------------------------------------', |
---|
642 | ,'--------------' |
---|
643 | WRITE(6,*) '$$$ Utilisation de l interpolation barycentrique ', |
---|
644 | , ' pour ST $$$ ' |
---|
645 | WRITE(6,*) '-------------------------------------------------', |
---|
646 | ,'--------------' |
---|
647 | CALL inter_barxy ( iml_phys,jml_phys -1,lon_rad,lat_rad , |
---|
648 | , var_ana, iml-1, jml-1, lon_in2, lat_in2, jml, tmp_var ) |
---|
649 | ELSE |
---|
650 | CALL grille_m(iml_phys, jml_phys, lon_rad, lat_rad, |
---|
651 | . var_ana, iml-1, jml, lon_in, lat_in, tmp_var ) |
---|
652 | ENDIF |
---|
653 | |
---|
654 | CALL gr_int_dyn(tmp_var, tsol, iml-1, jml) |
---|
655 | ! |
---|
656 | ! Soil moisture |
---|
657 | ! |
---|
658 | ALLOCATE(qsol(iml,jml)) |
---|
659 | CALL flinget(fid_phys, 'CDSW', iml_phys, jml_phys, |
---|
660 | . llm_tmp, ttm_tmp, 1, 1, var_ana) |
---|
661 | |
---|
662 | title='CDSW' |
---|
663 | CALL conf_dat2d(title,iml_phys, jml_phys, lon_ini, lat_ini, |
---|
664 | . lon_rad, lat_rad, var_ana, interbar ) |
---|
665 | |
---|
666 | IF ( interbar ) THEN |
---|
667 | WRITE(6,*) '-------------------------------------------------', |
---|
668 | ,'--------------' |
---|
669 | WRITE(6,*) '$$$ Utilisation de l interpolation barycentrique ', |
---|
670 | , ' pour CDSW $$$ ' |
---|
671 | WRITE(6,*) '-------------------------------------------------', |
---|
672 | ,'--------------' |
---|
673 | CALL inter_barxy ( iml_phys,jml_phys -1,lon_rad,lat_rad , |
---|
674 | , var_ana, iml-1, jml-1, lon_in2, lat_in2, jml, tmp_var ) |
---|
675 | ELSE |
---|
676 | CALL grille_m(iml_phys, jml_phys, lon_rad, lat_rad, |
---|
677 | . var_ana, iml-1, jml, lon_in, lat_in, tmp_var ) |
---|
678 | ENDIF |
---|
679 | c |
---|
680 | CALL gr_int_dyn(tmp_var, qsol, iml-1, jml) |
---|
681 | ! |
---|
682 | CALL flinclo(fid_phys) |
---|
683 | ! |
---|
684 | END SUBROUTINE start_init_phys |
---|
685 | ! |
---|
686 | ! |
---|
687 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
688 | ! |
---|
689 | ! |
---|
690 | SUBROUTINE startget_dyn(varname, iml, jml, lon_in, lat_in, |
---|
691 | . lml, pls, workvar, champ, val_exp,jml2, lon_in2, lat_in2 , |
---|
692 | , interbar ) |
---|
693 | ! |
---|
694 | ! ARGUMENTS |
---|
695 | ! |
---|
696 | CHARACTER*(*), INTENT(in) :: varname |
---|
697 | INTEGER, INTENT(in) :: iml, jml, lml, jml2 |
---|
698 | REAL, INTENT(in) :: lon_in(iml), lat_in(jml) |
---|
699 | REAL, INTENT(in) :: lon_in2(iml), lat_in2(jml2) |
---|
700 | REAL, INTENT(in) :: pls(iml, jml, lml) |
---|
701 | REAL, INTENT(in) :: workvar(iml, jml, lml) |
---|
702 | REAL, INTENT(inout) :: champ(iml, jml, lml) |
---|
703 | REAL, INTENT(in) :: val_exp |
---|
704 | LOGICAL interbar |
---|
705 | ! |
---|
706 | ! LOCAL |
---|
707 | ! |
---|
708 | INTEGER :: il, ij, ii |
---|
709 | REAL :: xppn, xpps |
---|
710 | ! |
---|
711 | #include "dimensions.h" |
---|
712 | #include "paramet.h" |
---|
713 | #include "comgeom2.h" |
---|
714 | #include "comconst.h" |
---|
715 | ! |
---|
716 | ! This routine only works if the variable does not exist or is constant |
---|
717 | ! |
---|
718 | IF ( MINVAL(champ(:,:,:)).EQ.MAXVAL(champ(:,:,:)) .AND. |
---|
719 | . MINVAL(champ(:,:,:)).EQ.val_exp ) THEN |
---|
720 | ! |
---|
721 | SELECTCASE(varname) |
---|
722 | CASE ('u') |
---|
723 | IF ( .NOT.ALLOCATED(psol_dyn)) THEN |
---|
724 | CALL start_init_dyn( iml, jml, lon_in, lat_in, jml2 , |
---|
725 | . lon_in2,lat_in2 , interbar ) |
---|
726 | ENDIF |
---|
727 | CALL start_inter_3d('U', iml, jml, lml, lon_in, |
---|
728 | . lat_in, jml2, lon_in2, lat_in2, pls, champ,interbar ) |
---|
729 | DO il=1,lml |
---|
730 | DO ij=1,jml |
---|
731 | DO ii=1,iml-1 |
---|
732 | champ(ii,ij,il) = champ(ii,ij,il) * cu(ii,ij) |
---|
733 | ENDDO |
---|
734 | champ(iml,ij, il) = champ(1,ij, il) |
---|
735 | ENDDO |
---|
736 | ENDDO |
---|
737 | CASE ('v') |
---|
738 | IF ( .NOT.ALLOCATED(psol_dyn)) THEN |
---|
739 | CALL start_init_dyn( iml, jml, lon_in, lat_in , jml2, |
---|
740 | . lon_in2, lat_in2 , interbar ) |
---|
741 | ENDIF |
---|
742 | CALL start_inter_3d('V', iml, jml, lml, lon_in, |
---|
743 | . lat_in, jml2, lon_in2, lat_in2, pls, champ, interbar ) |
---|
744 | DO il=1,lml |
---|
745 | DO ij=1,jml |
---|
746 | DO ii=1,iml-1 |
---|
747 | champ(ii,ij,il) = champ(ii,ij,il) * cv(ii,ij) |
---|
748 | ENDDO |
---|
749 | champ(iml,ij, il) = champ(1,ij, il) |
---|
750 | ENDDO |
---|
751 | ENDDO |
---|
752 | CASE ('t') |
---|
753 | IF ( .NOT.ALLOCATED(psol_dyn)) THEN |
---|
754 | CALL start_init_dyn( iml, jml, lon_in, lat_in, jml2 , |
---|
755 | . lon_in2, lat_in2 ,interbar ) |
---|
756 | ENDIF |
---|
757 | CALL start_inter_3d('TEMP', iml, jml, lml, lon_in, |
---|
758 | . lat_in, jml2, lon_in2, lat_in2, pls, champ, interbar ) |
---|
759 | |
---|
760 | CASE ('tpot') |
---|
761 | IF ( .NOT.ALLOCATED(psol_dyn)) THEN |
---|
762 | CALL start_init_dyn( iml, jml, lon_in, lat_in , jml2 , |
---|
763 | . lon_in2, lat_in2 , interbar ) |
---|
764 | ENDIF |
---|
765 | CALL start_inter_3d('TEMP', iml, jml, lml, lon_in, |
---|
766 | . lat_in, jml2, lon_in2, lat_in2, pls, champ, interbar ) |
---|
767 | IF ( MINVAL(workvar(:,:,:)) .NE. MAXVAL(workvar(:,:,:)) ) |
---|
768 | . THEN |
---|
769 | DO il=1,lml |
---|
770 | DO ij=1,jml |
---|
771 | DO ii=1,iml-1 |
---|
772 | champ(ii,ij,il) = champ(ii,ij,il) * cpp |
---|
773 | . / workvar(ii,ij,il) |
---|
774 | ENDDO |
---|
775 | champ(iml,ij,il) = champ(1,ij,il) |
---|
776 | ENDDO |
---|
777 | ENDDO |
---|
778 | DO il=1,lml |
---|
779 | xppn = SUM(aire(:,1)*champ(:,1,il))/apoln |
---|
780 | xpps = SUM(aire(:,jml)*champ(:,jml,il))/apols |
---|
781 | champ(:,1,il) = xppn |
---|
782 | champ(:,jml,il) = xpps |
---|
783 | ENDDO |
---|
784 | ELSE |
---|
785 | WRITE(*,*)'Could not compute potential temperature as the' |
---|
786 | WRITE(*,*)'Exner function is missing or constant.' |
---|
787 | STOP |
---|
788 | ENDIF |
---|
789 | CASE ('q') |
---|
790 | IF ( .NOT.ALLOCATED(psol_dyn)) THEN |
---|
791 | CALL start_init_dyn( iml, jml, lon_in, lat_in, jml2 , |
---|
792 | . lon_in2, lat_in2 , interbar ) |
---|
793 | ENDIF |
---|
794 | CALL start_inter_3d('R', iml, jml, lml, lon_in, lat_in, |
---|
795 | . jml2, lon_in2, lat_in2, pls, champ, interbar ) |
---|
796 | IF ( MINVAL(workvar(:,:,:)) .NE. MAXVAL(workvar(:,:,:)) ) |
---|
797 | . THEN |
---|
798 | DO il=1,lml |
---|
799 | DO ij=1,jml |
---|
800 | DO ii=1,iml-1 |
---|
801 | champ(ii,ij,il) = 0.01 * champ(ii,ij,il) * |
---|
802 | . workvar(ii,ij,il) |
---|
803 | ENDDO |
---|
804 | champ(iml,ij,il) = champ(1,ij,il) |
---|
805 | ENDDO |
---|
806 | ENDDO |
---|
807 | WHERE ( champ .LT. 0.) champ = 1.0E-10 |
---|
808 | DO il=1,lml |
---|
809 | xppn = SUM(aire(:,1)*champ(:,1,il))/apoln |
---|
810 | xpps = SUM(aire(:,jml)*champ(:,jml,il))/apols |
---|
811 | champ(:,1,il) = xppn |
---|
812 | champ(:,jml,il) = xpps |
---|
813 | ENDDO |
---|
814 | ELSE |
---|
815 | WRITE(*,*)'Could not compute specific humidity as the' |
---|
816 | WRITE(*,*)'saturated humidity is missing or constant.' |
---|
817 | STOP |
---|
818 | ENDIF |
---|
819 | CASE DEFAULT |
---|
820 | WRITE(*,*) 'startget_dyn' |
---|
821 | WRITE(*,*) 'No rule is present to extract variable ', |
---|
822 | . varname(:LEN_TRIM(varname)),' from any data set' |
---|
823 | STOP |
---|
824 | END SELECT |
---|
825 | ENDIF |
---|
826 | END SUBROUTINE startget_dyn |
---|
827 | ! |
---|
828 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
829 | ! |
---|
830 | SUBROUTINE start_init_dyn( iml, jml, lon_in, lat_in,jml2,lon_in2 , |
---|
831 | , lat_in2 , interbar ) |
---|
832 | ! |
---|
833 | INTEGER, INTENT(in) :: iml, jml, jml2 |
---|
834 | REAL, INTENT(in) :: lon_in(iml), lat_in(jml) |
---|
835 | REAL, INTENT(in) :: lon_in2(iml), lat_in2(jml2) |
---|
836 | LOGICAL interbar |
---|
837 | ! |
---|
838 | ! LOCAL |
---|
839 | ! |
---|
840 | REAL :: lev(1), date, dt |
---|
841 | INTEGER :: itau(1) |
---|
842 | INTEGER :: i, j |
---|
843 | integer :: iret |
---|
844 | ! |
---|
845 | CHARACTER*120 :: physfname |
---|
846 | LOGICAL :: check=.TRUE. |
---|
847 | ! |
---|
848 | REAL, ALLOCATABLE :: lon_rad(:), lat_rad(:) |
---|
849 | REAL, ALLOCATABLE :: lon_ini(:), lat_ini(:) |
---|
850 | REAL, ALLOCATABLE :: var_ana(:,:), tmp_var(:,:), z(:,:) |
---|
851 | REAL, ALLOCATABLE :: xppn(:), xpps(:) |
---|
852 | LOGICAL :: allo |
---|
853 | ! |
---|
854 | ! |
---|
855 | #include "dimensions.h" |
---|
856 | #include "paramet.h" |
---|
857 | #include "comgeom2.h" |
---|
858 | |
---|
859 | CHARACTER*25 title |
---|
860 | |
---|
861 | ! |
---|
862 | physfname = 'ECDYN.nc' |
---|
863 | ! |
---|
864 | IF ( check ) WRITE(*,*) 'Opening the surface analysis' |
---|
865 | ! |
---|
866 | CALL flininfo(physfname, iml_dyn, jml_dyn, llm_dyn, |
---|
867 | . ttm_dyn, fid_dyn) |
---|
868 | IF ( check ) WRITE(*,*) 'Values read: ', iml_dyn, jml_dyn, |
---|
869 | . llm_dyn, ttm_dyn |
---|
870 | ! |
---|
871 | ALLOCATE (lat_dyn(iml_dyn,jml_dyn), stat=iret) |
---|
872 | ALLOCATE (lon_dyn(iml_dyn,jml_dyn), stat=iret) |
---|
873 | ALLOCATE (levdyn_ini(llm_dyn), stat=iret) |
---|
874 | ! |
---|
875 | CALL flinopen(physfname, .FALSE., iml_dyn, jml_dyn, llm_dyn, |
---|
876 | . lon_dyn, lat_dyn, levdyn_ini, ttm_dyn, |
---|
877 | . itau, date, dt, fid_dyn) |
---|
878 | ! |
---|
879 | |
---|
880 | allo = allocated (var_ana) |
---|
881 | if (allo) then |
---|
882 | DEALLOCATE(var_ana, stat=iret) |
---|
883 | endif |
---|
884 | ALLOCATE(var_ana(iml_dyn, jml_dyn), stat=iret) |
---|
885 | |
---|
886 | allo = allocated (lon_rad) |
---|
887 | if (allo) then |
---|
888 | DEALLOCATE(lon_rad, stat=iret) |
---|
889 | endif |
---|
890 | |
---|
891 | ALLOCATE(lon_rad(iml_dyn), stat=iret) |
---|
892 | ALLOCATE(lon_ini(iml_dyn)) |
---|
893 | |
---|
894 | IF ( MAXVAL(lon_dyn(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN |
---|
895 | lon_ini(:) = lon_dyn(:,1) * 2.0 * ASIN(1.0) / 180.0 |
---|
896 | ELSE |
---|
897 | lon_ini(:) = lon_dyn(:,1) |
---|
898 | ENDIF |
---|
899 | |
---|
900 | ALLOCATE(lat_rad(jml_dyn)) |
---|
901 | ALLOCATE(lat_ini(jml_dyn)) |
---|
902 | |
---|
903 | IF ( MAXVAL(lat_dyn(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN |
---|
904 | lat_ini(:) = lat_dyn(1,:) * 2.0 * ASIN(1.0) / 180.0 |
---|
905 | ELSE |
---|
906 | lat_ini(:) = lat_dyn(1,:) |
---|
907 | ENDIF |
---|
908 | ! |
---|
909 | |
---|
910 | |
---|
911 | ALLOCATE(z(iml, jml)) |
---|
912 | ALLOCATE(tmp_var(iml-1,jml)) |
---|
913 | ! |
---|
914 | CALL flinget(fid_dyn, 'Z', iml_dyn, jml_dyn, 0, ttm_dyn, |
---|
915 | . 1, 1, var_ana) |
---|
916 | c |
---|
917 | title='Z' |
---|
918 | CALL conf_dat2d( title,iml_dyn, jml_dyn,lon_ini, lat_ini, |
---|
919 | . lon_rad, lat_rad, var_ana, interbar ) |
---|
920 | c |
---|
921 | IF ( interbar ) THEN |
---|
922 | WRITE(6,*) '-------------------------------------------------', |
---|
923 | ,'--------------' |
---|
924 | WRITE(6,*) '$$$ Utilisation de l interpolation barycentrique ', |
---|
925 | , ' pour Z $$$ ' |
---|
926 | WRITE(6,*) '-------------------------------------------------', |
---|
927 | ,'--------------' |
---|
928 | CALL inter_barxy ( iml_dyn,jml_dyn -1,lon_rad,lat_rad , |
---|
929 | , var_ana, iml-1, jml-1, lon_in2, lat_in2, jml, tmp_var) |
---|
930 | ELSE |
---|
931 | CALL grille_m(iml_dyn, jml_dyn , lon_rad, lat_rad, var_ana, |
---|
932 | . iml-1, jml, lon_in, lat_in, tmp_var) |
---|
933 | ENDIF |
---|
934 | |
---|
935 | CALL gr_int_dyn(tmp_var, z, iml-1, jml) |
---|
936 | ! |
---|
937 | ALLOCATE(psol_dyn(iml, jml)) |
---|
938 | ! |
---|
939 | CALL flinget(fid_dyn, 'SP', iml_dyn, jml_dyn, 0, ttm_dyn, |
---|
940 | . 1, 1, var_ana) |
---|
941 | |
---|
942 | title='SP' |
---|
943 | CALL conf_dat2d( title,iml_dyn, jml_dyn,lon_ini, lat_ini, |
---|
944 | . lon_rad, lat_rad, var_ana, interbar ) |
---|
945 | |
---|
946 | IF ( interbar ) THEN |
---|
947 | WRITE(6,*) '-------------------------------------------------', |
---|
948 | ,'--------------' |
---|
949 | WRITE(6,*) '$$$ Utilisation de l interpolation barycentrique ', |
---|
950 | , ' pour SP $$$ ' |
---|
951 | WRITE(6,*) '-------------------------------------------------', |
---|
952 | ,'--------------' |
---|
953 | CALL inter_barxy ( iml_dyn,jml_dyn -1,lon_rad,lat_rad , |
---|
954 | , var_ana, iml-1, jml-1, lon_in2, lat_in2, jml, tmp_var) |
---|
955 | ELSE |
---|
956 | CALL grille_m(iml_dyn, jml_dyn , lon_rad, lat_rad, var_ana, |
---|
957 | . iml-1, jml, lon_in, lat_in, tmp_var ) |
---|
958 | ENDIF |
---|
959 | |
---|
960 | CALL gr_int_dyn(tmp_var, psol_dyn, iml-1, jml) |
---|
961 | ! |
---|
962 | IF ( .NOT.ALLOCATED(tsol)) THEN |
---|
963 | ! These variables may have been allocated by the need to |
---|
964 | ! create a start field for them or by the varibale |
---|
965 | ! coming out of the restart file. In case we dor have it we will initialize it. |
---|
966 | ! |
---|
967 | CALL start_init_phys( iml, jml, lon_in, lat_in,jml2,lon_in2, |
---|
968 | . lat_in2 , interbar ) |
---|
969 | ELSE |
---|
970 | IF ( SIZE(tsol) .NE. SIZE(psol_dyn) ) THEN |
---|
971 | WRITE(*,*) 'start_init_dyn :' |
---|
972 | WRITE(*,*) 'The temperature field we have does not ', |
---|
973 | . 'have the right size' |
---|
974 | STOP |
---|
975 | ENDIF |
---|
976 | ENDIF |
---|
977 | IF ( .NOT.ALLOCATED(phis)) THEN |
---|
978 | ! |
---|
979 | ! These variables may have been allocated by the need to create a start field for them or by the varibale |
---|
980 | ! coming out of the restart file. In case we dor have it we will initialize it. |
---|
981 | ! |
---|
982 | CALL start_init_orog( iml, jml, lon_in, lat_in, jml2, lon_in2 , |
---|
983 | . lat_in2 , interbar ) |
---|
984 | ! |
---|
985 | ELSE |
---|
986 | ! |
---|
987 | IF (SIZE(phis) .NE. SIZE(psol_dyn)) THEN |
---|
988 | ! |
---|
989 | WRITE(*,*) 'start_init_dyn :' |
---|
990 | WRITE(*,*) 'The orography field we have does not ', |
---|
991 | . ' have the right size' |
---|
992 | STOP |
---|
993 | ENDIF |
---|
994 | ! |
---|
995 | ENDIF |
---|
996 | ! |
---|
997 | ! PSOL is computed in Pascals |
---|
998 | ! |
---|
999 | ! |
---|
1000 | DO j = 1, jml |
---|
1001 | DO i = 1, iml-1 |
---|
1002 | psol_dyn(i,j) = psol_dyn(i,j)*(1.0+(z(i,j)-phis(i,j)) |
---|
1003 | . /287.0/tsol(i,j)) |
---|
1004 | ENDDO |
---|
1005 | psol_dyn(iml,j) = psol_dyn(1,j) |
---|
1006 | ENDDO |
---|
1007 | ! |
---|
1008 | ! |
---|
1009 | ALLOCATE(xppn(iml-1)) |
---|
1010 | ALLOCATE(xpps(iml-1)) |
---|
1011 | ! |
---|
1012 | DO i = 1, iml-1 |
---|
1013 | xppn(i) = aire( i,1) * psol_dyn( i,1) |
---|
1014 | xpps(i) = aire( i,jml) * psol_dyn( i,jml) |
---|
1015 | ENDDO |
---|
1016 | ! |
---|
1017 | DO i = 1, iml |
---|
1018 | psol_dyn(i,1 ) = SUM(xppn)/apoln |
---|
1019 | psol_dyn(i,jml) = SUM(xpps)/apols |
---|
1020 | ENDDO |
---|
1021 | ! |
---|
1022 | RETURN |
---|
1023 | ! |
---|
1024 | END SUBROUTINE start_init_dyn |
---|
1025 | ! |
---|
1026 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
1027 | ! |
---|
1028 | SUBROUTINE start_inter_3d(varname, iml, jml, lml, lon_in, |
---|
1029 | . lat_in, jml2, lon_in2, lat_in2, pls_in, var3d, interbar ) |
---|
1030 | ! |
---|
1031 | ! This subroutine gets a variables from a 3D file and does the interpolations needed |
---|
1032 | ! |
---|
1033 | ! |
---|
1034 | ! ARGUMENTS |
---|
1035 | ! |
---|
1036 | CHARACTER*(*) :: varname |
---|
1037 | INTEGER :: iml, jml, lml, jml2 |
---|
1038 | REAL :: lon_in(iml), lat_in(jml), pls_in(iml, jml, lml) |
---|
1039 | REAL :: lon_in2(iml) , lat_in2(jml2) |
---|
1040 | REAL :: var3d(iml, jml, lml) |
---|
1041 | LOGICAL interbar |
---|
1042 | real chmin,chmax |
---|
1043 | ! |
---|
1044 | ! LOCAL |
---|
1045 | ! |
---|
1046 | CHARACTER*25 title |
---|
1047 | INTEGER :: ii, ij, il, jsort,i,j,l |
---|
1048 | REAL :: bx, by |
---|
1049 | REAL, ALLOCATABLE :: lon_rad(:), lat_rad(:) |
---|
1050 | REAL, ALLOCATABLE :: lon_ini(:), lat_ini(:) , lev_dyn(:) |
---|
1051 | REAL, ALLOCATABLE :: var_tmp2d(:,:), var_tmp3d(:,:,:) |
---|
1052 | REAL, ALLOCATABLE :: ax(:), ay(:), yder(:) |
---|
1053 | REAL, ALLOCATABLE :: varrr(:,:,:) |
---|
1054 | INTEGER, ALLOCATABLE :: lind(:) |
---|
1055 | ! |
---|
1056 | LOGICAL :: check = .TRUE. |
---|
1057 | ! |
---|
1058 | IF ( .NOT. ALLOCATED(var_ana3d)) THEN |
---|
1059 | ALLOCATE(var_ana3d(iml_dyn, jml_dyn, llm_dyn)) |
---|
1060 | ENDIF |
---|
1061 | ALLOCATE(varrr(iml_dyn, jml_dyn, llm_dyn)) |
---|
1062 | ! |
---|
1063 | ! |
---|
1064 | IF ( check) WRITE(*,*) 'Going into flinget to extract the 3D ', |
---|
1065 | . ' field.', fid_dyn |
---|
1066 | IF ( check) WRITE(*,*) fid_dyn, varname, iml_dyn, jml_dyn, |
---|
1067 | . llm_dyn,ttm_dyn |
---|
1068 | ! |
---|
1069 | CALL flinget(fid_dyn, varname, iml_dyn, jml_dyn, llm_dyn, |
---|
1070 | . ttm_dyn, 1, 1, var_ana3d) |
---|
1071 | ! |
---|
1072 | IF ( check) WRITE(*,*) 'Allocating space for the interpolation', |
---|
1073 | . iml, jml, llm_dyn |
---|
1074 | ! |
---|
1075 | ALLOCATE(lon_rad(iml_dyn)) |
---|
1076 | ALLOCATE(lon_ini(iml_dyn)) |
---|
1077 | |
---|
1078 | IF ( MAXVAL(lon_dyn(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN |
---|
1079 | lon_ini(:) = lon_dyn(:,1) * 2.0 * ASIN(1.0) / 180.0 |
---|
1080 | ELSE |
---|
1081 | lon_ini(:) = lon_dyn(:,1) |
---|
1082 | ENDIF |
---|
1083 | |
---|
1084 | ALLOCATE(lat_rad(jml_dyn)) |
---|
1085 | ALLOCATE(lat_ini(jml_dyn)) |
---|
1086 | |
---|
1087 | ALLOCATE(lev_dyn(llm_dyn)) |
---|
1088 | |
---|
1089 | IF ( MAXVAL(lat_dyn(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN |
---|
1090 | lat_ini(:) = lat_dyn(1,:) * 2.0 * ASIN(1.0) / 180.0 |
---|
1091 | ELSE |
---|
1092 | lat_ini(:) = lat_dyn(1,:) |
---|
1093 | ENDIF |
---|
1094 | ! |
---|
1095 | |
---|
1096 | CALL conf_dat3d ( varname,iml_dyn, jml_dyn, llm_dyn, lon_ini, |
---|
1097 | . lat_ini, levdyn_ini, lon_rad, lat_rad, lev_dyn, var_ana3d , |
---|
1098 | , interbar ) |
---|
1099 | |
---|
1100 | ALLOCATE(var_tmp2d(iml-1, jml)) |
---|
1101 | ALLOCATE(var_tmp3d(iml, jml, llm_dyn)) |
---|
1102 | ALLOCATE(ax(llm_dyn)) |
---|
1103 | ALLOCATE(ay(llm_dyn)) |
---|
1104 | ALLOCATE(yder(llm_dyn)) |
---|
1105 | ALLOCATE(lind(llm_dyn)) |
---|
1106 | ! |
---|
1107 | |
---|
1108 | DO il=1,llm_dyn |
---|
1109 | ! |
---|
1110 | IF( interbar ) THEN |
---|
1111 | IF( il.EQ.1 ) THEN |
---|
1112 | WRITE(6,*) '-------------------------------------------------', |
---|
1113 | ,'--------------' |
---|
1114 | WRITE(6,*) '$$$ Utilisation de l interpolation barycentrique ', |
---|
1115 | , ' pour ', varname |
---|
1116 | WRITE(6,*) '-------------------------------------------------', |
---|
1117 | ,'--------------' |
---|
1118 | ENDIF |
---|
1119 | CALL inter_barxy ( iml_dyn, jml_dyn -1,lon_rad, lat_rad, |
---|
1120 | , var_ana3d(:,:,il),iml-1, jml2, lon_in2, lat_in2,jml,var_tmp2d ) |
---|
1121 | ELSE |
---|
1122 | CALL grille_m(iml_dyn, jml_dyn, lon_rad, lat_rad, |
---|
1123 | . var_ana3d(:,:,il), iml-1, jml, lon_in, lat_in, var_tmp2d ) |
---|
1124 | ENDIF |
---|
1125 | ! |
---|
1126 | CALL gr_int_dyn(var_tmp2d, var_tmp3d(:,:,il), iml-1, jml) |
---|
1127 | ! |
---|
1128 | ENDDO |
---|
1129 | ! |
---|
1130 | DO il=1,llm_dyn |
---|
1131 | lind(il) = llm_dyn-il+1 |
---|
1132 | ENDDO |
---|
1133 | ! |
---|
1134 | c |
---|
1135 | c ... Pour l'interpolation verticale ,on interpole du haut de l'atmosphere |
---|
1136 | c vers le sol ... |
---|
1137 | c |
---|
1138 | DO ij=1,jml |
---|
1139 | DO ii=1,iml-1 |
---|
1140 | ! |
---|
1141 | ax(:) = lev_dyn(lind(:)) |
---|
1142 | ay(:) = var_tmp3d(ii, ij, lind(:)) |
---|
1143 | ! |
---|
1144 | |
---|
1145 | CALL SPLINE(ax, ay, llm_dyn, 1.e30, 1.e30, yder) |
---|
1146 | ! |
---|
1147 | DO il=1,lml |
---|
1148 | bx = pls_in(ii, ij, il) |
---|
1149 | CALL SPLINT(ax, ay, yder, llm_dyn, bx, by) |
---|
1150 | var3d(ii, ij, il) = by |
---|
1151 | ENDDO |
---|
1152 | ! |
---|
1153 | ENDDO |
---|
1154 | var3d(iml, ij, :) = var3d(1, ij, :) |
---|
1155 | ENDDO |
---|
1156 | |
---|
1157 | do il=1,lml |
---|
1158 | call minmax(iml*jml,var3d(1,1,il),chmin,chmax) |
---|
1159 | SELECTCASE(varname) |
---|
1160 | CASE('U') |
---|
1161 | WRITE(*,*) ' U min max l ',il,chmin,chmax |
---|
1162 | CASE('V') |
---|
1163 | WRITE(*,*) ' V min max l ',il,chmin,chmax |
---|
1164 | CASE('TEMP') |
---|
1165 | WRITE(*,*) ' TEMP min max l ',il,chmin,chmax |
---|
1166 | CASE('R') |
---|
1167 | WRITE(*,*) ' R min max l ',il,chmin,chmax |
---|
1168 | END SELECT |
---|
1169 | enddo |
---|
1170 | |
---|
1171 | DEALLOCATE(lon_rad) |
---|
1172 | DEALLOCATE(lat_rad) |
---|
1173 | DEALLOCATE(var_tmp2d) |
---|
1174 | DEALLOCATE(var_tmp3d) |
---|
1175 | DEALLOCATE(ax) |
---|
1176 | DEALLOCATE(ay) |
---|
1177 | DEALLOCATE(yder) |
---|
1178 | DEALLOCATE(lind) |
---|
1179 | |
---|
1180 | ! |
---|
1181 | RETURN |
---|
1182 | ! |
---|
1183 | END SUBROUTINE start_inter_3d |
---|
1184 | ! |
---|
1185 | END MODULE startvar |
---|