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