/[lmdze]/trunk/libf/dyn3d/inigeom.f90
ViewVC logotype

Contents of /trunk/libf/dyn3d/inigeom.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 70 - (show annotations)
Mon Jun 24 15:39:52 2013 UTC (10 years, 11 months ago) by guez
File size: 15939 byte(s)
In procedure, "addfi" access directly the module variable "dtphys"
instead of going through an argument.

In "conflx", do not create a local variable for temperature with
reversed order of vertical levels. Instead, give an actual argument
with reversed order in "physiq".

Changed names of variables "rmd" and "rmv" from module "suphec_m" to
"md" and "mv".

In "hgardfou", print only the first temperature out of range found.

1 module inigeom_m
2
3 IMPLICIT NONE
4
5 contains
6
7 SUBROUTINE inigeom
8
9 ! Auteur : P. Le Van
10
11 ! Calcul des élongations cuij1, ..., cuij4, cvij1, ..., cvij4 aux mêmes
12 ! endroits que les aires aireij1_2d, ..., aireij4_2d.
13
14 ! Choix entre une fonction "f(y)" à dérivée sinusoïdale ou à
15 ! dérivée tangente hyperbolique. Calcul des coefficients cu_2d,
16 ! cv_2d, 1. / cu_2d**2, 1. / cv_2d**2. Les coefficients cu_2d et cv_2d
17 ! permettent de passer des vitesses naturelles aux vitesses
18 ! covariantes et contravariantes, ou vice-versa.
19
20 ! On a :
21 ! u(covariant) = cu_2d * u(naturel), u(contravariant) = u(naturel) / cu_2d
22 ! v(covariant) = cv_2d * v(naturel), v(contravariant) = v(naturel) / cv_2d
23
24 ! On en tire :
25 ! u(covariant) = cu_2d * cu_2d * u(contravariant)
26 ! v(covariant) = cv_2d * cv_2d * v(contravariant)
27
28 ! On a l'application (x(X), y(Y)) avec - im / 2 + 1 <= X <= im / 2
29 ! et - jm / 2 <= Y <= jm / 2
30
31 ! x est la longitude du point en radians.
32 ! y est la latitude du point en radians.
33 !
34 ! On a : cu_2d(i, j) = rad * cos(y) * dx / dX
35 ! cv(j) = rad * dy / dY
36 ! aire_2d(i, j) = cu_2d(i, j) * cv(j)
37 !
38 ! y, dx / dX, dy / dY calculés aux points concernés. cv, bien que
39 ! dépendant de j uniquement, sera ici indicé aussi en i pour un
40 ! adressage plus facile en ij.
41
42 ! xprimu et xprimv sont respectivement les valeurs de dx / dX aux
43 ! points u et v. yprimu et yprimv sont respectivement les valeurs
44 ! de dy / dY aux points u et v. rlatu et rlatv sont respectivement
45 ! les valeurs de la latitude aux points u et v. cvu et cv_2d sont
46 ! respectivement les valeurs de cv_2d aux points u et v.
47
48 ! cu_2d, cuv, cuscal, cuz sont respectivement les valeurs de cu_2d
49 ! aux points u, v, scalaires, et z. Cf. "inigeom.txt".
50
51 USE comconst, ONLY : g, omeg, rad
52 USE comgeom, ONLY : airesurg_2d, aireu_2d, airev_2d, aire_2d, &
53 alpha1p2_2d, alpha1p4_2d, alpha1_2d, &
54 alpha2p3_2d, alpha2_2d, alpha3p4_2d, alpha3_2d, alpha4_2d, apoln, &
55 apols, constang_2d, cuscvugam_2d, cusurcvu_2d, cuvscvgam1_2d, &
56 cuvscvgam2_2d, cuvsurcv_2d, cu_2d, cvscuvgam_2d, cvsurcuv_2d, &
57 cvuscugam1_2d, cvuscugam2_2d, cvusurcu_2d, cv_2d, fext_2d, rlatu, &
58 rlatv, rlonu, rlonv, unsairez_2d, unsaire_2d, unsairz_gam_2d, &
59 unsair_gam1_2d, unsair_gam2_2d, unsapolnga1, unsapolnga2, &
60 unsapolsga1, unsapolsga2, unscu2_2d, unscv2_2d, xprimu, xprimv
61 USE comdissnew, ONLY : coefdis, nitergdiv, nitergrot, niterh
62 use conf_gcm_m, ONLY : fxyhypb, ysinus
63 USE dimens_m, ONLY : iim, jjm
64 use fxy_m, only: fxy
65 use fxyhyper_m, only: fxyhyper
66 use jumble, only: new_unit
67 use nr_util, only: pi
68 USE paramet_m, ONLY : iip1, jjp1
69 USE serre, ONLY : alphax, alphay, clat, clon, dzoomx, dzoomy, grossismx, &
70 grossismy, pxo, pyo, taux, tauy, transx, transy
71
72 ! Variables locales
73
74 INTEGER i, j, itmax, itmay, iter, unit
75 REAL cvu(iip1, jjp1), cuv(iip1, jjm)
76 REAL ai14, ai23, airez, un4rad2
77 REAL eps, x1, xo1, f, df, xdm, y1, yo1, ydm
78 REAL coslatm, coslatp, radclatm, radclatp
79 REAL, dimension(iip1, jjp1):: cuij1, cuij2, cuij3, cuij4 ! in m
80 REAL, dimension(iip1, jjp1):: cvij1, cvij2, cvij3, cvij4 ! in m
81 REAL rlatu1(jjm), yprimu1(jjm), rlatu2(jjm), yprimu2(jjm)
82 real yprimv(jjm), yprimu(jjp1)
83 REAL gamdi_gdiv, gamdi_grot, gamdi_h
84 REAL rlonm025(iip1), xprimm025(iip1), rlonp025(iip1), xprimp025(iip1)
85 real, dimension(iim + 1, jjm + 1):: aireij1_2d, aireij2_2d, aireij3_2d, &
86 aireij4_2d ! in m2
87 real airuscv2_2d(iim + 1, jjm)
88 real airvscu2_2d(iim + 1, jjm), aiuscv2gam_2d(iim + 1, jjm)
89 real aivscu2gam_2d(iim + 1, jjm)
90
91 !------------------------------------------------------------------
92
93 PRINT *, 'Call sequence information: inigeom'
94
95 IF (nitergdiv/=2) THEN
96 gamdi_gdiv = coefdis / (real(nitergdiv)-2.)
97 ELSE
98 gamdi_gdiv = 0.
99 END IF
100 IF (nitergrot/=2) THEN
101 gamdi_grot = coefdis / (real(nitergrot)-2.)
102 ELSE
103 gamdi_grot = 0.
104 END IF
105 IF (niterh/=2) THEN
106 gamdi_h = coefdis / (real(niterh)-2.)
107 ELSE
108 gamdi_h = 0.
109 END IF
110
111 print *, 'gamdi_gdiv = ', gamdi_gdiv
112 print *, "gamdi_grot = ", gamdi_grot
113 print *, "gamdi_h = ", gamdi_h
114
115 IF (.NOT. fxyhypb) THEN
116 IF (ysinus) THEN
117 print *, ' Inigeom, Y = Sinus (Latitude) '
118 ! utilisation de f(x, y) avec y = sinus de la latitude
119 CALL fxysinus(rlatu, yprimu, rlatv, yprimv, rlatu1, yprimu1, &
120 rlatu2, yprimu2, rlonu, xprimu, rlonv, xprimv, rlonm025, &
121 xprimm025, rlonp025, xprimp025)
122 ELSE
123 print *, 'Inigeom, Y = Latitude, der. sinusoid .'
124 ! utilisation de f(x, y) a tangente sinusoidale, y etant la latit
125
126 pxo = clon * pi / 180.
127 pyo = 2. * clat * pi / 180.
128
129 ! determination de transx (pour le zoom) par Newton-Raphson
130
131 itmax = 10
132 eps = .1E-7
133
134 xo1 = 0.
135 DO iter = 1, itmax
136 x1 = xo1
137 f = x1 + alphax * sin(x1-pxo)
138 df = 1. + alphax * cos(x1-pxo)
139 x1 = x1 - f / df
140 xdm = abs(x1-xo1)
141 IF (xdm<=eps) EXIT
142 xo1 = x1
143 END DO
144
145 transx = xo1
146
147 itmay = 10
148 eps = .1E-7
149
150 yo1 = 0.
151 DO iter = 1, itmay
152 y1 = yo1
153 f = y1 + alphay * sin(y1-pyo)
154 df = 1. + alphay * cos(y1-pyo)
155 y1 = y1 - f / df
156 ydm = abs(y1-yo1)
157 IF (ydm<=eps) EXIT
158 yo1 = y1
159 END DO
160
161 transy = yo1
162
163 CALL fxy(rlatu, yprimu, rlatv, yprimv, rlatu1, yprimu1, rlatu2, &
164 yprimu2, rlonu, xprimu, rlonv, xprimv, rlonm025, xprimm025, &
165 rlonp025, xprimp025)
166 END IF
167 ELSE
168 ! Utilisation de fxyhyper, f(x, y) à dérivée tangente hyperbolique
169 print *, 'Inigeom, Y = Latitude, dérivée tangente hyperbolique'
170 CALL fxyhyper(clat, grossismy, dzoomy, tauy, clon, grossismx, dzoomx, &
171 taux, rlatu, yprimu, rlatv, yprimv, rlatu1, yprimu1, rlatu2, &
172 yprimu2, rlonu, xprimu, rlonv, xprimv, rlonm025, xprimm025, &
173 rlonp025, xprimp025)
174 END IF
175
176 rlatu(1) = pi / 2.
177 rlatu(jjp1) = -rlatu(1)
178
179 ! Calcul aux pôles
180
181 yprimu(1) = 0.
182 yprimu(jjp1) = 0.
183
184 un4rad2 = 0.25 * rad * rad
185
186 ! Cf. "inigeom.txt". Calcul des quatre aires élémentaires
187 ! aireij1_2d, aireij2_2d, aireij3_2d, aireij4_2d qui entourent
188 ! chaque aire_2d(i, j), ainsi que les quatre élongations
189 ! élémentaires cuij et les quatre élongations cvij qui sont
190 ! calculées aux mêmes endroits que les aireij.
191
192 coslatm = cos(rlatu1(1))
193 radclatm = 0.5 * rad * coslatm
194
195 aireij1_2d(:iim, 1) = 0.
196 aireij2_2d(:iim, 1) = un4rad2 * coslatm * xprimp025(:iim) * yprimu1(1)
197 aireij3_2d(:iim, 1) = un4rad2 * coslatm * xprimm025(:iim) * yprimu1(1)
198 aireij4_2d(:iim, 1) = 0.
199
200 cuij1(:iim, 1) = 0.
201 cuij2(:iim, 1) = radclatm * xprimp025(:iim)
202 cuij3(:iim, 1) = radclatm * xprimm025(:iim)
203 cuij4(:iim, 1) = 0.
204
205 cvij1(:iim, 1) = 0.
206 cvij2(:iim, 1) = 0.5 * rad * yprimu1(1)
207 cvij3(:iim, 1) = cvij2(:iim, 1)
208 cvij4(:iim, 1) = 0.
209
210 do j = 2, jjm
211 coslatm = cos(rlatu1(j))
212 coslatp = cos(rlatu2(j-1))
213 radclatp = 0.5 * rad * coslatp
214 radclatm = 0.5 * rad * coslatm
215 ai14 = un4rad2 * coslatp * yprimu2(j-1)
216 ai23 = un4rad2 * coslatm * yprimu1(j)
217
218 aireij1_2d(:iim, j) = ai14 * xprimp025(:iim)
219 aireij2_2d(:iim, j) = ai23 * xprimp025(:iim)
220 aireij3_2d(:iim, j) = ai23 * xprimm025(:iim)
221 aireij4_2d(:iim, j) = ai14 * xprimm025(:iim)
222 cuij1(:iim, j) = radclatp * xprimp025(:iim)
223 cuij2(:iim, j) = radclatm * xprimp025(:iim)
224 cuij3(:iim, j) = radclatm * xprimm025(:iim)
225 cuij4(:iim, j) = radclatp * xprimm025(:iim)
226 cvij1(:iim, j) = 0.5 * rad * yprimu2(j-1)
227 cvij2(:iim, j) = 0.5 * rad * yprimu1(j)
228 cvij3(:iim, j) = cvij2(:iim, j)
229 cvij4(:iim, j) = cvij1(:iim, j)
230 end do
231
232 coslatp = cos(rlatu2(jjm))
233 radclatp = 0.5 * rad * coslatp
234
235 aireij1_2d(:iim, jjp1) = un4rad2 * coslatp * xprimp025(:iim) * yprimu2(jjm)
236 aireij2_2d(:iim, jjp1) = 0.
237 aireij3_2d(:iim, jjp1) = 0.
238 aireij4_2d(:iim, jjp1) = un4rad2 * coslatp * xprimm025(:iim) * yprimu2(jjm)
239
240 cuij1(:iim, jjp1) = radclatp * xprimp025(:iim)
241 cuij2(:iim, jjp1) = 0.
242 cuij3(:iim, jjp1) = 0.
243 cuij4(:iim, jjp1) = radclatp * xprimm025(:iim)
244
245 cvij1(:iim, jjp1) = 0.5 * rad * yprimu2(jjm)
246 cvij2(:iim, jjp1) = 0.
247 cvij3(:iim, jjp1) = 0.
248 cvij4(:iim, jjp1) = cvij1(:iim, jjp1)
249
250 ! Périodicité :
251
252 cvij1(iip1, :) = cvij1(1, :)
253 cvij2(iip1, :) = cvij2(1, :)
254 cvij3(iip1, :) = cvij3(1, :)
255 cvij4(iip1, :) = cvij4(1, :)
256
257 cuij1(iip1, :) = cuij1(1, :)
258 cuij2(iip1, :) = cuij2(1, :)
259 cuij3(iip1, :) = cuij3(1, :)
260 cuij4(iip1, :) = cuij4(1, :)
261
262 aireij1_2d(iip1, :) = aireij1_2d(1, :)
263 aireij2_2d(iip1, :) = aireij2_2d(1, :)
264 aireij3_2d(iip1, :) = aireij3_2d(1, :)
265 aireij4_2d(iip1, :) = aireij4_2d(1, :)
266
267 DO j = 1, jjp1
268 DO i = 1, iim
269 aire_2d(i, j) = aireij1_2d(i, j) + aireij2_2d(i, j) &
270 + aireij3_2d(i, j) + aireij4_2d(i, j)
271 alpha1_2d(i, j) = aireij1_2d(i, j) / aire_2d(i, j)
272 alpha2_2d(i, j) = aireij2_2d(i, j) / aire_2d(i, j)
273 alpha3_2d(i, j) = aireij3_2d(i, j) / aire_2d(i, j)
274 alpha4_2d(i, j) = aireij4_2d(i, j) / aire_2d(i, j)
275 alpha1p2_2d(i, j) = alpha1_2d(i, j) + alpha2_2d(i, j)
276 alpha1p4_2d(i, j) = alpha1_2d(i, j) + alpha4_2d(i, j)
277 alpha2p3_2d(i, j) = alpha2_2d(i, j) + alpha3_2d(i, j)
278 alpha3p4_2d(i, j) = alpha3_2d(i, j) + alpha4_2d(i, j)
279 END DO
280
281 aire_2d(iip1, j) = aire_2d(1, j)
282 alpha1_2d(iip1, j) = alpha1_2d(1, j)
283 alpha2_2d(iip1, j) = alpha2_2d(1, j)
284 alpha3_2d(iip1, j) = alpha3_2d(1, j)
285 alpha4_2d(iip1, j) = alpha4_2d(1, j)
286 alpha1p2_2d(iip1, j) = alpha1p2_2d(1, j)
287 alpha1p4_2d(iip1, j) = alpha1p4_2d(1, j)
288 alpha2p3_2d(iip1, j) = alpha2p3_2d(1, j)
289 alpha3p4_2d(iip1, j) = alpha3p4_2d(1, j)
290 END DO
291
292 DO j = 1, jjp1
293 DO i = 1, iim
294 aireu_2d(i, j) = aireij1_2d(i, j) + aireij2_2d(i, j) + &
295 aireij4_2d(i + 1, j) + aireij3_2d(i + 1, j)
296 unsaire_2d(i, j) = 1. / aire_2d(i, j)
297 unsair_gam1_2d(i, j) = unsaire_2d(i, j)**(-gamdi_gdiv)
298 unsair_gam2_2d(i, j) = unsaire_2d(i, j)**(-gamdi_h)
299 airesurg_2d(i, j) = aire_2d(i, j) / g
300 END DO
301 aireu_2d(iip1, j) = aireu_2d(1, j)
302 unsaire_2d(iip1, j) = unsaire_2d(1, j)
303 unsair_gam1_2d(iip1, j) = unsair_gam1_2d(1, j)
304 unsair_gam2_2d(iip1, j) = unsair_gam2_2d(1, j)
305 airesurg_2d(iip1, j) = airesurg_2d(1, j)
306 END DO
307
308 DO j = 1, jjm
309 DO i = 1, iim
310 airev_2d(i, j) = aireij2_2d(i, j) + aireij3_2d(i, j) + &
311 aireij1_2d(i, j + 1) + aireij4_2d(i, j + 1)
312 END DO
313 DO i = 1, iim
314 airez = aireij2_2d(i, j) + aireij1_2d(i, j + 1) &
315 + aireij3_2d(i + 1, j) + aireij4_2d(i + 1, j + 1)
316 unsairez_2d(i, j) = 1. / airez
317 unsairz_gam_2d(i, j) = unsairez_2d(i, j)**(-gamdi_grot)
318 fext_2d(i, j) = airez * sin(rlatv(j)) * 2. * omeg
319 END DO
320 airev_2d(iip1, j) = airev_2d(1, j)
321 unsairez_2d(iip1, j) = unsairez_2d(1, j)
322 fext_2d(iip1, j) = fext_2d(1, j)
323 unsairz_gam_2d(iip1, j) = unsairz_gam_2d(1, j)
324 END DO
325
326 ! Calcul des élongations cu_2d, cv_2d, cvu
327
328 DO j = 1, jjm
329 DO i = 1, iim
330 cv_2d(i, j) = 0.5 * &
331 (cvij2(i, j) + cvij3(i, j) + cvij1(i, j + 1) + cvij4(i, j + 1))
332 cvu(i, j) = 0.5 * (cvij1(i, j) + cvij4(i, j) + cvij2(i, j) &
333 + cvij3(i, j))
334 cuv(i, j) = 0.5 * (cuij2(i, j) + cuij3(i, j) + cuij1(i, j + 1) &
335 + cuij4(i, j + 1))
336 unscv2_2d(i, j) = 1. / cv_2d(i, j)**2
337 END DO
338 DO i = 1, iim
339 cuvsurcv_2d(i, j) = airev_2d(i, j) * unscv2_2d(i, j)
340 cvsurcuv_2d(i, j) = 1. / cuvsurcv_2d(i, j)
341 cuvscvgam1_2d(i, j) = cuvsurcv_2d(i, j)**(-gamdi_gdiv)
342 cuvscvgam2_2d(i, j) = cuvsurcv_2d(i, j)**(-gamdi_h)
343 cvscuvgam_2d(i, j) = cvsurcuv_2d(i, j)**(-gamdi_grot)
344 END DO
345 cv_2d(iip1, j) = cv_2d(1, j)
346 cvu(iip1, j) = cvu(1, j)
347 unscv2_2d(iip1, j) = unscv2_2d(1, j)
348 cuv(iip1, j) = cuv(1, j)
349 cuvsurcv_2d(iip1, j) = cuvsurcv_2d(1, j)
350 cvsurcuv_2d(iip1, j) = cvsurcuv_2d(1, j)
351 cuvscvgam1_2d(iip1, j) = cuvscvgam1_2d(1, j)
352 cuvscvgam2_2d(iip1, j) = cuvscvgam2_2d(1, j)
353 cvscuvgam_2d(iip1, j) = cvscuvgam_2d(1, j)
354 END DO
355
356 DO j = 2, jjm
357 DO i = 1, iim
358 cu_2d(i, j) = 0.5 * (cuij1(i, j) + cuij4(i + 1, j) + cuij2(i, j) &
359 + cuij3(i + 1, j))
360 unscu2_2d(i, j) = 1. / cu_2d(i, j)**2
361 cvusurcu_2d(i, j) = aireu_2d(i, j) * unscu2_2d(i, j)
362 cusurcvu_2d(i, j) = 1. / cvusurcu_2d(i, j)
363 cvuscugam1_2d(i, j) = cvusurcu_2d(i, j)**(-gamdi_gdiv)
364 cvuscugam2_2d(i, j) = cvusurcu_2d(i, j)**(-gamdi_h)
365 cuscvugam_2d(i, j) = cusurcvu_2d(i, j)**(-gamdi_grot)
366 END DO
367 cu_2d(iip1, j) = cu_2d(1, j)
368 unscu2_2d(iip1, j) = unscu2_2d(1, j)
369 cvusurcu_2d(iip1, j) = cvusurcu_2d(1, j)
370 cusurcvu_2d(iip1, j) = cusurcvu_2d(1, j)
371 cvuscugam1_2d(iip1, j) = cvuscugam1_2d(1, j)
372 cvuscugam2_2d(iip1, j) = cvuscugam2_2d(1, j)
373 cuscvugam_2d(iip1, j) = cuscvugam_2d(1, j)
374 END DO
375
376 ! Calcul aux pôles
377
378 cu_2d(:, 1) = 0.
379 unscu2_2d(:, 1) = 0.
380 cvu(:, 1) = 0.
381
382 cu_2d(:, jjp1) = 0.
383 unscu2_2d(:, jjp1) = 0.
384 cvu(:, jjp1) = 0.
385
386 DO j = 1, jjm
387 DO i = 1, iim
388 airvscu2_2d(i, j) = airev_2d(i, j) / (cuv(i, j) * cuv(i, j))
389 aivscu2gam_2d(i, j) = airvscu2_2d(i, j)**(-gamdi_grot)
390 END DO
391 airvscu2_2d(iip1, j) = airvscu2_2d(1, j)
392 aivscu2gam_2d(iip1, j) = aivscu2gam_2d(1, j)
393 END DO
394
395 DO j = 2, jjm
396 DO i = 1, iim
397 airuscv2_2d(i, j) = aireu_2d(i, j) / (cvu(i, j) * cvu(i, j))
398 aiuscv2gam_2d(i, j) = airuscv2_2d(i, j)**(-gamdi_grot)
399 END DO
400 airuscv2_2d(iip1, j) = airuscv2_2d(1, j)
401 aiuscv2gam_2d(iip1, j) = aiuscv2gam_2d(1, j)
402 END DO
403
404 ! Calcul des aires aux pôles :
405
406 apoln = sum(aire_2d(:iim, 1))
407 apols = sum(aire_2d(:iim, jjp1))
408 unsapolnga1 = 1. / (apoln**(-gamdi_gdiv))
409 unsapolsga1 = 1. / (apols**(-gamdi_gdiv))
410 unsapolnga2 = 1. / (apoln**(-gamdi_h))
411 unsapolsga2 = 1. / (apols**(-gamdi_h))
412
413 ! Changement F. Hourdin calcul conservatif pour fext_2d
414 ! constang_2d contient le produit a * cos (latitude) * omega
415
416 DO i = 1, iim
417 constang_2d(i, 1) = 0.
418 END DO
419 DO j = 1, jjm - 1
420 DO i = 1, iim
421 constang_2d(i, j + 1) = rad * omeg * cu_2d(i, j + 1) &
422 * cos(rlatu(j + 1))
423 END DO
424 END DO
425 DO i = 1, iim
426 constang_2d(i, jjp1) = 0.
427 END DO
428
429 ! Périodicité en longitude
430
431 DO j = 1, jjm
432 fext_2d(iip1, j) = fext_2d(1, j)
433 END DO
434 DO j = 1, jjp1
435 constang_2d(iip1, j) = constang_2d(1, j)
436 END DO
437
438 call new_unit(unit)
439 open(unit, file="longitude_latitude.txt", status="replace", action="write")
440 write(unit, fmt=*) '"longitudes at V points (degrees)"', rlonv * 180. / pi
441 write(unit, fmt=*) '"latitudes at V points (degrees)"', rlatv * 180. / pi
442 write(unit, fmt=*) '"longitudes at U points (degrees)"', rlonu * 180. / pi
443 write(unit, fmt=*) '"latitudes at U points (degrees)"', rlatu * 180. / pi
444 close(unit)
445
446 END SUBROUTINE inigeom
447
448 end module inigeom_m

  ViewVC Help
Powered by ViewVC 1.1.21