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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 39 - (show annotations)
Tue Jan 25 15:11:05 2011 UTC (13 years, 3 months ago) by guez
File size: 18315 byte(s)
"pi" comes from "nr_util". Removed subroutine "initialize" in module
"comconst".

Copied the content of "fxy_sin.h" into "fxysinus", instead of getting
it from an "include" line. Removed file "fxy_sin.h".

"ps" has rank 2 in "gcm" and "dynetat0".

Assumed-shape for argument "q" of "integrd".

1 module inigeom_m
2
3 IMPLICIT NONE
4
5 contains
6
7 SUBROUTINE inigeom
8
9 ! Auteur : P. Le Van
10 ! Version du 01/04/2001
11
12 ! Calcul des élongations cuij1, ..., cuij4, cvij1, ..., cvij4 aux mêmes
13 ! endroits que les aires aireij1_2d, ..., aireij4_2d.
14
15 ! Choix entre une fonction "f(y)" à dérivée sinusoïdale ou à dérivée
16 ! tangente hyperbolique
17 ! calcul des coefficients (cu_2d, cv_2d, 1./cu_2d**2, 1./cv_2d**2)
18
19 ! les coef. (cu_2d, cv_2d) permettent de passer des vitesses naturelles
20 ! aux vitesses covariantes et contravariantes, ou vice-versa
21
22 ! on a :
23 ! u (covariant) = cu_2d * u (naturel), u(contrav)= u(nat)/cu_2d
24 ! v (covariant) = cv_2d * v (naturel), v(contrav)= v(nat)/cv_2d
25
26 ! on en tire :
27 ! u(covariant) = cu_2d * cu_2d * u(contravariant)
28 ! v(covariant) = cv_2d * cv_2d * v(contravariant)
29
30 ! on a l'application (x(X), y(Y)) avec - im/2 +1 <= X <= im/2
31 ! et - jm/2 <= Y <= jm/2
32
33 ! x est la longitude du point en radians.
34 ! y est la latitude du point en radians.
35 !
36 ! on a : cu_2d(i, j) = rad * cos(y) * dx/dX
37 ! cv(j) = rad * dy/dY
38 ! aire_2d(i, j) = cu_2d(i, j) * cv(j)
39 !
40 ! y, dx/dX, dy/dY calcules aux points concernes
41 ! cv, bien que dependant de j uniquement, sera ici indice aussi en i
42 ! pour un adressage plus facile en ij.
43
44 ! aux points u et v,
45 ! xprimu et xprimv sont respectivement les valeurs de dx/dX
46 ! yprimu et yprimv sont respectivement les valeurs de dy/dY
47 ! rlatu et rlatv sont respectivement les valeurs de la latitude
48 ! cvu et cv_2d sont respectivement les valeurs de cv_2d
49
50 ! aux points u, v, scalaires, et z
51 ! cu_2d, cuv, cuscal, cuz sont respectivement les valeurs de cu_2d
52 ! Cf. "inigeom.txt".
53
54 USE comconst, ONLY : g, omeg, rad
55 USE comgeom, ONLY : airesurg_2d, aireu_2d, airev_2d, aire_2d, &
56 alpha1p2_2d, alpha1p4_2d, alpha1_2d, &
57 alpha2p3_2d, alpha2_2d, alpha3p4_2d, alpha3_2d, alpha4_2d, apoln, &
58 apols, constang_2d, cuscvugam_2d, cusurcvu_2d, cuvscvgam1_2d, &
59 cuvscvgam2_2d, cuvsurcv_2d, cu_2d, cvscuvgam_2d, cvsurcuv_2d, &
60 cvuscugam1_2d, cvuscugam2_2d, cvusurcu_2d, cv_2d, fext_2d, rlatu, &
61 rlatv, rlonu, rlonv, unsairez_2d, unsaire_2d, unsairz_gam_2d, &
62 unsair_gam1_2d, unsair_gam2_2d, unsapolnga1, unsapolnga2, &
63 unsapolsga1, unsapolsga2, unscu2_2d, unscv2_2d, xprimu, xprimv
64 USE comdissnew, ONLY : coefdis, nitergdiv, nitergrot, niterh
65 USE dimens_m, ONLY : iim, jjm
66 USE logic, ONLY : fxyhypb, ysinus
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
75 REAL cvu(iip1, jjp1), cuv(iip1, jjm)
76 REAL ai14, ai23, airez, rlatp, rlatm, xprm, xprp, un4rad2, yprp, yprm
77 REAL eps, x1, xo1, f, df, xdm, y1, yo1, ydm
78 REAL coslatm, coslatp, radclatm, radclatp
79 REAL cuij1(iip1, jjp1), cuij2(iip1, jjp1), cuij3(iip1, jjp1), &
80 cuij4(iip1, jjp1)
81 REAL cvij1(iip1, jjp1), cvij2(iip1, jjp1), cvij3(iip1, jjp1), &
82 cvij4(iip1, jjp1)
83 REAL rlonvv(iip1), rlatuu(jjp1)
84 REAL rlatu1(jjm), yprimu1(jjm), rlatu2(jjm), yprimu2(jjm), yprimv(jjm), &
85 yprimu(jjp1)
86 REAL gamdi_gdiv, gamdi_grot, gamdi_h
87
88 REAL rlonm025(iip1), xprimm025(iip1), rlonp025(iip1), xprimp025(iip1)
89 SAVE rlatu1, yprimu1, rlatu2, yprimu2, yprimv, yprimu
90 SAVE rlonm025, xprimm025, rlonp025, xprimp025
91
92 real aireij1_2d(iim + 1, jjm + 1)
93 real aireij2_2d(iim + 1, jjm + 1)
94 real aireij3_2d(iim + 1, jjm + 1), aireij4_2d(iim + 1, jjm + 1)
95 real airuscv2_2d(iim + 1, jjm)
96 real airvscu2_2d(iim + 1, jjm), aiuscv2gam_2d(iim + 1, jjm)
97 real aivscu2gam_2d(iim + 1, jjm)
98
99 !------------------------------------------------------------------
100
101 PRINT *, 'Call sequence information: inigeom'
102
103 IF (nitergdiv/=2) THEN
104 gamdi_gdiv = coefdis/(real(nitergdiv)-2.)
105 ELSE
106 gamdi_gdiv = 0.
107 END IF
108 IF (nitergrot/=2) THEN
109 gamdi_grot = coefdis/(real(nitergrot)-2.)
110 ELSE
111 gamdi_grot = 0.
112 END IF
113 IF (niterh/=2) THEN
114 gamdi_h = coefdis/(real(niterh)-2.)
115 ELSE
116 gamdi_h = 0.
117 END IF
118
119 print *, 'gamdi_gdiv = ', gamdi_gdiv
120 print *, "gamdi_grot = ", gamdi_grot
121 print *, "gamdi_h = ", gamdi_h
122
123 WRITE (6, 990)
124
125 IF (.NOT. fxyhypb) THEN
126 IF (ysinus) THEN
127 print *, ' Inigeom, Y = Sinus (Latitude) '
128 ! utilisation de f(x, y) avec y = sinus de la latitude
129 CALL fxysinus(rlatu, yprimu, rlatv, yprimv, rlatu1, yprimu1, &
130 rlatu2, yprimu2, rlonu, xprimu, rlonv, xprimv, rlonm025, &
131 xprimm025, rlonp025, xprimp025)
132 ELSE
133 print *, 'Inigeom, Y = Latitude, der. sinusoid .'
134 ! utilisation de f(x, y) a tangente sinusoidale, y etant la latit
135
136 pxo = clon*pi/180.
137 pyo = 2.*clat*pi/180.
138
139 ! determination de transx (pour le zoom) par Newton-Raphson
140
141 itmax = 10
142 eps = .1E-7
143
144 xo1 = 0.
145 DO iter = 1, itmax
146 x1 = xo1
147 f = x1 + alphax*sin(x1-pxo)
148 df = 1. + alphax*cos(x1-pxo)
149 x1 = x1 - f/df
150 xdm = abs(x1-xo1)
151 IF (xdm<=eps) EXIT
152 xo1 = x1
153 END DO
154
155 transx = xo1
156
157 itmay = 10
158 eps = .1E-7
159
160 yo1 = 0.
161 DO iter = 1, itmay
162 y1 = yo1
163 f = y1 + alphay*sin(y1-pyo)
164 df = 1. + alphay*cos(y1-pyo)
165 y1 = y1 - f/df
166 ydm = abs(y1-yo1)
167 IF (ydm<=eps) EXIT
168 yo1 = y1
169 END DO
170
171 transy = yo1
172
173 CALL fxy(rlatu, yprimu, rlatv, yprimv, rlatu1, yprimu1, rlatu2, &
174 yprimu2, rlonu, xprimu, rlonv, xprimv, rlonm025, xprimm025, &
175 rlonp025, xprimp025)
176 END IF
177 ELSE
178 ! Utilisation de fxyhyper, f(x, y) à dérivée tangente hyperbolique
179 print *, 'Inigeom, Y = Latitude, dérivée tangente hyperbolique'
180 CALL fxyhyper(clat, grossismy, dzoomy, tauy, clon, grossismx, dzoomx, &
181 taux, rlatu, yprimu, rlatv, yprimv, rlatu1, yprimu1, rlatu2, &
182 yprimu2, rlonu, xprimu, rlonv, xprimv, rlonm025, xprimm025, &
183 rlonp025, xprimp025)
184 END IF
185
186 rlatu(1) = asin(1.)
187 rlatu(jjp1) = -rlatu(1)
188
189 ! calcul aux poles
190
191 yprimu(1) = 0.
192 yprimu(jjp1) = 0.
193
194 un4rad2 = 0.25*rad*rad
195
196 ! calcul des aires (aire_2d, aireu_2d, airev_2d, 1./aire_2d, 1./airez)
197 ! - et de fext_2d, force de coriolis extensive
198
199 ! A 1 point scalaire P (i, j) de la grille, reguliere en (X, Y), sont
200 ! affectees 4 aires entourant P, calculees respectivement aux points
201 ! (i + 1/4, j - 1/4) : aireij1_2d (i, j)
202 ! (i + 1/4, j + 1/4) : aireij2_2d (i, j)
203 ! (i - 1/4, j + 1/4) : aireij3_2d (i, j)
204 ! (i - 1/4, j - 1/4) : aireij4_2d (i, j)
205
206 !,
207 ! Les cotes de chacun de ces 4 carres etant egaux a 1/2 suivant (X, Y).
208 ! Chaque aire centree en 1 point scalaire P(i, j) est egale a la somme
209 ! des 4 aires aireij1_2d, aireij2_2d, aireij3_2d, aireij4_2d qui sont
210 ! affectees au
211 ! point (i, j).
212 ! On definit en outre les coefficients alpha comme etant egaux a
213 ! (aireij / aire_2d), c.a.d par exp.
214 ! alpha1_2d(i, j)=aireij1_2d(i, j)/aire_2d(i, j)
215
216 ! De meme, toute aire centree en 1 point U est egale a la somme des
217 ! 4 aires aireij1_2d, aireij2_2d, aireij3_2d, aireij4_2d entourant
218 ! le point U.
219 ! Idem pour airev_2d, airez.
220
221 ! On a, pour chaque maille : dX = dY = 1
222
223 ! V
224
225 ! aireij4_2d . . aireij1_2d
226
227 ! U . . P . U
228
229 ! aireij3_2d . . aireij2_2d
230
231 ! V
232
233 ! Calcul des 4 aires elementaires aireij1_2d, aireij2_2d,
234 ! aireij3_2d, aireij4_2d
235 ! qui entourent chaque aire_2d(i, j), ainsi que les 4 elongations
236 ! elementaires
237 ! cuij et les 4 elongat. cvij qui sont calculees aux memes
238 ! endroits que les aireij.
239
240 ! do 35 : boucle sur les jjm + 1 latitudes
241
242 DO j = 1, jjp1
243
244 IF (j==1) THEN
245
246 yprm = yprimu1(j)
247 rlatm = rlatu1(j)
248
249 coslatm = cos(rlatm)
250 radclatm = 0.5*rad*coslatm
251
252 DO i = 1, iim
253 xprp = xprimp025(i)
254 xprm = xprimm025(i)
255 aireij2_2d(i, 1) = un4rad2*coslatm*xprp*yprm
256 aireij3_2d(i, 1) = un4rad2*coslatm*xprm*yprm
257 cuij2(i, 1) = radclatm*xprp
258 cuij3(i, 1) = radclatm*xprm
259 cvij2(i, 1) = 0.5*rad*yprm
260 cvij3(i, 1) = cvij2(i, 1)
261 END DO
262
263 DO i = 1, iim
264 aireij1_2d(i, 1) = 0.
265 aireij4_2d(i, 1) = 0.
266 cuij1(i, 1) = 0.
267 cuij4(i, 1) = 0.
268 cvij1(i, 1) = 0.
269 cvij4(i, 1) = 0.
270 END DO
271
272 END IF
273
274 IF (j==jjp1) THEN
275 yprp = yprimu2(j-1)
276 rlatp = rlatu2(j-1)
277
278 coslatp = cos(rlatp)
279 radclatp = 0.5*rad*coslatp
280
281 DO i = 1, iim
282 xprp = xprimp025(i)
283 xprm = xprimm025(i)
284 aireij1_2d(i, jjp1) = un4rad2*coslatp*xprp*yprp
285 aireij4_2d(i, jjp1) = un4rad2*coslatp*xprm*yprp
286 cuij1(i, jjp1) = radclatp*xprp
287 cuij4(i, jjp1) = radclatp*xprm
288 cvij1(i, jjp1) = 0.5*rad*yprp
289 cvij4(i, jjp1) = cvij1(i, jjp1)
290 END DO
291
292 DO i = 1, iim
293 aireij2_2d(i, jjp1) = 0.
294 aireij3_2d(i, jjp1) = 0.
295 cvij2(i, jjp1) = 0.
296 cvij3(i, jjp1) = 0.
297 cuij2(i, jjp1) = 0.
298 cuij3(i, jjp1) = 0.
299 END DO
300
301 END IF
302
303 IF (j>1 .AND. j<jjp1) THEN
304
305 rlatp = rlatu2(j-1)
306 yprp = yprimu2(j-1)
307 rlatm = rlatu1(j)
308 yprm = yprimu1(j)
309
310 coslatm = cos(rlatm)
311 coslatp = cos(rlatp)
312 radclatp = 0.5*rad*coslatp
313 radclatm = 0.5*rad*coslatm
314
315 DO i = 1, iim
316 xprp = xprimp025(i)
317 xprm = xprimm025(i)
318
319 ai14 = un4rad2*coslatp*yprp
320 ai23 = un4rad2*coslatm*yprm
321 aireij1_2d(i, j) = ai14*xprp
322 aireij2_2d(i, j) = ai23*xprp
323 aireij3_2d(i, j) = ai23*xprm
324 aireij4_2d(i, j) = ai14*xprm
325 cuij1(i, j) = radclatp*xprp
326 cuij2(i, j) = radclatm*xprp
327 cuij3(i, j) = radclatm*xprm
328 cuij4(i, j) = radclatp*xprm
329 cvij1(i, j) = 0.5*rad*yprp
330 cvij2(i, j) = 0.5*rad*yprm
331 cvij3(i, j) = cvij2(i, j)
332 cvij4(i, j) = cvij1(i, j)
333 END DO
334
335 END IF
336
337 ! periodicite
338
339 cvij1(iip1, j) = cvij1(1, j)
340 cvij2(iip1, j) = cvij2(1, j)
341 cvij3(iip1, j) = cvij3(1, j)
342 cvij4(iip1, j) = cvij4(1, j)
343 cuij1(iip1, j) = cuij1(1, j)
344 cuij2(iip1, j) = cuij2(1, j)
345 cuij3(iip1, j) = cuij3(1, j)
346 cuij4(iip1, j) = cuij4(1, j)
347 aireij1_2d(iip1, j) = aireij1_2d(1, j)
348 aireij2_2d(iip1, j) = aireij2_2d(1, j)
349 aireij3_2d(iip1, j) = aireij3_2d(1, j)
350 aireij4_2d(iip1, j) = aireij4_2d(1, j)
351
352 END DO
353
354 DO j = 1, jjp1
355 DO i = 1, iim
356 aire_2d(i, j) = aireij1_2d(i, j) + aireij2_2d(i, j) &
357 + aireij3_2d(i, j) + aireij4_2d(i, j)
358 alpha1_2d(i, j) = aireij1_2d(i, j)/aire_2d(i, j)
359 alpha2_2d(i, j) = aireij2_2d(i, j)/aire_2d(i, j)
360 alpha3_2d(i, j) = aireij3_2d(i, j)/aire_2d(i, j)
361 alpha4_2d(i, j) = aireij4_2d(i, j)/aire_2d(i, j)
362 alpha1p2_2d(i, j) = alpha1_2d(i, j) + alpha2_2d(i, j)
363 alpha1p4_2d(i, j) = alpha1_2d(i, j) + alpha4_2d(i, j)
364 alpha2p3_2d(i, j) = alpha2_2d(i, j) + alpha3_2d(i, j)
365 alpha3p4_2d(i, j) = alpha3_2d(i, j) + alpha4_2d(i, j)
366 END DO
367
368 aire_2d(iip1, j) = aire_2d(1, j)
369 alpha1_2d(iip1, j) = alpha1_2d(1, j)
370 alpha2_2d(iip1, j) = alpha2_2d(1, j)
371 alpha3_2d(iip1, j) = alpha3_2d(1, j)
372 alpha4_2d(iip1, j) = alpha4_2d(1, j)
373 alpha1p2_2d(iip1, j) = alpha1p2_2d(1, j)
374 alpha1p4_2d(iip1, j) = alpha1p4_2d(1, j)
375 alpha2p3_2d(iip1, j) = alpha2p3_2d(1, j)
376 alpha3p4_2d(iip1, j) = alpha3p4_2d(1, j)
377 END DO
378
379 DO j = 1, jjp1
380 DO i = 1, iim
381 aireu_2d(i, j) = aireij1_2d(i, j) + aireij2_2d(i, j) + &
382 aireij4_2d(i+1, j) + aireij3_2d(i+1, j)
383 unsaire_2d(i, j) = 1./aire_2d(i, j)
384 unsair_gam1_2d(i, j) = unsaire_2d(i, j)**(-gamdi_gdiv)
385 unsair_gam2_2d(i, j) = unsaire_2d(i, j)**(-gamdi_h)
386 airesurg_2d(i, j) = aire_2d(i, j)/g
387 END DO
388 aireu_2d(iip1, j) = aireu_2d(1, j)
389 unsaire_2d(iip1, j) = unsaire_2d(1, j)
390 unsair_gam1_2d(iip1, j) = unsair_gam1_2d(1, j)
391 unsair_gam2_2d(iip1, j) = unsair_gam2_2d(1, j)
392 airesurg_2d(iip1, j) = airesurg_2d(1, j)
393 END DO
394
395 DO j = 1, jjm
396
397 DO i = 1, iim
398 airev_2d(i, j) = aireij2_2d(i, j) + aireij3_2d(i, j) + &
399 aireij1_2d(i, j+1) + aireij4_2d(i, j+1)
400 END DO
401 DO i = 1, iim
402 airez = aireij2_2d(i, j) + aireij1_2d(i, j+1) + aireij3_2d(i+1, j) &
403 + aireij4_2d(i+1, j+1)
404 unsairez_2d(i, j) = 1./airez
405 unsairz_gam_2d(i, j) = unsairez_2d(i, j)**(-gamdi_grot)
406 fext_2d(i, j) = airez*sin(rlatv(j))*2.*omeg
407 END DO
408 airev_2d(iip1, j) = airev_2d(1, j)
409 unsairez_2d(iip1, j) = unsairez_2d(1, j)
410 fext_2d(iip1, j) = fext_2d(1, j)
411 unsairz_gam_2d(iip1, j) = unsairz_gam_2d(1, j)
412
413 END DO
414
415 ! Calcul des elongations cu_2d, cv_2d, cvu
416
417 DO j = 1, jjm
418 DO i = 1, iim
419 cv_2d(i, j) = 0.5 * &
420 (cvij2(i, j) + cvij3(i, j) + cvij1(i, j+1) + cvij4(i, j+1))
421 cvu(i, j) = 0.5*(cvij1(i, j)+cvij4(i, j)+cvij2(i, j)+cvij3(i, j))
422 cuv(i, j) = 0.5*(cuij2(i, j)+cuij3(i, j)+cuij1(i, j+1)+cuij4(i, j+1))
423 unscv2_2d(i, j) = 1./(cv_2d(i, j)*cv_2d(i, j))
424 END DO
425 DO i = 1, iim
426 cuvsurcv_2d(i, j) = airev_2d(i, j)*unscv2_2d(i, j)
427 cvsurcuv_2d(i, j) = 1./cuvsurcv_2d(i, j)
428 cuvscvgam1_2d(i, j) = cuvsurcv_2d(i, j)**(-gamdi_gdiv)
429 cuvscvgam2_2d(i, j) = cuvsurcv_2d(i, j)**(-gamdi_h)
430 cvscuvgam_2d(i, j) = cvsurcuv_2d(i, j)**(-gamdi_grot)
431 END DO
432 cv_2d(iip1, j) = cv_2d(1, j)
433 cvu(iip1, j) = cvu(1, j)
434 unscv2_2d(iip1, j) = unscv2_2d(1, j)
435 cuv(iip1, j) = cuv(1, j)
436 cuvsurcv_2d(iip1, j) = cuvsurcv_2d(1, j)
437 cvsurcuv_2d(iip1, j) = cvsurcuv_2d(1, j)
438 cuvscvgam1_2d(iip1, j) = cuvscvgam1_2d(1, j)
439 cuvscvgam2_2d(iip1, j) = cuvscvgam2_2d(1, j)
440 cvscuvgam_2d(iip1, j) = cvscuvgam_2d(1, j)
441 END DO
442
443 DO j = 2, jjm
444 DO i = 1, iim
445 cu_2d(i, j) = 0.5 * (cuij1(i, j) + cuij4(i+1, j) + cuij2(i, j) &
446 + cuij3(i+1, j))
447 unscu2_2d(i, j) = 1./(cu_2d(i, j)*cu_2d(i, j))
448 cvusurcu_2d(i, j) = aireu_2d(i, j)*unscu2_2d(i, j)
449 cusurcvu_2d(i, j) = 1./cvusurcu_2d(i, j)
450 cvuscugam1_2d(i, j) = cvusurcu_2d(i, j)**(-gamdi_gdiv)
451 cvuscugam2_2d(i, j) = cvusurcu_2d(i, j)**(-gamdi_h)
452 cuscvugam_2d(i, j) = cusurcvu_2d(i, j)**(-gamdi_grot)
453 END DO
454 cu_2d(iip1, j) = cu_2d(1, j)
455 unscu2_2d(iip1, j) = unscu2_2d(1, j)
456 cvusurcu_2d(iip1, j) = cvusurcu_2d(1, j)
457 cusurcvu_2d(iip1, j) = cusurcvu_2d(1, j)
458 cvuscugam1_2d(iip1, j) = cvuscugam1_2d(1, j)
459 cvuscugam2_2d(iip1, j) = cvuscugam2_2d(1, j)
460 cuscvugam_2d(iip1, j) = cuscvugam_2d(1, j)
461 END DO
462
463 ! calcul aux poles
464
465 DO i = 1, iip1
466 cu_2d(i, 1) = 0.
467 unscu2_2d(i, 1) = 0.
468 cvu(i, 1) = 0.
469
470 cu_2d(i, jjp1) = 0.
471 unscu2_2d(i, jjp1) = 0.
472 cvu(i, jjp1) = 0.
473 END DO
474
475 DO j = 1, jjm
476 DO i = 1, iim
477 airvscu2_2d(i, j) = airev_2d(i, j)/(cuv(i, j)*cuv(i, j))
478 aivscu2gam_2d(i, j) = airvscu2_2d(i, j)**(-gamdi_grot)
479 END DO
480 airvscu2_2d(iip1, j) = airvscu2_2d(1, j)
481 aivscu2gam_2d(iip1, j) = aivscu2gam_2d(1, j)
482 END DO
483
484 DO j = 2, jjm
485 DO i = 1, iim
486 airuscv2_2d(i, j) = aireu_2d(i, j)/(cvu(i, j)*cvu(i, j))
487 aiuscv2gam_2d(i, j) = airuscv2_2d(i, j)**(-gamdi_grot)
488 END DO
489 airuscv2_2d(iip1, j) = airuscv2_2d(1, j)
490 aiuscv2gam_2d(iip1, j) = aiuscv2gam_2d(1, j)
491 END DO
492
493 ! calcul des aires aux poles :
494
495 apoln = sum(aire_2d(:iim, 1))
496 apols = sum(aire_2d(:iim, jjp1))
497 unsapolnga1 = 1./(apoln**(-gamdi_gdiv))
498 unsapolsga1 = 1./(apols**(-gamdi_gdiv))
499 unsapolnga2 = 1./(apoln**(-gamdi_h))
500 unsapolsga2 = 1./(apols**(-gamdi_h))
501
502 ! changement F. Hourdin calcul conservatif pour fext_2d
503 ! constang_2d contient le produit a * cos (latitude) * omega
504
505 DO i = 1, iim
506 constang_2d(i, 1) = 0.
507 END DO
508 DO j = 1, jjm - 1
509 DO i = 1, iim
510 constang_2d(i, j+1) = rad*omeg*cu_2d(i, j+1)*cos(rlatu(j+1))
511 END DO
512 END DO
513 DO i = 1, iim
514 constang_2d(i, jjp1) = 0.
515 END DO
516
517 ! periodicite en longitude
518
519 DO j = 1, jjm
520 fext_2d(iip1, j) = fext_2d(1, j)
521 END DO
522 DO j = 1, jjp1
523 constang_2d(iip1, j) = constang_2d(1, j)
524 END DO
525
526 ! fin du changement
527
528 print *, ' Coordonnees de la grille '
529 WRITE (6, 995)
530
531 print *, ' LONGITUDES aux pts. V (degres) '
532 WRITE (6, 995)
533 DO i = 1, iip1
534 rlonvv(i) = rlonv(i)*180./pi
535 END DO
536 WRITE (6, 400) rlonvv
537
538 WRITE (6, 995)
539 print *, ' LATITUDES aux pts. V (degres) '
540 WRITE (6, 995)
541 DO i = 1, jjm
542 rlatuu(i) = rlatv(i)*180./pi
543 END DO
544 WRITE (6, 400) (rlatuu(i), i=1, jjm)
545
546 DO i = 1, iip1
547 rlonvv(i) = rlonu(i)*180./pi
548 END DO
549 WRITE (6, 995)
550 print *, ' LONGITUDES aux pts. U (degres) '
551 WRITE (6, 995)
552 WRITE (6, 400) rlonvv
553 WRITE (6, 995)
554
555 print *, ' LATITUDES aux pts. U (degres) '
556 WRITE (6, 995)
557 DO i = 1, jjp1
558 rlatuu(i) = rlatu(i)*180./pi
559 END DO
560 WRITE (6, 400) (rlatuu(i), i=1, jjp1)
561 WRITE (6, 995)
562
563 400 FORMAT (1X, 8F8.2)
564 990 FORMAT (//)
565 995 FORMAT (/)
566
567 END SUBROUTINE inigeom
568
569 end module inigeom_m

  ViewVC Help
Powered by ViewVC 1.1.21