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