1 |
module dynetat0_m |
2 |
|
3 |
use dimensions, only: iim, jjm |
4 |
|
5 |
IMPLICIT NONE |
6 |
|
7 |
private iim, jjm, principal_cshift, invert_zoom_x, funcd |
8 |
|
9 |
INTEGER day_ini |
10 |
! day number at the beginning of the run, based at value 1 on |
11 |
! January 1st of annee_ref |
12 |
|
13 |
integer:: day_ref = 1 ! jour de l'ann\'ee de l'\'etat initial |
14 |
! (= 350 si 20 d\'ecembre par exemple) |
15 |
|
16 |
integer:: annee_ref = 1998 ! Annee de l'etat initial (avec 4 chiffres) |
17 |
|
18 |
REAL, protected:: clon ! longitude of the center of the zoom, in rad |
19 |
real, protected:: clat ! latitude of the center of the zoom, in rad |
20 |
|
21 |
real, protected:: grossismx, grossismy |
22 |
! facteurs de grossissement du zoom, selon la longitude et la latitude |
23 |
! = 2 si 2 fois, = 3 si 3 fois, etc. |
24 |
|
25 |
real, protected:: dzoomx, dzoomy |
26 |
! extensions en longitude et latitude de la zone du zoom (fractions |
27 |
! de la zone totale) |
28 |
|
29 |
real, protected:: taux, tauy |
30 |
! raideur de la transition de l'int\'erieur \`a l'ext\'erieur du zoom |
31 |
|
32 |
real rlatu(jjm + 1) |
33 |
! latitudes of points of the "scalar" and "u" grid, in rad |
34 |
|
35 |
real rlatv(jjm) |
36 |
! latitudes of points of the "v" grid, in rad, in decreasing order |
37 |
|
38 |
real rlonu(iim + 1) ! longitudes of points of the "u" grid, in rad |
39 |
|
40 |
real rlonv(iim + 1) |
41 |
! longitudes of points of the "scalar" and "v" grid, in rad |
42 |
|
43 |
real, protected:: xprimu(iim + 1), xprimv(iim + 1) |
44 |
! 2 pi / iim * (derivative of the longitudinal zoom function)(rlon[uv]) |
45 |
|
46 |
REAL, protected:: xprimm025(iim + 1), xprimp025(iim + 1) |
47 |
REAL, protected:: rlatu1(jjm), rlatu2(jjm), yprimu1(jjm), yprimu2(jjm) |
48 |
REAL ang0, etot0, ptot0, ztot0, stot0 |
49 |
INTEGER, PARAMETER, private:: nmax = 30000 |
50 |
DOUBLE PRECISION, private:: abs_y |
51 |
|
52 |
save |
53 |
|
54 |
contains |
55 |
|
56 |
SUBROUTINE dynetat0(vcov, ucov, teta, q, masse, ps, phis) |
57 |
|
58 |
! From dynetat0.F, version 1.2, 2004/06/22 11:45:30 |
59 |
! Authors: P. Le Van, L. Fairhead |
60 |
! This procedure reads the initial state of the atmosphere. |
61 |
|
62 |
use comconst, only: dtvr |
63 |
use conf_gcm_m, only: raz_date |
64 |
use dimensions, only: iim, jjm, llm, nqmx |
65 |
use disvert_m, only: pa |
66 |
use iniadvtrac_m, only: tname |
67 |
use netcdf, only: NF90_NOWRITE, NF90_NOERR |
68 |
use netcdf95, only: NF95_GET_VAR, nf95_open, nf95_inq_varid, NF95_CLOSE, & |
69 |
NF95_Gw_VAR |
70 |
use nr_util, only: assert |
71 |
use temps, only: itau_dyn |
72 |
use unit_nml_m, only: unit_nml |
73 |
|
74 |
REAL, intent(out):: vcov(: , :, :) ! (iim + 1, jjm, llm) |
75 |
REAL, intent(out):: ucov(:, :, :) ! (iim + 1, jjm + 1, llm) |
76 |
REAL, intent(out):: teta(:, :, :) ! (iim + 1, jjm + 1, llm) |
77 |
REAL, intent(out):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx) |
78 |
REAL, intent(out):: masse(:, :, :) ! (iim + 1, jjm + 1, llm) |
79 |
REAL, intent(out):: ps(:, :) ! (iim + 1, jjm + 1) in Pa |
80 |
REAL, intent(out):: phis(:, :) ! (iim + 1, jjm + 1) |
81 |
|
82 |
! Local variables: |
83 |
INTEGER iq |
84 |
REAL, allocatable:: tab_cntrl(:) ! tableau des param\`etres du run |
85 |
INTEGER ierr, ncid, varid |
86 |
|
87 |
namelist /dynetat0_nml/ day_ref, annee_ref |
88 |
|
89 |
!----------------------------------------------------------------------- |
90 |
|
91 |
print *, "Call sequence information: dynetat0" |
92 |
|
93 |
call assert((/size(ucov, 1), size(vcov, 1), size(masse, 1), size(ps, 1), & |
94 |
size(phis, 1), size(q, 1), size(teta, 1)/) == iim + 1, "dynetat0 iim") |
95 |
call assert((/size(ucov, 2), size(vcov, 2) + 1, size(masse, 2), & |
96 |
size(ps, 2), size(phis, 2), size(q, 2), size(teta, 2)/) == jjm + 1, & |
97 |
"dynetat0 jjm") |
98 |
call assert((/size(vcov, 3), size(ucov, 3), size(teta, 3), size(q, 3), & |
99 |
size(masse, 3)/) == llm, "dynetat0 llm") |
100 |
call assert(size(q, 4) == nqmx, "dynetat0 q nqmx") |
101 |
|
102 |
! Fichier \'etat initial : |
103 |
call nf95_open("start.nc", NF90_NOWRITE, ncid) |
104 |
|
105 |
call nf95_inq_varid(ncid, "controle", varid) |
106 |
call NF95_Gw_VAR(ncid, varid, tab_cntrl) |
107 |
|
108 |
call assert(int(tab_cntrl(1)) == iim, "dynetat0 tab_cntrl iim") |
109 |
call assert(int(tab_cntrl(2)) == jjm, "dynetat0 tab_cntrl jjm") |
110 |
call assert(int(tab_cntrl(3)) == llm, "dynetat0 tab_cntrl llm") |
111 |
|
112 |
IF (dtvr /= tab_cntrl(12)) THEN |
113 |
print *, 'Warning: the time steps from day_step and "start.nc" ' // & |
114 |
'are different.' |
115 |
print *, 'dtvr from day_step: ', dtvr |
116 |
print *, 'dtvr from "start.nc": ', tab_cntrl(12) |
117 |
print *, 'Using the value from day_step.' |
118 |
ENDIF |
119 |
|
120 |
etot0 = tab_cntrl(13) |
121 |
ptot0 = tab_cntrl(14) |
122 |
ztot0 = tab_cntrl(15) |
123 |
stot0 = tab_cntrl(16) |
124 |
ang0 = tab_cntrl(17) |
125 |
pa = tab_cntrl(18) |
126 |
|
127 |
clon = tab_cntrl(20) |
128 |
clat = tab_cntrl(21) |
129 |
grossismx = tab_cntrl(22) |
130 |
grossismy = tab_cntrl(23) |
131 |
dzoomx = tab_cntrl(25) |
132 |
dzoomy = tab_cntrl(26) |
133 |
taux = tab_cntrl(28) |
134 |
tauy = tab_cntrl(29) |
135 |
|
136 |
print *, "Enter namelist 'dynetat0_nml'." |
137 |
read(unit=*, nml=dynetat0_nml) |
138 |
write(unit_nml, nml=dynetat0_nml) |
139 |
|
140 |
if (raz_date) then |
141 |
print *, 'Resetting the date, using the namelist.' |
142 |
day_ini = day_ref |
143 |
itau_dyn = 0 |
144 |
else |
145 |
day_ref = tab_cntrl(4) |
146 |
annee_ref = tab_cntrl(5) |
147 |
itau_dyn = tab_cntrl(31) |
148 |
day_ini = tab_cntrl(30) |
149 |
end if |
150 |
|
151 |
print *, "day_ini = ", day_ini |
152 |
|
153 |
call NF95_INQ_VARID (ncid, "rlonu", varid) |
154 |
call NF95_GET_VAR(ncid, varid, rlonu) |
155 |
|
156 |
call NF95_INQ_VARID (ncid, "rlatu", varid) |
157 |
call NF95_GET_VAR(ncid, varid, rlatu) |
158 |
|
159 |
call NF95_INQ_VARID (ncid, "rlonv", varid) |
160 |
call NF95_GET_VAR(ncid, varid, rlonv) |
161 |
|
162 |
call NF95_INQ_VARID (ncid, "rlatv", varid) |
163 |
call NF95_GET_VAR(ncid, varid, rlatv) |
164 |
|
165 |
CALL nf95_inq_varid(ncid, 'xprimu', varid) |
166 |
CALL nf95_get_var(ncid, varid, xprimu) |
167 |
|
168 |
CALL nf95_inq_varid(ncid, 'xprimv', varid) |
169 |
CALL nf95_get_var(ncid, varid, xprimv) |
170 |
|
171 |
CALL nf95_inq_varid(ncid, 'xprimm025', varid) |
172 |
CALL nf95_get_var(ncid, varid, xprimm025) |
173 |
|
174 |
CALL nf95_inq_varid(ncid, 'xprimp025', varid) |
175 |
CALL nf95_get_var(ncid, varid, xprimp025) |
176 |
|
177 |
call NF95_INQ_VARID (ncid, "rlatu1", varid) |
178 |
call NF95_GET_VAR(ncid, varid, rlatu1) |
179 |
|
180 |
call NF95_INQ_VARID (ncid, "rlatu2", varid) |
181 |
call NF95_GET_VAR(ncid, varid, rlatu2) |
182 |
|
183 |
CALL nf95_inq_varid(ncid, 'yprimu1', varid) |
184 |
CALL nf95_get_var(ncid, varid, yprimu1) |
185 |
|
186 |
CALL nf95_inq_varid(ncid, 'yprimu2', varid) |
187 |
CALL nf95_get_var(ncid, varid, yprimu2) |
188 |
|
189 |
call NF95_INQ_VARID (ncid, "phis", varid) |
190 |
call NF95_GET_VAR(ncid, varid, phis) |
191 |
|
192 |
call NF95_INQ_VARID (ncid, "ucov", varid) |
193 |
call NF95_GET_VAR(ncid, varid, ucov) |
194 |
|
195 |
call NF95_INQ_VARID (ncid, "vcov", varid) |
196 |
call NF95_GET_VAR(ncid, varid, vcov) |
197 |
|
198 |
call NF95_INQ_VARID (ncid, "teta", varid) |
199 |
call NF95_GET_VAR(ncid, varid, teta) |
200 |
|
201 |
DO iq = 1, nqmx |
202 |
call NF95_INQ_VARID(ncid, tname(iq), varid, ierr) |
203 |
IF (ierr == NF90_NOERR) THEN |
204 |
call NF95_GET_VAR(ncid, varid, q(:, :, :, iq)) |
205 |
ELSE |
206 |
PRINT *, 'dynetat0: "' // tname(iq) // '" not found, ' // & |
207 |
"setting it to zero..." |
208 |
q(:, :, :, iq) = 0. |
209 |
ENDIF |
210 |
ENDDO |
211 |
|
212 |
call NF95_INQ_VARID (ncid, "masse", varid) |
213 |
call NF95_GET_VAR(ncid, varid, masse) |
214 |
|
215 |
call NF95_INQ_VARID (ncid, "ps", varid) |
216 |
call NF95_GET_VAR(ncid, varid, ps) |
217 |
! Check that there is a single value at each pole: |
218 |
call assert(ps(1, 1) == ps(2:, 1), "dynetat0 ps north pole") |
219 |
call assert(ps(1, jjm + 1) == ps(2:, jjm + 1), "dynetat0 ps south pole") |
220 |
|
221 |
call NF95_CLOSE(ncid) |
222 |
|
223 |
END SUBROUTINE dynetat0 |
224 |
|
225 |
!******************************************************************** |
226 |
|
227 |
subroutine read_serre |
228 |
|
229 |
use unit_nml_m, only: unit_nml |
230 |
use nr_util, only: assert, pi |
231 |
|
232 |
REAL:: clon_deg = 0. ! longitude of the center of the zoom, in degrees |
233 |
real:: clat_deg = 0. ! latitude of the center of the zoom, in degrees |
234 |
|
235 |
namelist /serre_nml/ clon_deg, clat_deg, grossismx, grossismy, dzoomx, & |
236 |
dzoomy, taux, tauy |
237 |
|
238 |
!------------------------------------------------- |
239 |
|
240 |
! Default values: |
241 |
grossismx = 1. |
242 |
grossismy = 1. |
243 |
dzoomx = 0.2 |
244 |
dzoomy = 0.2 |
245 |
taux = 3. |
246 |
tauy = 3. |
247 |
|
248 |
print *, "Enter namelist 'serre_nml'." |
249 |
read(unit=*, nml=serre_nml) |
250 |
write(unit_nml, nml=serre_nml) |
251 |
|
252 |
call assert(grossismx >= 1. .and. grossismy >= 1., "read_serre grossism") |
253 |
call assert(dzoomx > 0., dzoomx < 1., dzoomy < 1., & |
254 |
"read_serre dzoomx dzoomy") |
255 |
clon = clon_deg / 180. * pi |
256 |
clat = clat_deg / 180. * pi |
257 |
|
258 |
end subroutine read_serre |
259 |
|
260 |
!******************************************************************** |
261 |
|
262 |
SUBROUTINE fyhyp |
263 |
|
264 |
! From LMDZ4/libf/dyn3d/fyhyp.F, version 1.2, 2005/06/03 09:11:32 |
265 |
|
266 |
! Author: P. Le Van, from analysis by R. Sadourny |
267 |
|
268 |
! Define rlatu, rlatv, rlatu2, yprimu2, rlatu1, yprimu1, using |
269 |
! clat, grossismy, dzoomy, tauy. |
270 |
|
271 |
! Calcule les latitudes et dérivées dans la grille du GCM pour une |
272 |
! fonction f(y) à dérivée tangente hyperbolique. |
273 |
|
274 |
! Il vaut mieux avoir : grossismy * dzoom < pi / 2 |
275 |
|
276 |
use coefpoly_m, only: coefpoly, a0, a1, a2, a3 |
277 |
USE dimensions, only: jjm |
278 |
use heavyside_m, only: heavyside |
279 |
|
280 |
! Local: |
281 |
|
282 |
INTEGER, PARAMETER:: nmax=30000, nmax2=2*nmax |
283 |
REAL dzoom ! distance totale de la zone du zoom (en radians) |
284 |
DOUBLE PRECISION ylat(jjm + 1), yprim(jjm + 1) |
285 |
DOUBLE PRECISION yuv |
286 |
DOUBLE PRECISION, save:: yt(0:nmax2) |
287 |
DOUBLE PRECISION fhyp(0:nmax2), beta |
288 |
DOUBLE PRECISION, save:: ytprim(0:nmax2) |
289 |
DOUBLE PRECISION fxm(0:nmax2) |
290 |
DOUBLE PRECISION, save:: yf(0:nmax2) |
291 |
DOUBLE PRECISION yypr(0:nmax2) |
292 |
DOUBLE PRECISION yvrai(jjm + 1), yprimm(jjm + 1), ylatt(jjm + 1) |
293 |
DOUBLE PRECISION pi, pis2, epsilon, pisjm |
294 |
DOUBLE PRECISION yo1, yi, ylon2, ymoy, yprimin |
295 |
DOUBLE PRECISION yfi, yf1, ffdy |
296 |
DOUBLE PRECISION ypn |
297 |
DOUBLE PRECISION, save::deply, y00 |
298 |
|
299 |
INTEGER i, j, it, ik, iter, jlat, jjpn |
300 |
INTEGER, save:: jpn |
301 |
DOUBLE PRECISION yi2, heavyy0, heavyy0m |
302 |
DOUBLE PRECISION fa(0:nmax2), fb(0:nmax2) |
303 |
REAL y0min, y0max |
304 |
|
305 |
!------------------------------------------------------------------- |
306 |
|
307 |
print *, "Call sequence information: fyhyp" |
308 |
|
309 |
pi = 2.*asin(1.) |
310 |
pis2 = pi/2. |
311 |
pisjm = pi/real(jjm) |
312 |
epsilon = 1e-3 |
313 |
dzoom = dzoomy*pi |
314 |
|
315 |
DO i = 0, nmax2 |
316 |
yt(i) = -pis2 + real(i)*pi/nmax2 |
317 |
END DO |
318 |
|
319 |
heavyy0m = heavyside(-clat) |
320 |
heavyy0 = heavyside(clat) |
321 |
y0min = 2.*clat*heavyy0m - pis2 |
322 |
y0max = 2.*clat*heavyy0 + pis2 |
323 |
|
324 |
fa = 999.999 |
325 |
fb = 999.999 |
326 |
|
327 |
DO i = 0, nmax2 |
328 |
IF (yt(i)<clat) THEN |
329 |
fa(i) = tauy*(yt(i)-clat + dzoom/2.) |
330 |
fb(i) = (yt(i)-2.*clat*heavyy0m + pis2)*(clat-yt(i)) |
331 |
ELSE IF (yt(i)>clat) THEN |
332 |
fa(i) = tauy*(clat-yt(i) + dzoom/2.) |
333 |
fb(i) = (2.*clat*heavyy0-yt(i) + pis2)*(yt(i)-clat) |
334 |
END IF |
335 |
|
336 |
IF (200.*fb(i)<-fa(i)) THEN |
337 |
fhyp(i) = -1. |
338 |
ELSE IF (200.*fb(i)<fa(i)) THEN |
339 |
fhyp(i) = 1. |
340 |
ELSE |
341 |
fhyp(i) = tanh(fa(i)/fb(i)) |
342 |
END IF |
343 |
|
344 |
IF (yt(i)==clat) fhyp(i) = 1. |
345 |
IF (yt(i)==y0min .OR. yt(i)==y0max) fhyp(i) = -1. |
346 |
END DO |
347 |
|
348 |
! Calcul de beta |
349 |
|
350 |
ffdy = 0. |
351 |
|
352 |
DO i = 1, nmax2 |
353 |
ymoy = 0.5*(yt(i-1) + yt(i)) |
354 |
IF (ymoy<clat) THEN |
355 |
fa(i) = tauy*(ymoy-clat + dzoom/2.) |
356 |
fb(i) = (ymoy-2.*clat*heavyy0m + pis2)*(clat-ymoy) |
357 |
ELSE IF (ymoy>clat) THEN |
358 |
fa(i) = tauy*(clat-ymoy + dzoom/2.) |
359 |
fb(i) = (2.*clat*heavyy0-ymoy + pis2)*(ymoy-clat) |
360 |
END IF |
361 |
|
362 |
IF (200.*fb(i)<-fa(i)) THEN |
363 |
fxm(i) = -1. |
364 |
ELSE IF (200.*fb(i)<fa(i)) THEN |
365 |
fxm(i) = 1. |
366 |
ELSE |
367 |
fxm(i) = tanh(fa(i)/fb(i)) |
368 |
END IF |
369 |
IF (ymoy==clat) fxm(i) = 1. |
370 |
IF (ymoy==y0min .OR. yt(i)==y0max) fxm(i) = -1. |
371 |
ffdy = ffdy + fxm(i)*(yt(i)-yt(i-1)) |
372 |
END DO |
373 |
|
374 |
beta = (grossismy*ffdy-pi)/(ffdy-pi) |
375 |
|
376 |
IF (2. * beta - grossismy <= 0.) THEN |
377 |
print *, 'Attention ! La valeur beta calculee dans la routine fyhyp ' & |
378 |
// 'est mauvaise. Modifier les valeurs de grossismy, tauy ou ' & |
379 |
// 'dzoomy et relancer.' |
380 |
STOP 1 |
381 |
END IF |
382 |
|
383 |
! calcul de Ytprim |
384 |
|
385 |
DO i = 0, nmax2 |
386 |
ytprim(i) = beta + (grossismy-beta)*fhyp(i) |
387 |
END DO |
388 |
|
389 |
! Calcul de Yf |
390 |
|
391 |
yf(0) = -pis2 |
392 |
DO i = 1, nmax2 |
393 |
yypr(i) = beta + (grossismy-beta)*fxm(i) |
394 |
END DO |
395 |
|
396 |
DO i = 1, nmax2 |
397 |
yf(i) = yf(i-1) + yypr(i)*(yt(i)-yt(i-1)) |
398 |
END DO |
399 |
|
400 |
! yuv = 0. si calcul des latitudes aux pts. U |
401 |
! yuv = 0.5 si calcul des latitudes aux pts. V |
402 |
|
403 |
loop_ik: DO ik = 1, 4 |
404 |
IF (ik==1) THEN |
405 |
yuv = 0. |
406 |
jlat = jjm + 1 |
407 |
ELSE IF (ik==2) THEN |
408 |
yuv = 0.5 |
409 |
jlat = jjm |
410 |
ELSE IF (ik==3) THEN |
411 |
yuv = 0.25 |
412 |
jlat = jjm |
413 |
ELSE IF (ik==4) THEN |
414 |
yuv = 0.75 |
415 |
jlat = jjm |
416 |
END IF |
417 |
|
418 |
yo1 = 0. |
419 |
DO j = 1, jlat |
420 |
yo1 = 0. |
421 |
ylon2 = -pis2 + pisjm*(real(j) + yuv-1.) |
422 |
yfi = ylon2 |
423 |
|
424 |
it = nmax2 |
425 |
DO while (it >= 1 .and. yfi < yf(it)) |
426 |
it = it - 1 |
427 |
END DO |
428 |
|
429 |
yi = yt(it) |
430 |
IF (it==nmax2) THEN |
431 |
it = nmax2 - 1 |
432 |
yf(it + 1) = pis2 |
433 |
END IF |
434 |
|
435 |
! Interpolation entre yi(it) et yi(it + 1) pour avoir Y(yi) |
436 |
! et Y'(yi) |
437 |
|
438 |
CALL coefpoly(yf(it), yf(it + 1), ytprim(it), ytprim(it + 1), & |
439 |
yt(it), yt(it + 1)) |
440 |
|
441 |
yf1 = yf(it) |
442 |
yprimin = a1 + 2.*a2*yi + 3.*a3*yi*yi |
443 |
|
444 |
iter = 1 |
445 |
DO |
446 |
yi = yi - (yf1-yfi)/yprimin |
447 |
IF (abs(yi-yo1)<=epsilon .or. iter == 300) exit |
448 |
yo1 = yi |
449 |
yi2 = yi*yi |
450 |
yf1 = a0 + a1*yi + a2*yi2 + a3*yi2*yi |
451 |
yprimin = a1 + 2.*a2*yi + 3.*a3*yi2 |
452 |
END DO |
453 |
if (abs(yi-yo1) > epsilon) then |
454 |
print *, 'Pas de solution.', j, ylon2 |
455 |
STOP 1 |
456 |
end if |
457 |
|
458 |
yprimin = a1 + 2.*a2*yi + 3.*a3*yi*yi |
459 |
yprim(j) = pi/(jjm*yprimin) |
460 |
yvrai(j) = yi |
461 |
END DO |
462 |
|
463 |
DO j = 1, jlat - 1 |
464 |
IF (yvrai(j + 1)<yvrai(j)) THEN |
465 |
print *, 'Problème avec rlat(', j + 1, ') plus petit que rlat(', & |
466 |
j, ')' |
467 |
STOP 1 |
468 |
END IF |
469 |
END DO |
470 |
|
471 |
print *, 'Reorganisation des latitudes pour avoir entre - pi/2 et pi/2' |
472 |
|
473 |
IF (ik==1) THEN |
474 |
ypn = pis2 |
475 |
DO j = jjm + 1, 1, -1 |
476 |
IF (yvrai(j)<=ypn) exit |
477 |
END DO |
478 |
|
479 |
jpn = j |
480 |
y00 = yvrai(jpn) |
481 |
deply = pis2 - y00 |
482 |
END IF |
483 |
|
484 |
DO j = 1, jjm + 1 - jpn |
485 |
ylatt(j) = -pis2 - y00 + yvrai(jpn + j-1) |
486 |
yprimm(j) = yprim(jpn + j-1) |
487 |
END DO |
488 |
|
489 |
jjpn = jpn |
490 |
IF (jlat==jjm) jjpn = jpn - 1 |
491 |
|
492 |
DO j = 1, jjpn |
493 |
ylatt(j + jjm + 1-jpn) = yvrai(j) + deply |
494 |
yprimm(j + jjm + 1-jpn) = yprim(j) |
495 |
END DO |
496 |
|
497 |
! Fin de la reorganisation |
498 |
|
499 |
DO j = 1, jlat |
500 |
ylat(j) = ylatt(jlat + 1-j) |
501 |
yprim(j) = yprimm(jlat + 1-j) |
502 |
END DO |
503 |
|
504 |
DO j = 1, jlat |
505 |
yvrai(j) = ylat(j)*180./pi |
506 |
END DO |
507 |
|
508 |
IF (ik==1) THEN |
509 |
DO j = 1, jjm + 1 |
510 |
rlatu(j) = ylat(j) |
511 |
END DO |
512 |
ELSE IF (ik==2) THEN |
513 |
DO j = 1, jjm |
514 |
rlatv(j) = ylat(j) |
515 |
END DO |
516 |
ELSE IF (ik==3) THEN |
517 |
DO j = 1, jjm |
518 |
rlatu2(j) = ylat(j) |
519 |
yprimu2(j) = yprim(j) |
520 |
END DO |
521 |
ELSE IF (ik==4) THEN |
522 |
DO j = 1, jjm |
523 |
rlatu1(j) = ylat(j) |
524 |
yprimu1(j) = yprim(j) |
525 |
END DO |
526 |
END IF |
527 |
END DO loop_ik |
528 |
|
529 |
DO j = 1, jjm |
530 |
ylat(j) = rlatu(j) - rlatu(j + 1) |
531 |
END DO |
532 |
|
533 |
DO j = 1, jjm |
534 |
IF (rlatu1(j) <= rlatu2(j)) THEN |
535 |
print *, 'Attention ! rlatu1 < rlatu2 ', rlatu1(j), rlatu2(j), j |
536 |
STOP 13 |
537 |
ENDIF |
538 |
|
539 |
IF (rlatu2(j) <= rlatu(j+1)) THEN |
540 |
print *, 'Attention ! rlatu2 < rlatup1 ', rlatu2(j), rlatu(j+1), j |
541 |
STOP 14 |
542 |
ENDIF |
543 |
|
544 |
IF (rlatu(j) <= rlatu1(j)) THEN |
545 |
print *, ' Attention ! rlatu < rlatu1 ', rlatu(j), rlatu1(j), j |
546 |
STOP 15 |
547 |
ENDIF |
548 |
|
549 |
IF (rlatv(j) <= rlatu2(j)) THEN |
550 |
print *, ' Attention ! rlatv < rlatu2 ', rlatv(j), rlatu2(j), j |
551 |
STOP 16 |
552 |
ENDIF |
553 |
|
554 |
IF (rlatv(j) >= rlatu1(j)) THEN |
555 |
print *, ' Attention ! rlatv > rlatu1 ', rlatv(j), rlatu1(j), j |
556 |
STOP 17 |
557 |
ENDIF |
558 |
|
559 |
IF (rlatv(j) >= rlatu(j)) THEN |
560 |
print *, ' Attention ! rlatv > rlatu ', rlatv(j), rlatu(j), j |
561 |
STOP 18 |
562 |
ENDIF |
563 |
ENDDO |
564 |
|
565 |
print *, 'Latitudes' |
566 |
print 3, minval(ylat(:jjm)) *180d0/pi, maxval(ylat(:jjm))*180d0/pi |
567 |
|
568 |
3 Format(1x, ' Au centre du zoom, la longueur de la maille est', & |
569 |
' d environ ', f0.2, ' degres ', /, & |
570 |
' alors que la maille en dehors de la zone du zoom est ', & |
571 |
"d'environ ", f0.2, ' degres ') |
572 |
|
573 |
rlatu(1) = pi / 2. |
574 |
rlatu(jjm + 1) = -rlatu(1) |
575 |
|
576 |
END SUBROUTINE fyhyp |
577 |
|
578 |
!******************************************************************** |
579 |
|
580 |
SUBROUTINE fxhyp |
581 |
|
582 |
! From LMDZ4/libf/dyn3d/fxhyp.F, version 1.2, 2005/06/03 09:11:32 |
583 |
! Author: P. Le Van, from formulas by R. Sadourny |
584 |
|
585 |
! Compute xprimm025, rlonv, xprimv, rlonu, xprimu, xprimp025, |
586 |
! using clon, grossismx, dzoomx, taux. |
587 |
|
588 |
! Calcule les longitudes et dérivées dans la grille du GCM pour |
589 |
! une fonction x_f(\tilde x) à dérivée tangente hyperbolique. |
590 |
|
591 |
! Il vaut mieux avoir : grossismx \times delta < pi |
592 |
|
593 |
! Le premier point scalaire pour une grille regulière (grossismx = |
594 |
! 1) avec clon = 0 est à - 180 degrés. |
595 |
|
596 |
USE dimensions, ONLY: iim |
597 |
use nr_util, only: pi, pi_d, twopi, twopi_d, arth |
598 |
use tanh_cautious_m, only: tanh_cautious |
599 |
|
600 |
! Local: |
601 |
real rlonm025(iim + 1), rlonp025(iim + 1), d_rlonv(iim) |
602 |
REAL delta, h |
603 |
DOUBLE PRECISION, dimension(0:nmax):: xtild, fhyp, G, Xf, ffdx |
604 |
DOUBLE PRECISION beta |
605 |
INTEGER i, is2 |
606 |
DOUBLE PRECISION xmoy(nmax), fxm(nmax) |
607 |
|
608 |
!---------------------------------------------------------------------- |
609 |
|
610 |
print *, "Call sequence information: fxhyp" |
611 |
|
612 |
if (grossismx == 1.) then |
613 |
h = twopi / iim |
614 |
|
615 |
xprimm025(:iim) = h |
616 |
xprimp025(:iim) = h |
617 |
xprimv(:iim) = h |
618 |
xprimu(:iim) = h |
619 |
|
620 |
rlonv(:iim) = arth(- pi + clon, h, iim) |
621 |
rlonm025(:iim) = rlonv(:iim) - 0.25 * h |
622 |
rlonp025(:iim) = rlonv(:iim) + 0.25 * h |
623 |
rlonu(:iim) = rlonv(:iim) + 0.5 * h |
624 |
else |
625 |
delta = dzoomx * twopi_d |
626 |
xtild = arth(0d0, pi_d / nmax, nmax + 1) |
627 |
forall (i = 1:nmax) xmoy(i) = 0.5d0 * (xtild(i-1) + xtild(i)) |
628 |
|
629 |
! Compute fhyp: |
630 |
fhyp(1:nmax - 1) = tanh_cautious(taux * (delta / 2d0 & |
631 |
- xtild(1:nmax - 1)), xtild(1:nmax - 1) & |
632 |
* (pi_d - xtild(1:nmax - 1))) |
633 |
fhyp(0) = 1d0 |
634 |
fhyp(nmax) = -1d0 |
635 |
|
636 |
fxm = tanh_cautious(taux * (delta / 2d0 - xmoy), xmoy * (pi_d - xmoy)) |
637 |
|
638 |
! Compute \int_0 ^{\tilde x} F: |
639 |
|
640 |
ffdx(0) = 0d0 |
641 |
|
642 |
DO i = 1, nmax |
643 |
ffdx(i) = ffdx(i - 1) + fxm(i) * (xtild(i) - xtild(i-1)) |
644 |
END DO |
645 |
|
646 |
print *, "ffdx(nmax) = ", ffdx(nmax) |
647 |
beta = (pi_d - grossismx * ffdx(nmax)) / (pi_d - ffdx(nmax)) |
648 |
print *, "beta = ", beta |
649 |
|
650 |
IF (2d0 * beta - grossismx <= 0d0) THEN |
651 |
print *, 'Bad choice of grossismx, taux, dzoomx.' |
652 |
print *, 'Decrease dzoomx or grossismx.' |
653 |
STOP 1 |
654 |
END IF |
655 |
|
656 |
G = beta + (grossismx - beta) * fhyp |
657 |
|
658 |
Xf(:nmax - 1) = beta * xtild(:nmax - 1) + (grossismx - beta) & |
659 |
* ffdx(:nmax - 1) |
660 |
Xf(nmax) = pi_d |
661 |
|
662 |
call invert_zoom_x(beta, xf, xtild, G, rlonm025(:iim), xprimm025(:iim), & |
663 |
xuv = - 0.25d0) |
664 |
call invert_zoom_x(beta, xf, xtild, G, rlonv(:iim), xprimv(:iim), & |
665 |
xuv = 0d0) |
666 |
call invert_zoom_x(beta, xf, xtild, G, rlonu(:iim), xprimu(:iim), & |
667 |
xuv = 0.5d0) |
668 |
call invert_zoom_x(beta, xf, xtild, G, rlonp025(:iim), xprimp025(:iim), & |
669 |
xuv = 0.25d0) |
670 |
end if |
671 |
|
672 |
is2 = 0 |
673 |
|
674 |
IF (MINval(rlonm025(:iim)) < - pi - 0.1 & |
675 |
.or. MAXval(rlonm025(:iim)) > pi + 0.1) THEN |
676 |
IF (clon <= 0.) THEN |
677 |
is2 = 1 |
678 |
|
679 |
do while (rlonm025(is2) < - pi .and. is2 < iim) |
680 |
is2 = is2 + 1 |
681 |
end do |
682 |
|
683 |
if (rlonm025(is2) < - pi) then |
684 |
print *, 'Rlonm025 plus petit que - pi !' |
685 |
STOP 1 |
686 |
end if |
687 |
ELSE |
688 |
is2 = iim |
689 |
|
690 |
do while (rlonm025(is2) > pi .and. is2 > 1) |
691 |
is2 = is2 - 1 |
692 |
end do |
693 |
|
694 |
if (rlonm025(is2) > pi) then |
695 |
print *, 'Rlonm025 plus grand que pi !' |
696 |
STOP 1 |
697 |
end if |
698 |
END IF |
699 |
END IF |
700 |
|
701 |
call principal_cshift(is2, rlonm025, xprimm025) |
702 |
call principal_cshift(is2, rlonv, xprimv) |
703 |
call principal_cshift(is2, rlonu, xprimu) |
704 |
call principal_cshift(is2, rlonp025, xprimp025) |
705 |
|
706 |
forall (i = 1: iim) d_rlonv(i) = rlonv(i + 1) - rlonv(i) |
707 |
print *, "Minimum longitude step:", MINval(d_rlonv) * 180. / pi, "degrees" |
708 |
print *, "Maximum longitude step:", MAXval(d_rlonv) * 180. / pi, "degrees" |
709 |
|
710 |
! Check that rlonm025 <= rlonv <= rlonp025 <= rlonu: |
711 |
DO i = 1, iim + 1 |
712 |
IF (rlonp025(i) < rlonv(i)) THEN |
713 |
print *, 'rlonp025(', i, ') = ', rlonp025(i) |
714 |
print *, "< rlonv(", i, ") = ", rlonv(i) |
715 |
STOP 1 |
716 |
END IF |
717 |
|
718 |
IF (rlonv(i) < rlonm025(i)) THEN |
719 |
print *, 'rlonv(', i, ') = ', rlonv(i) |
720 |
print *, "< rlonm025(", i, ") = ", rlonm025(i) |
721 |
STOP 1 |
722 |
END IF |
723 |
|
724 |
IF (rlonp025(i) > rlonu(i)) THEN |
725 |
print *, 'rlonp025(', i, ') = ', rlonp025(i) |
726 |
print *, "> rlonu(", i, ") = ", rlonu(i) |
727 |
STOP 1 |
728 |
END IF |
729 |
END DO |
730 |
|
731 |
END SUBROUTINE fxhyp |
732 |
|
733 |
!******************************************************************** |
734 |
|
735 |
subroutine principal_cshift(is2, xlon, xprimm) |
736 |
|
737 |
! Add or subtract 2 pi so that xlon is near [-pi, pi], then cshift |
738 |
! so that xlon is in ascending order. Make the same cshift on |
739 |
! xprimm. Use clon. |
740 |
|
741 |
USE dimensions, ONLY: iim |
742 |
use nr_util, only: twopi |
743 |
|
744 |
integer, intent(in):: is2 |
745 |
real, intent(inout):: xlon(:), xprimm(:) ! (iim + 1) |
746 |
|
747 |
!----------------------------------------------------- |
748 |
|
749 |
if (is2 /= 0) then |
750 |
IF (clon <= 0.) THEN |
751 |
IF (is2 /= 1) THEN |
752 |
xlon(:is2 - 1) = xlon(:is2 - 1) + twopi |
753 |
xlon(:iim) = cshift(xlon(:iim), shift = is2 - 1) |
754 |
xprimm(:iim) = cshift(xprimm(:iim), shift = is2 - 1) |
755 |
END IF |
756 |
else |
757 |
xlon(is2 + 1:iim) = xlon(is2 + 1:iim) - twopi |
758 |
xlon(:iim) = cshift(xlon(:iim), shift = is2) |
759 |
xprimm(:iim) = cshift(xprimm(:iim), shift = is2) |
760 |
end IF |
761 |
end if |
762 |
|
763 |
xlon(iim + 1) = xlon(1) + twopi |
764 |
xprimm(iim + 1) = xprimm(1) |
765 |
|
766 |
end subroutine principal_cshift |
767 |
|
768 |
!********************************************************************** |
769 |
|
770 |
subroutine invert_zoom_x(beta, xf, xtild, G, xlon, xprim, xuv) |
771 |
|
772 |
! Using clon and grossismx. |
773 |
|
774 |
use coefpoly_m, only: coefpoly, a1, a2, a3 |
775 |
USE dimensions, ONLY: iim |
776 |
use nr_util, only: pi_d, twopi_d |
777 |
use numer_rec_95, only: hunt, rtsafe |
778 |
|
779 |
DOUBLE PRECISION, intent(in):: beta, Xf(0:), xtild(0:), G(0:) ! (0:nmax) |
780 |
|
781 |
real, intent(out):: xlon(:), xprim(:) ! (iim) |
782 |
|
783 |
DOUBLE PRECISION, intent(in):: xuv |
784 |
! between - 0.25 and 0.5 |
785 |
! 0. si calcul aux points scalaires |
786 |
! 0.5 si calcul aux points U |
787 |
|
788 |
! Local: |
789 |
DOUBLE PRECISION Y |
790 |
DOUBLE PRECISION h ! step of the uniform grid |
791 |
integer i, it |
792 |
|
793 |
DOUBLE PRECISION xvrai(iim), Gvrai(iim) |
794 |
! intermediary variables because xlon and xprim are single precision |
795 |
|
796 |
!------------------------------------------------------------------ |
797 |
|
798 |
print *, "Call sequence information: invert_zoom_x" |
799 |
it = 0 ! initial guess |
800 |
h = twopi_d / iim |
801 |
|
802 |
DO i = 1, iim |
803 |
Y = - pi_d + (i + xuv - 0.75d0) * h |
804 |
! - pi <= y < pi |
805 |
abs_y = abs(y) |
806 |
|
807 |
! Distinguish boundaries in order to avoid roundoff error. |
808 |
! funcd should be exactly equal to 0 at xtild(it) or xtild(it + |
809 |
! 1) and could be very small with the wrong sign so rtsafe |
810 |
! would fail. |
811 |
if (abs_y == 0d0) then |
812 |
xvrai(i) = 0d0 |
813 |
gvrai(i) = grossismx |
814 |
else if (abs_y == pi_d) then |
815 |
xvrai(i) = pi_d |
816 |
gvrai(i) = 2d0 * beta - grossismx |
817 |
else |
818 |
call hunt(xf, abs_y, it, my_lbound = 0) |
819 |
! {0 <= it <= nmax - 1} |
820 |
|
821 |
! Calcul de xvrai(i) et Gvrai(i) |
822 |
CALL coefpoly(Xf(it), Xf(it + 1), G(it), G(it + 1), xtild(it), & |
823 |
xtild(it + 1)) |
824 |
xvrai(i) = rtsafe(funcd, xtild(it), xtild(it + 1), xacc = 1d-6) |
825 |
Gvrai(i) = a1 + xvrai(i) * (2d0 * a2 + xvrai(i) * 3d0 * a3) |
826 |
end if |
827 |
|
828 |
if (y < 0d0) xvrai(i) = - xvrai(i) |
829 |
end DO |
830 |
|
831 |
DO i = 1, iim -1 |
832 |
IF (xvrai(i + 1) < xvrai(i)) THEN |
833 |
print *, 'xvrai(', i + 1, ') < xvrai(', i, ')' |
834 |
STOP 1 |
835 |
END IF |
836 |
END DO |
837 |
|
838 |
xlon = xvrai + clon |
839 |
xprim = h / Gvrai |
840 |
|
841 |
end subroutine invert_zoom_x |
842 |
|
843 |
!********************************************************************** |
844 |
|
845 |
SUBROUTINE funcd(x, fval, fderiv) |
846 |
|
847 |
use coefpoly_m, only: a0, a1, a2, a3 |
848 |
|
849 |
DOUBLE PRECISION, INTENT(IN):: x |
850 |
DOUBLE PRECISION, INTENT(OUT):: fval, fderiv |
851 |
|
852 |
fval = a0 + x * (a1 + x * (a2 + x * a3)) - abs_y |
853 |
fderiv = a1 + x * (2d0 * a2 + x * 3d0 * a3) |
854 |
|
855 |
END SUBROUTINE funcd |
856 |
|
857 |
end module dynetat0_m |