1 | subroutine groupe_p(pext,pbaru,pbarv,pbarum,pbarvm,wm) |
---|
2 | USE parallel_lmdz |
---|
3 | implicit none |
---|
4 | |
---|
5 | c sous-programme servant a fitlrer les champs de flux de masse aux |
---|
6 | c poles en "regroupant" les mailles 2 par 2 puis 4 par 4 etc. au fur |
---|
7 | c et a mesure qu'on se rapproche du pole. |
---|
8 | c |
---|
9 | c en entree: pext, pbaru et pbarv |
---|
10 | c |
---|
11 | c en sortie: pbarum,pbarvm et wm. |
---|
12 | c |
---|
13 | c remarque, le wm est recalcule a partir des pbaru pbarv et on n'a donc |
---|
14 | c pas besoin de w en entree. |
---|
15 | |
---|
16 | !----------------------------------------------------------------------- |
---|
17 | ! INCLUDE 'dimensions.h' |
---|
18 | ! |
---|
19 | ! dimensions.h contient les dimensions du modele |
---|
20 | ! ndm est tel que iim=2**ndm |
---|
21 | !----------------------------------------------------------------------- |
---|
22 | |
---|
23 | INTEGER iim,jjm,llm,ndm |
---|
24 | |
---|
25 | PARAMETER (iim= 128,jjm=96,llm=64,ndm=1) |
---|
26 | |
---|
27 | !----------------------------------------------------------------------- |
---|
28 | ! |
---|
29 | ! $Header$ |
---|
30 | ! |
---|
31 | ! |
---|
32 | ! ATTENTION!!!!: ce fichier include est compatible format fixe/format libre |
---|
33 | ! veillez n'utiliser que des ! pour les commentaires |
---|
34 | ! et bien positionner les & des lignes de continuation |
---|
35 | ! (les placer en colonne 6 et en colonne 73) |
---|
36 | ! |
---|
37 | ! |
---|
38 | !----------------------------------------------------------------------- |
---|
39 | ! INCLUDE 'paramet.h' |
---|
40 | |
---|
41 | INTEGER iip1,iip2,iip3,jjp1,llmp1,llmp2,llmm1 |
---|
42 | INTEGER kftd,ip1jm,ip1jmp1,ip1jmi1,ijp1llm |
---|
43 | INTEGER ijmllm,mvar |
---|
44 | INTEGER jcfil,jcfllm |
---|
45 | |
---|
46 | PARAMETER( iip1= iim+1,iip2=iim+2,iip3=iim+3 & |
---|
47 | & ,jjp1=jjm+1-1/jjm) |
---|
48 | PARAMETER( llmp1 = llm+1, llmp2 = llm+2, llmm1 = llm-1 ) |
---|
49 | PARAMETER( kftd = iim/2 -ndm ) |
---|
50 | PARAMETER( ip1jm = iip1*jjm, ip1jmp1= iip1*jjp1 ) |
---|
51 | PARAMETER( ip1jmi1= ip1jm - iip1 ) |
---|
52 | PARAMETER( ijp1llm= ip1jmp1 * llm, ijmllm= ip1jm * llm ) |
---|
53 | PARAMETER( mvar= ip1jmp1*( 2*llm+1) + ijmllm ) |
---|
54 | PARAMETER( jcfil=jjm/2+5, jcfllm=jcfil*llm ) |
---|
55 | |
---|
56 | !----------------------------------------------------------------------- |
---|
57 | ! |
---|
58 | ! $Id: comconst.h 1437 2010-09-30 08:29:10Z emillour $ |
---|
59 | ! |
---|
60 | !----------------------------------------------------------------------- |
---|
61 | ! INCLUDE comconst.h |
---|
62 | |
---|
63 | COMMON/comconsti/im,jm,lllm,imp1,jmp1,lllmm1,lllmp1,lcl, & |
---|
64 | & iflag_top_bound,mode_top_bound |
---|
65 | COMMON/comconstr/dtvr,daysec, & |
---|
66 | & pi,dtphys,dtdiss,rad,r,kappa,cotot,unsim,g,omeg & |
---|
67 | & ,dissip_fac_mid,dissip_fac_up,dissip_deltaz,dissip_hdelta & |
---|
68 | & ,dissip_pupstart ,tau_top_bound, & |
---|
69 | & daylen,molmass, ihf |
---|
70 | COMMON/cpdetvenus/cpp,nu_venus,t0_venus |
---|
71 | |
---|
72 | INTEGER im,jm,lllm,imp1,jmp1,lllmm1,lllmp1,lcl |
---|
73 | REAL dtvr ! dynamical time step (in s) |
---|
74 | REAL daysec !length (in s) of a standard day |
---|
75 | REAL pi ! something like 3.14159.... |
---|
76 | REAL dtphys ! (s) time step for the physics |
---|
77 | REAL dtdiss ! (s) time step for the dissipation |
---|
78 | REAL rad ! (m) radius of the planet |
---|
79 | REAL r ! Reduced Gas constant r=R/mu |
---|
80 | ! with R=8.31.. J.K-1.mol-1, mu: mol mass of atmosphere (kg/mol) |
---|
81 | REAL cpp ! Cp |
---|
82 | REAL kappa ! kappa=R/Cp |
---|
83 | REAL cotot |
---|
84 | REAL unsim ! = 1./iim |
---|
85 | REAL g ! (m/s2) gravity |
---|
86 | REAL omeg ! (rad/s) rotation rate of the planet |
---|
87 | ! Dissipation factors, for Earth model: |
---|
88 | REAL dissip_factz,dissip_zref !dissip_deltaz |
---|
89 | ! Dissipation factors, for other planets: |
---|
90 | REAL dissip_fac_mid,dissip_fac_up,dissip_deltaz,dissip_hdelta |
---|
91 | REAL dissip_pupstart |
---|
92 | INTEGER iflag_top_bound,mode_top_bound |
---|
93 | REAL tau_top_bound |
---|
94 | REAL daylen ! length of solar day, in 'standard' day length |
---|
95 | REAL molmass ! (g/mol) molar mass of the atmosphere |
---|
96 | |
---|
97 | REAL nu_venus,t0_venus ! coeffs needed for Cp(T), Venus atmosphere |
---|
98 | REAL ihf ! (W/m2) intrinsic heat flux for giant planets |
---|
99 | |
---|
100 | |
---|
101 | !----------------------------------------------------------------------- |
---|
102 | ! |
---|
103 | ! $Header$ |
---|
104 | ! |
---|
105 | !CDK comgeom2 |
---|
106 | COMMON/comgeom/ & |
---|
107 | & cu(iip1,jjp1),cv(iip1,jjm),unscu2(iip1,jjp1),unscv2(iip1,jjm) , & |
---|
108 | & aire(iip1,jjp1),airesurg(iip1,jjp1),aireu(iip1,jjp1) , & |
---|
109 | & airev(iip1,jjm),unsaire(iip1,jjp1),apoln,apols , & |
---|
110 | & unsairez(iip1,jjm),airuscv2(iip1,jjm),airvscu2(iip1,jjm) , & |
---|
111 | & aireij1(iip1,jjp1),aireij2(iip1,jjp1),aireij3(iip1,jjp1) , & |
---|
112 | & aireij4(iip1,jjp1),alpha1(iip1,jjp1),alpha2(iip1,jjp1) , & |
---|
113 | & alpha3(iip1,jjp1),alpha4(iip1,jjp1),alpha1p2(iip1,jjp1) , & |
---|
114 | & alpha1p4(iip1,jjp1),alpha2p3(iip1,jjp1),alpha3p4(iip1,jjp1) , & |
---|
115 | & fext(iip1,jjm),constang(iip1,jjp1), rlatu(jjp1),rlatv(jjm), & |
---|
116 | & rlonu(iip1),rlonv(iip1),cuvsurcv(iip1,jjm),cvsurcuv(iip1,jjm) , & |
---|
117 | & cvusurcu(iip1,jjp1),cusurcvu(iip1,jjp1) , & |
---|
118 | & cuvscvgam1(iip1,jjm),cuvscvgam2(iip1,jjm),cvuscugam1(iip1,jjp1), & |
---|
119 | & cvuscugam2(iip1,jjp1),cvscuvgam(iip1,jjm),cuscvugam(iip1,jjp1) , & |
---|
120 | & unsapolnga1,unsapolnga2,unsapolsga1,unsapolsga2 , & |
---|
121 | & unsair_gam1(iip1,jjp1),unsair_gam2(iip1,jjp1) , & |
---|
122 | & unsairz_gam(iip1,jjm),aivscu2gam(iip1,jjm),aiuscv2gam(iip1,jjm) & |
---|
123 | & , xprimu(iip1),xprimv(iip1) |
---|
124 | |
---|
125 | |
---|
126 | REAL & |
---|
127 | & cu,cv,unscu2,unscv2,aire,airesurg,aireu,airev,apoln,apols,unsaire & |
---|
128 | & ,unsairez,airuscv2,airvscu2,aireij1,aireij2,aireij3,aireij4 , & |
---|
129 | & alpha1,alpha2,alpha3,alpha4,alpha1p2,alpha1p4,alpha2p3,alpha3p4 , & |
---|
130 | & fext,constang,rlatu,rlatv,rlonu,rlonv,cuvscvgam1,cuvscvgam2 , & |
---|
131 | & cvuscugam1,cvuscugam2,cvscuvgam,cuscvugam,unsapolnga1 , & |
---|
132 | & unsapolnga2,unsapolsga1,unsapolsga2,unsair_gam1,unsair_gam2 , & |
---|
133 | & unsairz_gam,aivscu2gam,aiuscv2gam,cuvsurcv,cvsurcuv,cvusurcu , & |
---|
134 | & cusurcvu,xprimu,xprimv |
---|
135 | ! |
---|
136 | ! $Id: comvert.h 1654 2012-09-24 15:07:18Z aslmd $ |
---|
137 | ! |
---|
138 | !----------------------------------------------------------------------- |
---|
139 | ! INCLUDE 'comvert.h' |
---|
140 | |
---|
141 | COMMON/comvertr/ap(llm+1),bp(llm+1),presnivs(llm),dpres(llm), & |
---|
142 | & pa,preff,nivsigs(llm),nivsig(llm+1), & |
---|
143 | & aps(llm),bps(llm),scaleheight,pseudoalt(llm) |
---|
144 | |
---|
145 | common/comverti/disvert_type, pressure_exner |
---|
146 | |
---|
147 | real ap ! hybrid pressure contribution at interlayers |
---|
148 | real bp ! hybrid sigma contribution at interlayer |
---|
149 | real presnivs ! (reference) pressure at mid-layers |
---|
150 | real dpres |
---|
151 | real pa ! reference pressure (Pa) at which hybrid coordinates |
---|
152 | ! become purely pressure |
---|
153 | real preff ! reference surface pressure (Pa) |
---|
154 | real nivsigs |
---|
155 | real nivsig |
---|
156 | real aps ! hybrid pressure contribution at mid-layers |
---|
157 | real bps ! hybrid sigma contribution at mid-layers |
---|
158 | real scaleheight ! atmospheric (reference) scale height (km) |
---|
159 | real pseudoalt ! pseudo-altitude of model levels (km), based on presnivs(), |
---|
160 | ! preff and scaleheight |
---|
161 | |
---|
162 | integer disvert_type ! type of vertical discretization: |
---|
163 | ! 1: Earth (default for planet_type==earth), |
---|
164 | ! automatic generation |
---|
165 | ! 2: Planets (default for planet_type!=earth), |
---|
166 | ! using 'z2sig.def' (or 'esasig.def) file |
---|
167 | |
---|
168 | logical pressure_exner |
---|
169 | ! compute pressure inside layers using Exner function, else use mean |
---|
170 | ! of pressure values at interfaces |
---|
171 | |
---|
172 | !----------------------------------------------------------------------- |
---|
173 | |
---|
174 | integer ngroup |
---|
175 | parameter (ngroup=3) |
---|
176 | |
---|
177 | |
---|
178 | real pbaru(iip1,jjp1,llm),pbarv(iip1,jjm,llm) |
---|
179 | real pext(iip1,jjp1,llm) |
---|
180 | |
---|
181 | real pbarum(iip1,jjp1,llm),pbarvm(iip1,jjm,llm) |
---|
182 | real wm(iip1,jjp1,llm) |
---|
183 | |
---|
184 | real,save :: zconvm(iip1,jjp1,llm) |
---|
185 | real,save :: zconvmm(iip1,jjp1,llm) |
---|
186 | |
---|
187 | real uu |
---|
188 | |
---|
189 | integer i,j,l |
---|
190 | |
---|
191 | logical firstcall,groupe_ok |
---|
192 | save firstcall,groupe_ok |
---|
193 | c$OMP THREADPRIVATE(firstcall,groupe_ok) |
---|
194 | |
---|
195 | data firstcall/.true./ |
---|
196 | data groupe_ok/.true./ |
---|
197 | |
---|
198 | integer ijb,ije,jjb,jje |
---|
199 | |
---|
200 | if (iim==1) then |
---|
201 | groupe_ok=.false. |
---|
202 | endif |
---|
203 | |
---|
204 | if (firstcall) then |
---|
205 | if (groupe_ok) then |
---|
206 | if(mod(iim,2**ngroup).ne.0) stop'probleme du nombre de point' |
---|
207 | endif |
---|
208 | firstcall=.false. |
---|
209 | endif |
---|
210 | |
---|
211 | c Champs 1D |
---|
212 | |
---|
213 | call convflu_p(pbaru,pbarv,llm,zconvm) |
---|
214 | |
---|
215 | c |
---|
216 | c call scopy(ijp1llm,zconvm,1,zconvmm,1) |
---|
217 | c call scopy(ijmllm,pbarv,1,pbarvm,1) |
---|
218 | |
---|
219 | jjb=jj_begin |
---|
220 | jje=jj_end |
---|
221 | |
---|
222 | c$OMP DO SCHEDULE(STATIC,OMP_CHUNK) |
---|
223 | do l=1,llm |
---|
224 | zconvmm(:,jjb:jje,l)=zconvm(:,jjb:jje,l) |
---|
225 | enddo |
---|
226 | c$OMP END DO NOWAIT |
---|
227 | |
---|
228 | if (groupe_ok) then |
---|
229 | call groupeun_p(jjp1,llm,jjb,jje,zconvmm) |
---|
230 | endif |
---|
231 | |
---|
232 | jjb=jj_begin-1 |
---|
233 | jje=jj_end |
---|
234 | if (pole_nord) jjb=jj_begin |
---|
235 | if (pole_sud) jje=jj_end-1 |
---|
236 | c$OMP DO SCHEDULE(STATIC,OMP_CHUNK) |
---|
237 | do l=1,llm |
---|
238 | pbarvm(:,jjb:jje,l)=pbarv(:,jjb:jje,l) |
---|
239 | enddo |
---|
240 | c$OMP END DO NOWAIT |
---|
241 | |
---|
242 | if (groupe_ok) then |
---|
243 | call groupeun_p(jjm,llm,jjb,jje,pbarvm) |
---|
244 | endif |
---|
245 | |
---|
246 | c Champs 3D |
---|
247 | |
---|
248 | jjb=jj_begin |
---|
249 | jje=jj_end |
---|
250 | if (pole_nord) jjb=jj_begin+1 |
---|
251 | if (pole_sud) jje=jj_end-1 |
---|
252 | |
---|
253 | c$OMP DO SCHEDULE(STATIC,OMP_CHUNK) |
---|
254 | do l=1,llm |
---|
255 | do j=jjb,jje |
---|
256 | uu=pbaru(iim,j,l) |
---|
257 | do i=1,iim |
---|
258 | uu=uu+pbarvm(i,j,l)-pbarvm(i,j-1,l)-zconvmm(i,j,l) |
---|
259 | pbarum(i,j,l)=uu |
---|
260 | c zconvm(i,j,l ) = xflu(i-1,j,l)-xflu(i,j,l)+ |
---|
261 | c * yflu(i,j,l)-yflu(i,j-1,l) |
---|
262 | enddo |
---|
263 | pbarum(iip1,j,l)=pbarum(1,j,l) |
---|
264 | enddo |
---|
265 | enddo |
---|
266 | c$OMP END DO NOWAIT |
---|
267 | |
---|
268 | c integration de la convergence de masse de haut en bas ...... |
---|
269 | |
---|
270 | jjb=jj_begin |
---|
271 | jje=jj_end |
---|
272 | |
---|
273 | c$OMP BARRIER |
---|
274 | c$OMP MASTER |
---|
275 | do l = llm-1,1,-1 |
---|
276 | do j=jjb,jje |
---|
277 | do i=1,iip1 |
---|
278 | zconvmm(i,j,l)=zconvmm(i,j,l)+zconvmm(i,j,l+1) |
---|
279 | enddo |
---|
280 | enddo |
---|
281 | enddo |
---|
282 | |
---|
283 | if (.not. pole_sud) then |
---|
284 | zconvmm(:,jj_end+1,:)=0 |
---|
285 | cym wm(:,jj_end+1,:)=0 |
---|
286 | endif |
---|
287 | |
---|
288 | c$OMP END MASTER |
---|
289 | c$OMP BARRIER |
---|
290 | |
---|
291 | CALL vitvert_p(zconvmm(1,1,1),wm(1,1,1)) |
---|
292 | |
---|
293 | return |
---|
294 | end |
---|
295 | |
---|