1 | ! |
---|
2 | ! $Id: $ |
---|
3 | ! |
---|
4 | |
---|
5 | c#define DEBUG_IO |
---|
6 | |
---|
7 | SUBROUTINE caldyn_p |
---|
8 | $ (itau,ucov,vcov,teta,ps,masse,pk,pkf,tsurpk,phis , |
---|
9 | $ phi,conser,du,dv,dteta,dp,w,pbaru,pbarv,time ) |
---|
10 | USE parallel_lmdz |
---|
11 | USE Write_Field_p |
---|
12 | |
---|
13 | IMPLICIT NONE |
---|
14 | |
---|
15 | !======================================================================= |
---|
16 | ! |
---|
17 | ! Auteur : P. Le Van |
---|
18 | ! |
---|
19 | ! Objet: |
---|
20 | ! ------ |
---|
21 | ! |
---|
22 | ! Calcul des tendances dynamiques. |
---|
23 | ! |
---|
24 | ! Modif 04/93 F.Forget |
---|
25 | !======================================================================= |
---|
26 | |
---|
27 | !----------------------------------------------------------------------- |
---|
28 | ! 0. Declarations: |
---|
29 | ! ---------------- |
---|
30 | |
---|
31 | !----------------------------------------------------------------------- |
---|
32 | ! INCLUDE 'dimensions.h' |
---|
33 | ! |
---|
34 | ! dimensions.h contient les dimensions du modele |
---|
35 | ! ndm est tel que iim=2**ndm |
---|
36 | !----------------------------------------------------------------------- |
---|
37 | |
---|
38 | INTEGER iim,jjm,llm,ndm |
---|
39 | |
---|
40 | PARAMETER (iim= 128,jjm=96,llm=64,ndm=1) |
---|
41 | |
---|
42 | !----------------------------------------------------------------------- |
---|
43 | ! |
---|
44 | ! $Header$ |
---|
45 | ! |
---|
46 | ! |
---|
47 | ! ATTENTION!!!!: ce fichier include est compatible format fixe/format libre |
---|
48 | ! veillez n'utiliser que des ! pour les commentaires |
---|
49 | ! et bien positionner les & des lignes de continuation |
---|
50 | ! (les placer en colonne 6 et en colonne 73) |
---|
51 | ! |
---|
52 | ! |
---|
53 | !----------------------------------------------------------------------- |
---|
54 | ! INCLUDE 'paramet.h' |
---|
55 | |
---|
56 | INTEGER iip1,iip2,iip3,jjp1,llmp1,llmp2,llmm1 |
---|
57 | INTEGER kftd,ip1jm,ip1jmp1,ip1jmi1,ijp1llm |
---|
58 | INTEGER ijmllm,mvar |
---|
59 | INTEGER jcfil,jcfllm |
---|
60 | |
---|
61 | PARAMETER( iip1= iim+1,iip2=iim+2,iip3=iim+3 & |
---|
62 | & ,jjp1=jjm+1-1/jjm) |
---|
63 | PARAMETER( llmp1 = llm+1, llmp2 = llm+2, llmm1 = llm-1 ) |
---|
64 | PARAMETER( kftd = iim/2 -ndm ) |
---|
65 | PARAMETER( ip1jm = iip1*jjm, ip1jmp1= iip1*jjp1 ) |
---|
66 | PARAMETER( ip1jmi1= ip1jm - iip1 ) |
---|
67 | PARAMETER( ijp1llm= ip1jmp1 * llm, ijmllm= ip1jm * llm ) |
---|
68 | PARAMETER( mvar= ip1jmp1*( 2*llm+1) + ijmllm ) |
---|
69 | PARAMETER( jcfil=jjm/2+5, jcfllm=jcfil*llm ) |
---|
70 | |
---|
71 | !----------------------------------------------------------------------- |
---|
72 | ! |
---|
73 | ! $Id: comconst.h 1437 2010-09-30 08:29:10Z emillour $ |
---|
74 | ! |
---|
75 | !----------------------------------------------------------------------- |
---|
76 | ! INCLUDE comconst.h |
---|
77 | |
---|
78 | COMMON/comconsti/im,jm,lllm,imp1,jmp1,lllmm1,lllmp1,lcl, & |
---|
79 | & iflag_top_bound,mode_top_bound |
---|
80 | COMMON/comconstr/dtvr,daysec, & |
---|
81 | & pi,dtphys,dtdiss,rad,r,kappa,cotot,unsim,g,omeg & |
---|
82 | & ,dissip_fac_mid,dissip_fac_up,dissip_deltaz,dissip_hdelta & |
---|
83 | & ,dissip_pupstart ,tau_top_bound, & |
---|
84 | & daylen,molmass, ihf |
---|
85 | COMMON/cpdetvenus/cpp,nu_venus,t0_venus |
---|
86 | |
---|
87 | INTEGER im,jm,lllm,imp1,jmp1,lllmm1,lllmp1,lcl |
---|
88 | REAL dtvr ! dynamical time step (in s) |
---|
89 | REAL daysec !length (in s) of a standard day |
---|
90 | REAL pi ! something like 3.14159.... |
---|
91 | REAL dtphys ! (s) time step for the physics |
---|
92 | REAL dtdiss ! (s) time step for the dissipation |
---|
93 | REAL rad ! (m) radius of the planet |
---|
94 | REAL r ! Reduced Gas constant r=R/mu |
---|
95 | ! with R=8.31.. J.K-1.mol-1, mu: mol mass of atmosphere (kg/mol) |
---|
96 | REAL cpp ! Cp |
---|
97 | REAL kappa ! kappa=R/Cp |
---|
98 | REAL cotot |
---|
99 | REAL unsim ! = 1./iim |
---|
100 | REAL g ! (m/s2) gravity |
---|
101 | REAL omeg ! (rad/s) rotation rate of the planet |
---|
102 | ! Dissipation factors, for Earth model: |
---|
103 | REAL dissip_factz,dissip_zref !dissip_deltaz |
---|
104 | ! Dissipation factors, for other planets: |
---|
105 | REAL dissip_fac_mid,dissip_fac_up,dissip_deltaz,dissip_hdelta |
---|
106 | REAL dissip_pupstart |
---|
107 | INTEGER iflag_top_bound,mode_top_bound |
---|
108 | REAL tau_top_bound |
---|
109 | REAL daylen ! length of solar day, in 'standard' day length |
---|
110 | REAL molmass ! (g/mol) molar mass of the atmosphere |
---|
111 | |
---|
112 | REAL nu_venus,t0_venus ! coeffs needed for Cp(T), Venus atmosphere |
---|
113 | REAL ihf ! (W/m2) intrinsic heat flux for giant planets |
---|
114 | |
---|
115 | |
---|
116 | !----------------------------------------------------------------------- |
---|
117 | ! |
---|
118 | ! $Id: comvert.h 1654 2012-09-24 15:07:18Z aslmd $ |
---|
119 | ! |
---|
120 | !----------------------------------------------------------------------- |
---|
121 | ! INCLUDE 'comvert.h' |
---|
122 | |
---|
123 | COMMON/comvertr/ap(llm+1),bp(llm+1),presnivs(llm),dpres(llm), & |
---|
124 | & pa,preff,nivsigs(llm),nivsig(llm+1), & |
---|
125 | & aps(llm),bps(llm),scaleheight,pseudoalt(llm) |
---|
126 | |
---|
127 | common/comverti/disvert_type, pressure_exner |
---|
128 | |
---|
129 | real ap ! hybrid pressure contribution at interlayers |
---|
130 | real bp ! hybrid sigma contribution at interlayer |
---|
131 | real presnivs ! (reference) pressure at mid-layers |
---|
132 | real dpres |
---|
133 | real pa ! reference pressure (Pa) at which hybrid coordinates |
---|
134 | ! become purely pressure |
---|
135 | real preff ! reference surface pressure (Pa) |
---|
136 | real nivsigs |
---|
137 | real nivsig |
---|
138 | real aps ! hybrid pressure contribution at mid-layers |
---|
139 | real bps ! hybrid sigma contribution at mid-layers |
---|
140 | real scaleheight ! atmospheric (reference) scale height (km) |
---|
141 | real pseudoalt ! pseudo-altitude of model levels (km), based on presnivs(), |
---|
142 | ! preff and scaleheight |
---|
143 | |
---|
144 | integer disvert_type ! type of vertical discretization: |
---|
145 | ! 1: Earth (default for planet_type==earth), |
---|
146 | ! automatic generation |
---|
147 | ! 2: Planets (default for planet_type!=earth), |
---|
148 | ! using 'z2sig.def' (or 'esasig.def) file |
---|
149 | |
---|
150 | logical pressure_exner |
---|
151 | ! compute pressure inside layers using Exner function, else use mean |
---|
152 | ! of pressure values at interfaces |
---|
153 | |
---|
154 | !----------------------------------------------------------------------- |
---|
155 | ! |
---|
156 | ! $Header$ |
---|
157 | ! |
---|
158 | !CDK comgeom |
---|
159 | COMMON/comgeom/ & |
---|
160 | & cu(ip1jmp1),cv(ip1jm),unscu2(ip1jmp1),unscv2(ip1jm), & |
---|
161 | & aire(ip1jmp1),airesurg(ip1jmp1),aireu(ip1jmp1), & |
---|
162 | & airev(ip1jm),unsaire(ip1jmp1),apoln,apols, & |
---|
163 | & unsairez(ip1jm),airuscv2(ip1jm),airvscu2(ip1jm), & |
---|
164 | & aireij1(ip1jmp1),aireij2(ip1jmp1),aireij3(ip1jmp1), & |
---|
165 | & aireij4(ip1jmp1),alpha1(ip1jmp1),alpha2(ip1jmp1), & |
---|
166 | & alpha3(ip1jmp1),alpha4(ip1jmp1),alpha1p2(ip1jmp1), & |
---|
167 | & alpha1p4(ip1jmp1),alpha2p3(ip1jmp1),alpha3p4(ip1jmp1), & |
---|
168 | & fext(ip1jm),constang(ip1jmp1),rlatu(jjp1),rlatv(jjm), & |
---|
169 | & rlonu(iip1),rlonv(iip1),cuvsurcv(ip1jm),cvsurcuv(ip1jm), & |
---|
170 | & cvusurcu(ip1jmp1),cusurcvu(ip1jmp1),cuvscvgam1(ip1jm), & |
---|
171 | & cuvscvgam2(ip1jm),cvuscugam1(ip1jmp1), & |
---|
172 | & cvuscugam2(ip1jmp1),cvscuvgam(ip1jm),cuscvugam(ip1jmp1), & |
---|
173 | & unsapolnga1,unsapolnga2,unsapolsga1,unsapolsga2, & |
---|
174 | & unsair_gam1(ip1jmp1),unsair_gam2(ip1jmp1),unsairz_gam(ip1jm), & |
---|
175 | & aivscu2gam(ip1jm),aiuscv2gam(ip1jm),xprimu(iip1),xprimv(iip1) |
---|
176 | |
---|
177 | ! |
---|
178 | REAL & |
---|
179 | & cu,cv,unscu2,unscv2,aire,airesurg,aireu,airev,unsaire,apoln ,& |
---|
180 | & apols,unsairez,airuscv2,airvscu2,aireij1,aireij2,aireij3,aireij4,& |
---|
181 | & alpha1,alpha2,alpha3,alpha4,alpha1p2,alpha1p4,alpha2p3,alpha3p4 ,& |
---|
182 | & fext,constang,rlatu,rlatv,rlonu,rlonv,cuvscvgam1,cuvscvgam2 ,& |
---|
183 | & cvuscugam1,cvuscugam2,cvscuvgam,cuscvugam,unsapolnga1,unsapolnga2& |
---|
184 | & ,unsapolsga1,unsapolsga2,unsair_gam1,unsair_gam2,unsairz_gam ,& |
---|
185 | & aivscu2gam ,aiuscv2gam,cuvsurcv,cvsurcuv,cvusurcu,cusurcvu,xprimu& |
---|
186 | & , xprimv |
---|
187 | ! |
---|
188 | |
---|
189 | ! Arguments: |
---|
190 | ! ---------- |
---|
191 | |
---|
192 | LOGICAL,INTENT(IN) :: conser ! triggers printing some diagnostics |
---|
193 | INTEGER,INTENT(IN) :: itau ! time step index |
---|
194 | REAL,INTENT(IN) :: vcov(ip1jm,llm) ! covariant meridional wind |
---|
195 | REAL,INTENT(IN) :: ucov(ip1jmp1,llm) ! covariant zonal wind |
---|
196 | REAL,INTENT(IN) :: teta(ip1jmp1,llm) ! potential temperature |
---|
197 | REAL,INTENT(IN) :: ps(ip1jmp1) ! surface pressure |
---|
198 | REAL,INTENT(IN) :: phis(ip1jmp1) ! geopotential at the surface |
---|
199 | REAL,INTENT(IN) :: pk(ip1jmp1,llm) ! Exner at mid-layer |
---|
200 | REAL,INTENT(IN) :: pkf(ip1jmp1,llm) ! filtered Exner |
---|
201 | REAL,INTENT(IN) :: tsurpk(ip1jmp1,llm) ! cpp * temperature / pk |
---|
202 | REAL,INTENT(IN) :: phi(ip1jmp1,llm) ! geopotential |
---|
203 | REAL,INTENT(OUT) :: masse(ip1jmp1,llm) ! air mass |
---|
204 | REAL,INTENT(OUT) :: dv(ip1jm,llm) ! tendency on vcov |
---|
205 | REAL,INTENT(OUT) :: du(ip1jmp1,llm) ! tendency on ucov |
---|
206 | REAL,INTENT(OUT) :: dteta(ip1jmp1,llm) ! tenddency on teta |
---|
207 | REAL,INTENT(OUT) :: dp(ip1jmp1) ! tendency on ps |
---|
208 | REAL,INTENT(OUT) :: w(ip1jmp1,llm) ! vertical velocity |
---|
209 | REAL,INTENT(OUT) :: pbaru(ip1jmp1,llm) ! mass flux in the zonal direction |
---|
210 | REAL,INTENT(OUT) :: pbarv(ip1jm,llm) ! mass flux in the meridional direction |
---|
211 | REAL,INTENT(IN) :: time ! current time |
---|
212 | |
---|
213 | ! Local: |
---|
214 | ! ------ |
---|
215 | |
---|
216 | REAL,SAVE :: vcont(ip1jm,llm),ucont(ip1jmp1,llm) |
---|
217 | REAL,SAVE :: ang(ip1jmp1,llm) |
---|
218 | REAL,SAVE :: p(ip1jmp1,llmp1) |
---|
219 | REAL,SAVE :: massebx(ip1jmp1,llm),masseby(ip1jm,llm) |
---|
220 | REAL,SAVE :: psexbarxy(ip1jm) |
---|
221 | REAL,SAVE :: vorpot(ip1jm,llm) |
---|
222 | REAL,SAVE :: ecin(ip1jmp1,llm) |
---|
223 | REAL,SAVE :: bern(ip1jmp1,llm) |
---|
224 | REAL,SAVE :: massebxy(ip1jm,llm) |
---|
225 | REAL,SAVE :: convm(ip1jmp1,llm) |
---|
226 | ! REAL,SAVE :: temp(ip1jmp1,llm) |
---|
227 | INTEGER ij,l,ijb,ije,ierr |
---|
228 | |
---|
229 | !----------------------------------------------------------------------- |
---|
230 | ! Compute dynamical tendencies: |
---|
231 | !-------------------------------- |
---|
232 | |
---|
233 | ! compute contravariant winds ucont() and vcont |
---|
234 | CALL covcont_p ( llm , ucov , vcov , ucont, vcont ) |
---|
235 | ! compute pressure p() |
---|
236 | CALL pression_p ( ip1jmp1, ap , bp , ps , p ) |
---|
237 | cym CALL psextbar ( ps , psexbarxy ) |
---|
238 | c$OMP BARRIER |
---|
239 | ! compute mass in each atmospheric mesh: masse() |
---|
240 | CALL massdair_p ( p , masse ) |
---|
241 | ! compute X and Y-averages of mass, massebx() and masseby() |
---|
242 | CALL massbar_p ( masse, massebx , masseby ) |
---|
243 | ! compute XY-average of mass, massebxy() |
---|
244 | call massbarxy_p( masse, massebxy ) |
---|
245 | ! compute mass fluxes pbaru() and pbarv() |
---|
246 | CALL flumass_p ( massebx, masseby , vcont, ucont ,pbaru, pbarv ) |
---|
247 | ! compute dteta() , horizontal converging flux of theta |
---|
248 | CALL dteta1_p ( teta , pbaru , pbarv, dteta ) |
---|
249 | ! compute convm(), horizontal converging flux of mass |
---|
250 | CALL convmas1_p ( pbaru, pbarv , convm ) |
---|
251 | c$OMP BARRIER |
---|
252 | CALL convmas2_p ( convm ) |
---|
253 | c$OMP BARRIER |
---|
254 | |
---|
255 | c$OMP BARRIER |
---|
256 | c$OMP MASTER |
---|
257 | ijb=ij_begin |
---|
258 | ije=ij_end |
---|
259 | ! compute pressure variation due to mass convergence |
---|
260 | DO ij =ijb, ije |
---|
261 | dp( ij ) = convm( ij,1 ) / airesurg( ij ) |
---|
262 | ENDDO |
---|
263 | c$OMP END MASTER |
---|
264 | c$OMP BARRIER |
---|
265 | c$OMP FLUSH |
---|
266 | |
---|
267 | ! compute vertical velocity w() |
---|
268 | CALL vitvert_p ( convm , w ) |
---|
269 | ! compute potential vorticity vorpot() |
---|
270 | CALL tourpot_p ( vcov , ucov , massebxy , vorpot ) |
---|
271 | ! compute rotation induced du() and dv() |
---|
272 | CALL dudv1_p ( vorpot , pbaru , pbarv , du , dv ) |
---|
273 | |
---|
274 | |
---|
275 | ! compute kinetic energy ecin() |
---|
276 | CALL enercin_p ( vcov , ucov , vcont , ucont , ecin ) |
---|
277 | ! compute Bernouilli function bern() |
---|
278 | CALL bernoui_p ( ip1jmp1, llm , phi , ecin , bern ) |
---|
279 | ! compute and add du() and dv() contributions from Bernouilli and pressure |
---|
280 | CALL dudv2_p ( tsurpk , pkf , bern , du , dv ) |
---|
281 | |
---|
282 | |
---|
283 | ijb=ij_begin-iip1 |
---|
284 | ije=ij_end+iip1 |
---|
285 | |
---|
286 | if (pole_nord) ijb=ij_begin |
---|
287 | if (pole_sud) ije=ij_end |
---|
288 | |
---|
289 | c$OMP DO SCHEDULE(STATIC,OMP_CHUNK) |
---|
290 | DO l=1,llm |
---|
291 | DO ij=ijb,ije |
---|
292 | ang(ij,l) = ucov(ij,l) + constang(ij) |
---|
293 | ENDDO |
---|
294 | ENDDO |
---|
295 | c$OMP END DO |
---|
296 | |
---|
297 | ! compute vertical advection contributions to du(), dv() and dteta() |
---|
298 | CALL advect_new_p(ang,vcov,teta,w,massebx,masseby,du,dv,dteta) |
---|
299 | |
---|
300 | C WARNING probleme de peridocite de dv sur les PC/1. Pb d'arrondi |
---|
301 | C probablement. Observe sur le code compile avec pgf90 3.0-1 |
---|
302 | ijb=ij_begin |
---|
303 | ije=ij_end |
---|
304 | if (pole_sud) ije=ij_end-iip1 |
---|
305 | |
---|
306 | c$OMP DO SCHEDULE(STATIC,OMP_CHUNK) |
---|
307 | DO l = 1, llm |
---|
308 | DO ij = ijb, ije, iip1 |
---|
309 | IF( dv(ij,l).NE.dv(ij+iim,l) ) THEN |
---|
310 | c PRINT *,'!!!ATTENTION!!! probleme de periodicite sur vcov', |
---|
311 | c , ' dans caldyn' |
---|
312 | c PRINT *,' l, ij = ', l, ij, ij+iim,dv(ij+iim,l),dv(ij,l) |
---|
313 | dv(ij+iim,l) = dv(ij,l) |
---|
314 | endif |
---|
315 | enddo |
---|
316 | enddo |
---|
317 | c$OMP END DO NOWAIT |
---|
318 | c----------------------------------------------------------------------- |
---|
319 | c Sorties eventuelles des variables de controle: |
---|
320 | c ---------------------------------------------- |
---|
321 | |
---|
322 | IF( conser ) THEN |
---|
323 | c ym ---> exige communication collective ( aussi dans advect) |
---|
324 | CALL sortvarc |
---|
325 | $ (itau,ucov,tsurpk,ps,masse,pk,phis,vorpot,phi,bern,dp,time,vcov) |
---|
326 | |
---|
327 | ENDIF |
---|
328 | |
---|
329 | END |
---|