1 | MODULE topmons_mod |
---|
2 | |
---|
3 | IMPLICIT NONE |
---|
4 | |
---|
5 | ! sub-grid scale mountain mesh fraction |
---|
6 | REAL, SAVE, ALLOCATABLE :: alpha_hmons(:) |
---|
7 | |
---|
8 | !$OMP THREADPRIVATE(alpha_hmons) |
---|
9 | |
---|
10 | CONTAINS |
---|
11 | |
---|
12 | !======================================================================= |
---|
13 | ! ENTRAINMENT OF DUST ABOVE THE SUB-GRID SCALE TOPOGRAPHY |
---|
14 | !======================================================================= |
---|
15 | ! calculation of the vertical wind velocity of topdust tracers |
---|
16 | ! transport scheme of topdust tracers |
---|
17 | ! detrainement of topdust into background dust |
---|
18 | ! ----------------------------------------------------------------------- |
---|
19 | ! Authors: M. Vals; C. Wang; T. Bertrand; F. Forget |
---|
20 | ! Institution: Laboratoire de Meteorologie Dynamique (LMD) Paris, France |
---|
21 | ! ----------------------------------------------------------------------- |
---|
22 | |
---|
23 | SUBROUTINE topmons(ngrid,nlayer,nq,ptime,ptimestep, & |
---|
24 | pq,pdq,pt,pdt,pplev,pplay,pzlev, & |
---|
25 | pzlay,pdtsw,pdtlw, & |
---|
26 | ! input for radiative transfer |
---|
27 | icount,zday,zls,tsurf,igout,aerosol, & |
---|
28 | tauscaling,dust_rad_adjust, & |
---|
29 | ! input sub-grid scale rocket dust storm |
---|
30 | totstormfract,clearatm, & |
---|
31 | ! input sub-grid scale cloud |
---|
32 | clearsky,totcloudfrac, & |
---|
33 | ! input sub-grid scale mountain |
---|
34 | nohmons,hsummit, & |
---|
35 | ! output |
---|
36 | pdqtop,wfin,dsodust,dsords,dsotop, & |
---|
37 | tau_pref_scenario,tau_pref_gcm) |
---|
38 | |
---|
39 | USE tracer_mod, only: igcm_topdust_mass,igcm_topdust_number & |
---|
40 | ,igcm_dust_mass,igcm_dust_number & |
---|
41 | ,rho_dust |
---|
42 | USE comcstfi_h, only: r,g,cpp,rcp |
---|
43 | USE dimradmars_mod, only: albedo,naerkind |
---|
44 | USE comsaison_h, only: dist_sol,mu0,fract |
---|
45 | USE surfdat_h, only: emis,co2ice,hmons,summit |
---|
46 | USE callradite_mod, only: callradite |
---|
47 | |
---|
48 | IMPLICIT NONE |
---|
49 | |
---|
50 | !-------------------------------------------------------- |
---|
51 | ! Input Variables |
---|
52 | !-------------------------------------------------------- |
---|
53 | |
---|
54 | INTEGER, INTENT(IN) :: ngrid ! number of horizontal grid points |
---|
55 | INTEGER, INTENT(IN) :: nlayer ! number of vertical grid points |
---|
56 | INTEGER, INTENT(IN) :: nq ! number of tracer species |
---|
57 | REAL, INTENT(IN) :: ptime |
---|
58 | REAL, INTENT(IN) :: ptimestep |
---|
59 | |
---|
60 | REAL, INTENT(IN) :: pq(ngrid,nlayer,nq) ! advected field nq |
---|
61 | REAL, INTENT(IN) :: pdq(ngrid,nlayer,nq)! tendancy field pq |
---|
62 | REAL, INTENT(IN) :: pt(ngrid,nlayer) ! temperature at mid-layer (K) |
---|
63 | REAL, INTENT(IN) :: pdt(ngrid,nlayer) ! tendancy temperature at mid-layer (K) |
---|
64 | |
---|
65 | REAL, INTENT(IN) :: pplay(ngrid,nlayer) ! pressure at middle of the layers |
---|
66 | REAL, INTENT(IN) :: pplev(ngrid,nlayer+1) ! pressure at intermediate levels |
---|
67 | REAL, INTENT(IN) :: pzlay(ngrid,nlayer) ! altitude at the middle of the layers |
---|
68 | REAL, INTENT(IN) :: pzlev(ngrid,nlayer+1) ! altitude at layer boundaries |
---|
69 | |
---|
70 | REAL, INTENT(IN) :: pdtsw(ngrid,nlayer) ! (K/s) env |
---|
71 | REAL, INTENT(IN) :: pdtlw(ngrid,nlayer) ! (K/s) env |
---|
72 | |
---|
73 | ! input for second radiative transfer |
---|
74 | INTEGER, INTENT(INOUT) :: icount |
---|
75 | REAL, INTENT(IN) :: zday |
---|
76 | REAL, INTENT(IN) :: zls |
---|
77 | REAL, INTENT(IN) :: tsurf(ngrid) |
---|
78 | INTEGER, INTENT(IN) :: igout |
---|
79 | REAL, INTENT(INOUT) :: tauscaling(ngrid) |
---|
80 | REAL,INTENT(OUT) :: dust_rad_adjust(ngrid) |
---|
81 | ! input sub-grid scale rocket dust storm |
---|
82 | LOGICAL, INTENT(IN) :: clearatm |
---|
83 | REAL, INTENT(IN) :: totstormfract(ngrid) |
---|
84 | ! input sub-grid scale water ice clouds |
---|
85 | LOGICAL, INTENT(IN) :: clearsky |
---|
86 | REAL, INTENT(IN) :: totcloudfrac(ngrid) |
---|
87 | ! input sub-grid scale mountain |
---|
88 | LOGICAL, INTENT(IN) :: nohmons |
---|
89 | REAL, INTENT(IN) :: hsummit(ngrid) |
---|
90 | |
---|
91 | !-------------------------------------------------------- |
---|
92 | ! Output Variables |
---|
93 | !-------------------------------------------------------- |
---|
94 | |
---|
95 | REAL, INTENT(OUT) :: pdqtop(ngrid,nlayer,nq) ! tendancy field for dust |
---|
96 | REAL, INTENT(OUT) :: wfin(ngrid,nlayer+1) ! final vertical wind velocity: combination of both wup and wrad |
---|
97 | ! wfin < 0 means upward wind (in pressure levels/s) |
---|
98 | |
---|
99 | ! output for second radiative transfer |
---|
100 | REAL, INTENT(OUT) :: aerosol(ngrid,nlayer,naerkind) |
---|
101 | REAL, INTENT(OUT) :: dsodust(ngrid,nlayer) |
---|
102 | REAL, INTENT(OUT) :: dsords(ngrid,nlayer) |
---|
103 | REAL, INTENT(OUT) :: dsotop(ngrid,nlayer) |
---|
104 | REAL,INTENT(OUT) :: tau_pref_scenario(ngrid) ! prescribed dust column |
---|
105 | ! visible opacity at odpref |
---|
106 | REAL,INTENT(OUT) :: tau_pref_gcm(ngrid) ! dust column visible opacity at |
---|
107 | ! odpref in the GCM |
---|
108 | |
---|
109 | !-------------------------------------------------------- |
---|
110 | ! Local variables |
---|
111 | !-------------------------------------------------------- |
---|
112 | LOGICAL,SAVE :: firstcall=.true. |
---|
113 | |
---|
114 | !$OMP THREADPRIVATE(firstcall) |
---|
115 | |
---|
116 | INTEGER l,ig,tsub,iq,ll |
---|
117 | REAL zq0(ngrid,nlayer,nq) ! initial tracers |
---|
118 | REAL zq(ngrid,nlayer,nq) ! updated tracers |
---|
119 | REAL,PARAMETER:: mu0lim=0.001 ! solar zenith angle limit to determine the daytime |
---|
120 | |
---|
121 | ! Variables for radiative transfer |
---|
122 | REAL fluxsurf_lw1(ngrid) |
---|
123 | REAL fluxsurf_sw1(ngrid,2) |
---|
124 | REAL fluxtop_lw1(ngrid) |
---|
125 | REAL fluxtop_sw1(ngrid,2) |
---|
126 | REAL tau(ngrid,naerkind) |
---|
127 | REAL taucloudtes(ngrid) |
---|
128 | REAL rdust(ngrid,nlayer) |
---|
129 | REAL rstormdust(ngrid,nlayer) |
---|
130 | REAL rtopdust(ngrid,nlayer) |
---|
131 | REAL rice(ngrid,nlayer) |
---|
132 | REAL nuice(ngrid,nlayer) |
---|
133 | DOUBLE PRECISION riceco2(ngrid,nlayer) |
---|
134 | REAL nuiceco2(ngrid,nlayer) |
---|
135 | |
---|
136 | ! Temperature profile |
---|
137 | REAL zt(ngrid,nlayer) ! actual temperature at mid-layer (K) |
---|
138 | REAL ztlev(ngrid,nlayer) ! temperature at intermediate levels l+1/2 without last level |
---|
139 | REAL zdtlw1(ngrid,nlayer) ! (K/s) topdust |
---|
140 | REAL zdtsw1(ngrid,nlayer) ! (K/s) topdust |
---|
141 | REAL zdtlw1_lev(ngrid,nlayer),zdtsw1_lev(ngrid,nlayer) ! rad. heating rate at intermediate levels l+1/2 |
---|
142 | REAL zdtlw_lev(ngrid,nlayer),zdtsw_lev(ngrid,nlayer) ! rad. heating rate at intermediate levels l+1/2 |
---|
143 | REAL zdtvert(ngrid,nlayer) ! dT/dz , lapse rate |
---|
144 | REAL deltahr(ngrid,nlayer+1) ! extra heating between the concentrated and non-concentrated dust |
---|
145 | |
---|
146 | REAL newzt(ngrid,nlayer) ! recalculated temperature with extra radiative heating |
---|
147 | REAL t_top(ngrid,nlayer) ! estimated temperature above the mountain |
---|
148 | REAL t_top_lev(ngrid,nlayer)! estimated temperature above the mountain at the levels |
---|
149 | REAL dt_top(ngrid) ! temperature difference between the summit and the environment |
---|
150 | INTEGER lmons(ngrid) ! layer of the mountain's slope |
---|
151 | INTEGER lsummit(ngrid) ! layer of the mountain's summit |
---|
152 | REAL t_env(ngrid,nlayer) ! recalculated temperature of the environment air "next" to the mountain |
---|
153 | |
---|
154 | ! Vertical velocity profile |
---|
155 | REAL wrad(ngrid,nlayer+1) ! vertical wind velocity resulting from radiative heating |
---|
156 | REAL wup(ngrid,nlayer+1) ! vertical wind velocity calculated above the sub-grid scale mountain |
---|
157 | REAL wup2(ngrid,nlayer+1) ! square of the vertical velocity calculated above the sub-grid scale mountain |
---|
158 | REAL w0(ngrid) ! vertical wind velocity convergence at the top of the sub-grid scale mountain |
---|
159 | INTEGER lwmax(ngrid) ! level of the maximum vertical velocity |
---|
160 | |
---|
161 | REAL,PARAMETER:: wmax=10. ! maximum vertical velocity resulting from radiative heating (m/s) |
---|
162 | REAL,PARAMETER:: secu=3. ! coefficient on wfin to avoid dust crossing many layers during subtimestep |
---|
163 | REAL,PARAMETER:: k0=0.25 !max 10m/s: k0=1. !max 3m/s: k0=0.25 |
---|
164 | REAL,PARAMETER:: k1=5.e-4 !max 10m/s: k1=5.e-4 !max 3m/s: k1=5.e-4 |
---|
165 | REAL,PARAMETER:: k2=5.e-3 !max 10m/s: k2=30.e-3 !max 3m/s: k2=5.e-3 |
---|
166 | |
---|
167 | ! Entrainment |
---|
168 | REAL entr(ngrid) ! entrainment flux |
---|
169 | REAL rho(ngrid,nlayer),rhobarz(ngrid,nlayer+1) ! density, density at the levels |
---|
170 | REAL masse(ngrid,nlayer) ! air mass |
---|
171 | REAL masse_pbl(ngrid) ! total air mass within the PBL |
---|
172 | REAL masse_pbl_dust_mass(ngrid),masse_pbl_dust_number(ngrid) ! total air mass + dust mass/number within the PBL (kg/m2) |
---|
173 | REAL dqm_mass_pbl(ngrid),dqm_number_pbl(ngrid) ! total mass of the entrained dust mass/number from the PBL (kg/m2) |
---|
174 | |
---|
175 | ! Van Leer |
---|
176 | REAL zq_dust_mass(ngrid,nlayer),zq_dust_number(ngrid,nlayer) ! mixing ratio background dust (kg/kg) |
---|
177 | REAL zq_topdust_mass(ngrid,nlayer),zq_topdust_number(ngrid,nlayer) ! mixing ratio topdust (kg/kg) |
---|
178 | REAL mr_dust_mass(ngrid,nlayer),mr_dust_number(ngrid,nlayer) ! mixing ratio background dust (kg/kg) |
---|
179 | REAL mr_topdust_mass(ngrid,nlayer),mr_topdust_number(ngrid,nlayer) ! mixing ratio background dust + concentrated topdust (kg/kg) |
---|
180 | |
---|
181 | REAL w(ngrid,nlayer) ! vertical flux above the sub-grid scale mountain |
---|
182 | REAL wqmass(ngrid,nlayer+1) ! flux mass |
---|
183 | REAL wqnumber(ngrid,nlayer+1) ! flux number |
---|
184 | REAL qbar_mass_pbl(ngrid) ! total dust mass mixing ratio within the PBL for Van Leer |
---|
185 | REAL qbar_number_pbl(ngrid) ! total dust number mixing ratio within the PBL for Van Leer |
---|
186 | REAL masse_col(nlayer) ! masse of the atmosphere in one column |
---|
187 | |
---|
188 | REAL dqvl_dust_mass(ngrid,nlayer) ! tendancy pdq dust mass after vertical transport |
---|
189 | REAL dqvl_dust_number(ngrid,nlayer) ! tendancy pdq dust number after vertical transport |
---|
190 | REAL dqvl_topdust_mass(ngrid,nlayer) ! tendancy pdq topdust mass after vertical transport |
---|
191 | REAL dqvl_topdust_number(ngrid,nlayer) ! tendancy pdq topdust number after vertical transport |
---|
192 | |
---|
193 | INTEGER nsubtimestep(ngrid) ! number of subtimestep when calling Van Leer scheme |
---|
194 | REAL subtimestep(ngrid) ! ptimestep/nsubtimestep |
---|
195 | REAL dtmax ! considered minimum time interval needed for the dust to cross one vertical layer |
---|
196 | |
---|
197 | ! Detrainment |
---|
198 | REAL coefdetrain(ngrid,nlayer) ! coefficient for detrainment : % of stormdust detrained |
---|
199 | REAL dqdet_topdust_mass(ngrid,nlayer) ! tendancy pdq topdust mass after detrainment only |
---|
200 | REAL dqdet_topdust_number(ngrid,nlayer) ! tendancy pdq topdust number after detrainment only |
---|
201 | |
---|
202 | ! Diagnostics |
---|
203 | REAL zlaywmax(ngrid) |
---|
204 | REAL detrain_rate(ngrid,nlayer) |
---|
205 | |
---|
206 | ! ********************************************************************** |
---|
207 | ! ********************************************************************** |
---|
208 | ! Parametrization of the entrainment by slope wind above the sub-grid |
---|
209 | ! scale topography |
---|
210 | ! ********************************************************************** |
---|
211 | ! ********************************************************************** |
---|
212 | !! 1. Radiative transfer in topdust |
---|
213 | !! 2. Compute vertical velocity for topdust |
---|
214 | !! 3. Vertical transport |
---|
215 | !! 4. Detrainment |
---|
216 | ! ********************************************************************** |
---|
217 | ! ********************************************************************** |
---|
218 | |
---|
219 | ! ********************************************************************** |
---|
220 | ! 0. Initializations |
---|
221 | ! ********************************************************************** |
---|
222 | aerosol(:,:,:)=0. |
---|
223 | pdqtop(:,:,:) = 0. |
---|
224 | dqvl_topdust_mass(:,:)=0. |
---|
225 | dqvl_topdust_number(:,:)=0. |
---|
226 | dqvl_dust_mass(:,:)=0. |
---|
227 | dqvl_dust_number(:,:)=0. |
---|
228 | dqdet_topdust_mass(:,:)=0. |
---|
229 | dqdet_topdust_number(:,:)=0. |
---|
230 | wup(:,:)=0. |
---|
231 | wup2(:,:)=0. |
---|
232 | wrad(:,:)=0. |
---|
233 | w0(:)=0. |
---|
234 | wfin(:,:)=0. |
---|
235 | w(:,:)=0. |
---|
236 | zq(:,:,:) = 0. |
---|
237 | zq0(:,:,:) = 0. |
---|
238 | entr(:) = 0. |
---|
239 | mr_dust_mass(:,:) = 0. |
---|
240 | mr_dust_number(:,:) = 0. |
---|
241 | mr_topdust_mass(:,:) = 0. |
---|
242 | mr_topdust_number(:,:) = 0. |
---|
243 | t_top(:,:) = 0. |
---|
244 | dt_top(:) = 0. |
---|
245 | newzt(:,:) = 0. |
---|
246 | lmons(:) = 1 |
---|
247 | lsummit(:) =1 |
---|
248 | lwmax(:) = 1 |
---|
249 | deltahr(:,:) = 0. |
---|
250 | zdtvert(:,:) = 0. |
---|
251 | wqmass(:,:) = 0. |
---|
252 | wqnumber(:,:) = 0. |
---|
253 | coefdetrain(:,:) = 1. |
---|
254 | dqm_mass_pbl(:)=0. |
---|
255 | dqm_number_pbl(:)=0. |
---|
256 | nsubtimestep(:)=1 |
---|
257 | masse_pbl(:)=0. |
---|
258 | detrain_rate(:,:) = 0. |
---|
259 | |
---|
260 | !! Update the initial temperature |
---|
261 | zt(1:ngrid,1:nlayer)=pt(1:ngrid,1:nlayer)+pdt(1:ngrid,1:nlayer)*ptimestep |
---|
262 | newzt(1:ngrid,1:nlayer)=zt(1:ngrid,1:nlayer) |
---|
263 | t_env(1:ngrid,1:nlayer)=zt(1:ngrid,1:nlayer) |
---|
264 | t_top(1:ngrid,1:nlayer)=zt(1:ngrid,1:nlayer) |
---|
265 | |
---|
266 | !! Update the initial mixing ratios |
---|
267 | zq0(1:ngrid,1:nlayer,1:nq)=pq(1:ngrid,1:nlayer,1:nq)+pdq(1:ngrid,1:nlayer,1:nq)*ptimestep ! update of the background dust after rocket dust storm scheme |
---|
268 | zq(1:ngrid,1:nlayer,1:nq)=zq0(1:ngrid,1:nlayer,1:nq) |
---|
269 | zq_dust_mass(1:ngrid,1:nlayer)=zq(1:ngrid,1:nlayer,igcm_dust_mass) |
---|
270 | zq_dust_number(1:ngrid,1:nlayer)=zq(1:ngrid,1:nlayer,igcm_dust_number) |
---|
271 | zq_topdust_mass(1:ngrid,1:nlayer)=zq(1:ngrid,1:nlayer,igcm_topdust_mass) |
---|
272 | zq_topdust_number(1:ngrid,1:nlayer)=zq(1:ngrid,1:nlayer,igcm_topdust_number) |
---|
273 | |
---|
274 | ! ---------------------------------------------------------------------- |
---|
275 | ! 1. Radiative heating |
---|
276 | ! ---------------------------------------------------------------------- |
---|
277 | |
---|
278 | ! ********************************************************************* |
---|
279 | ! 1.1. Call the second radiative transfer for topdust, obtain the extra heating |
---|
280 | ! ********************************************************************* |
---|
281 | CALL callradite(icount,ngrid,nlayer,nq,zday,zls,zq,albedo, & |
---|
282 | emis,mu0,pplev,pplay,pt,tsurf,fract,dist_sol,igout, & |
---|
283 | zdtlw1,zdtsw1,fluxsurf_lw1,fluxsurf_sw1,fluxtop_lw1, & |
---|
284 | fluxtop_sw1,tau_pref_scenario,tau_pref_gcm, & |
---|
285 | tau,aerosol,dsodust,tauscaling,dust_rad_adjust, & |
---|
286 | taucloudtes,rdust,rice,nuice,riceco2,nuiceco2,co2ice,rstormdust,rtopdust, & |
---|
287 | totstormfract,clearatm,dsords,dsotop,alpha_hmons,nohmons,& |
---|
288 | clearsky,totcloudfrac) |
---|
289 | ! ********************************************************************** |
---|
290 | ! 1.2. Compute radiative vertical velocity for entrained topdust |
---|
291 | ! ********************************************************************** |
---|
292 | DO ig=1,ngrid |
---|
293 | IF ( (mu0(ig) .gt. mu0lim) .and. (alpha_hmons(ig) .gt. 0.) ) THEN |
---|
294 | !! ********************************************************************** |
---|
295 | !! Temperature profile above the mountain and in the close environment |
---|
296 | call t_topmons(nlayer,summit(ig),hsummit(ig),hmons(ig),zt(ig,:),pzlay(ig,:), & |
---|
297 | t_top(ig,:),dt_top(ig),t_env(ig,:),lsummit(ig),lmons(ig)) |
---|
298 | !! ********************************************************************** |
---|
299 | !! Calculation of the extra heating : computing heating rates |
---|
300 | !! gradient at boundaries of each layer |
---|
301 | zdtlw1_lev(ig,1)=0. |
---|
302 | zdtsw1_lev(ig,1)=0. |
---|
303 | zdtlw_lev(ig,1)=0. |
---|
304 | zdtsw_lev(ig,1)=0. |
---|
305 | ztlev(ig,1)=zt(ig,1) |
---|
306 | t_top_lev(ig,1)=tsurf(ig) |
---|
307 | |
---|
308 | DO l=1,nlayer-1 |
---|
309 | !! Calculation for the background dust AND the concentrated topdust |
---|
310 | zdtlw1_lev(ig,l+1)=(zdtlw1(ig,l)*(pzlay(ig,l+1)-pzlev(ig,l+1))+ & |
---|
311 | zdtlw1(ig,l+1)*(pzlev(ig,l+1)-pzlay(ig,l))) / & |
---|
312 | (pzlay(ig,l+1)-pzlay(ig,l)) |
---|
313 | |
---|
314 | zdtsw1_lev(ig,l+1)=(zdtsw1(ig,l)*(pzlay(ig,l+1)-pzlev(ig,l+1))+ & |
---|
315 | zdtsw1(ig,l+1)*(pzlev(ig,l+1)-pzlay(ig,l))) / & |
---|
316 | (pzlay(ig,l+1)-pzlay(ig,l)) |
---|
317 | !! Calculation for the background dust only |
---|
318 | zdtlw_lev(ig,l+1)=(pdtlw(ig,l)*(pzlay(ig,l+1)-pzlev(ig,l+1))+ & |
---|
319 | pdtlw(ig,l+1)*(pzlev(ig,l+1)-pzlay(ig,l))) / & |
---|
320 | (pzlay(ig,l+1)-pzlay(ig,l)) |
---|
321 | |
---|
322 | zdtsw_lev(ig,l+1)=(pdtsw(ig,l)*(pzlay(ig,l+1)-pzlev(ig,l+1))+ & |
---|
323 | pdtsw(ig,l+1)*(pzlev(ig,l+1)-pzlay(ig,l))) / & |
---|
324 | (pzlay(ig,l+1)-pzlay(ig,l)) |
---|
325 | !! Temperature at the levels |
---|
326 | ztlev(ig,l+1)=(zt(ig,l)*(pzlay(ig,l+1)-pzlev(ig,l+1))+ & |
---|
327 | zt(ig,l+1)*(pzlev(ig,l+1)-pzlay(ig,l))) / & |
---|
328 | (pzlay(ig,l+1)-pzlay(ig,l)) |
---|
329 | !! Temperature above the mountain at the levels |
---|
330 | t_top_lev(ig,l+1)=(t_top(ig,l)*(pzlay(ig,l+1)-pzlev(ig,l+1))+ & |
---|
331 | t_top(ig,l+1)*(pzlev(ig,l+1)-pzlay(ig,l))) / & |
---|
332 | (pzlay(ig,l+1)-pzlay(ig,l)) |
---|
333 | ENDDO ! DO l=1,nlayer-1 |
---|
334 | !! ********************************************************************** |
---|
335 | !! Vertical velocity profile for extra heating balanced by adiabatic cooling |
---|
336 | !! Gradient at boundaries of each layer, start from surface |
---|
337 | zdtvert(ig,1)=0. |
---|
338 | DO l=1,nlayer-1 |
---|
339 | zdtvert(ig,l+1)=(t_top_lev(ig,l+1)-t_top_lev(ig,l))/(pzlay(ig,l+1)-pzlay(ig,l)) |
---|
340 | ENDDO |
---|
341 | !! Extra heating balanced by adiabatic cooling |
---|
342 | DO l=1,nlayer |
---|
343 | deltahr(ig,l)=(zdtlw1_lev(ig,l)+zdtsw1_lev(ig,l)) & |
---|
344 | -(zdtlw_lev(ig,l)+zdtsw_lev(ig,l)) |
---|
345 | wrad(ig,l)=-deltahr(ig,l)/(g/cpp+ & |
---|
346 | max(zdtvert(ig,l),-0.99*g/cpp)) |
---|
347 | !! Limit of the vertical wind velocity in case of lapse rate close to adiabatic |
---|
348 | wrad(ig,l)=max(wrad(ig,l),-wmax) |
---|
349 | wrad(ig,l)=min(wrad(ig,l),wmax) |
---|
350 | !! ATTENTION !! to simplify here, no downward velocity calculated |
---|
351 | if (wrad(ig,l).gt.0.) then |
---|
352 | wrad(ig,l)=0. |
---|
353 | endif |
---|
354 | ENDDO |
---|
355 | ENDIF ! IF ((mu0(ig) .gt. mu0lim) .and. alpha_hmons(ig) .gt. 0.) |
---|
356 | ENDDO ! DO ig=1,ngrid |
---|
357 | |
---|
358 | ! ---------------------------------------------------------------------- |
---|
359 | ! 2. Vertical velocity profile |
---|
360 | ! ---------------------------------------------------------------------- |
---|
361 | |
---|
362 | ! ********************************************************************** |
---|
363 | ! 2.1 Compute total vertical velocity for entrained topdust: buoyancy + radiative |
---|
364 | ! ********************************************************************** |
---|
365 | DO ig=1,ngrid |
---|
366 | IF ( (mu0(ig) .gt. mu0lim) .and. (alpha_hmons(ig) .gt. 0.) ) THEN |
---|
367 | !! Positive buoyancy: negative vertical velocity entrains UP |
---|
368 | IF (dt_top(ig) .gt. 0.) THEN |
---|
369 | !! Convergence of vertical winds at the top of the mountain |
---|
370 | w0(ig)=-k0*g*sqrt((dt_top(ig)/t_env(ig,lsummit(ig)))) |
---|
371 | !! Case the vertical velocity at the top of the mountain > wrad |
---|
372 | IF (w0(ig).lt.wrad(ig,lsummit(ig)+1)) then |
---|
373 | wup(ig,lsummit(ig)+1)=w0(ig) |
---|
374 | wfin(ig,lsummit(ig)+1)=w0(ig) |
---|
375 | newzt(ig,lsummit(ig))= t_top(ig,lsummit(ig)) |
---|
376 | !! Temperature profile from the top of the summit to the top of the atmosphere |
---|
377 | DO l=lsummit(ig)+1,nlayer-1 |
---|
378 | !ATTENTION: test without wrad, T is only decreasing because of the adiabatic gradient |
---|
379 | !newzt(ig,l)=newzt(ig,l-1)-(g/cpp)* & |
---|
380 | ! (pzlay(ig,l)-pzlay(ig,l-1)) |
---|
381 | !! Case of superadiabatic layer |
---|
382 | IF (zdtvert(ig,l) .lt. -g/cpp) then |
---|
383 | newzt(ig,l)=t_top(ig,l) |
---|
384 | !! New temperature due to radiative heating balanced by adiabatic cooling |
---|
385 | ELSE |
---|
386 | newzt(ig,l)=newzt(ig,l-1)+ & |
---|
387 | ( (deltahr(ig,l)/(-wfin(ig,l)))-g/cpp )* & |
---|
388 | ( pzlay(ig,l)-pzlay(ig,l-1) ) |
---|
389 | ENDIF ! (zdtvert(ig,l) .lt. -g/cpp) |
---|
390 | !! Vertical velocity profile due to the presence of the mountain with the new temperature |
---|
391 | !! Square of the vertical velocity profile with the new temperature |
---|
392 | wup2(ig,l+1)=(1.-2*k1*(pzlev(ig,l+1)-pzlev(ig,l)))*& |
---|
393 | (wup(ig,l)**2)+2.*k2*g*(pzlev(ig,l+1) & |
---|
394 | -pzlev(ig,l))*((newzt(ig,l) & |
---|
395 | -t_env(ig,l))/t_env(ig,l)) |
---|
396 | !! Vertical velocity profile if W2 > 0 |
---|
397 | IF (wup2(ig,l+1) .gt. 0.) THEN |
---|
398 | wup(ig,l+1)=-sqrt(wup2(ig,l+1)) |
---|
399 | wfin(ig,l+1)=wup(ig,l+1) |
---|
400 | IF (wup(ig,l+1) .gt. wrad(ig,l+1)) then |
---|
401 | DO ll=l,nlayer-1 |
---|
402 | newzt(ig,ll)=zt(ig,ll) |
---|
403 | wfin(ig,ll+1)=wrad(ig,ll+1) |
---|
404 | ENDDO |
---|
405 | EXIT |
---|
406 | ENDIF |
---|
407 | !! Vertical velocity profile if W2 < 0 |
---|
408 | ELSE |
---|
409 | DO ll=l,nlayer-1 |
---|
410 | newzt(ig,ll)=zt(ig,ll) |
---|
411 | wfin(ig,ll+1)=wrad(ig,ll+1) |
---|
412 | ENDDO |
---|
413 | EXIT |
---|
414 | ENDIF ! IF (wup2(ig,l+1) .gt. 0.) |
---|
415 | ENDDO ! DO l=lsummit(ig)+1,nlayer-1 |
---|
416 | !! Case the vertical velocity at the top of the mountain < wrad |
---|
417 | ELSE |
---|
418 | DO ll=lsummit(ig),nlayer-1 |
---|
419 | newzt(ig,ll)=zt(ig,ll) |
---|
420 | wfin(ig,ll+1)=wrad(ig,ll+1) |
---|
421 | ENDDO |
---|
422 | ENDIF !(w0(ig).lt.wrad(ig,lsummit(ig)+1)) |
---|
423 | !! Positive buoyancy: positive vertical velocity entrains DOWN |
---|
424 | ELSE |
---|
425 | DO l=lsummit(ig)+1,nlayer |
---|
426 | wfin(ig,l)=wrad(ig,l) |
---|
427 | ENDDO |
---|
428 | ENDIF ! IF (dt_top(ig) .gt. 0.) |
---|
429 | ! ********************************************************************** |
---|
430 | ! 2.2 Find the level lwmax of the maximum vertical velocity caused by the |
---|
431 | ! mountain presence (wup), where to entrain the dust from the PBL |
---|
432 | ! ********************************************************************** |
---|
433 | IF (minloc(wup(ig,:),dim=1).lt.nlayer) THEN |
---|
434 | lwmax(ig)=minloc(wup(ig,:),dim=1) |
---|
435 | ENDIF |
---|
436 | zlaywmax(ig)=pzlay(ig,lwmax(ig)) ! wup(lwmax) acts on level lwmax, going to layer lwmax |
---|
437 | ! ********************************************************************** |
---|
438 | ! 2.3 Compute the subtimestep to conserve the mass in the Van Leer transport |
---|
439 | ! ********************************************************************** |
---|
440 | !! Calculation of the subtimesteps |
---|
441 | dtmax=ptimestep |
---|
442 | DO l=lwmax(ig),nlayer |
---|
443 | IF (wfin(ig,l).lt.0.) THEN |
---|
444 | dtmax=min(dtmax,(pzlev(ig,l)-pzlev(ig,l-1))/ & |
---|
445 | (secu*abs(wfin(ig,l)))) |
---|
446 | ELSE IF (wfin(ig,l).gt.0.) then |
---|
447 | dtmax=min(dtmax,(pzlev(ig,l+1)-pzlev(ig,l))/ & |
---|
448 | (secu*abs(wfin(ig,l)))) |
---|
449 | ENDIF |
---|
450 | ENDDO |
---|
451 | nsubtimestep(ig)= int(ptimestep/dtmax) |
---|
452 | subtimestep(ig)=ptimestep/float(nsubtimestep(ig)) |
---|
453 | !! Mass flux generated by wup in kg/m2 |
---|
454 | DO l=1,nlayer!+1 |
---|
455 | w(ig,l)=wfin(ig,l)*pplev(ig,l)/(r*ztlev(ig,l)) & |
---|
456 | *subtimestep(ig) |
---|
457 | ENDDO ! l=1,nlayer |
---|
458 | |
---|
459 | ENDIF ! ((mu0(ig) .gt. mu0lim)... |
---|
460 | ENDDO ! DO ig=1,ngrid |
---|
461 | |
---|
462 | ! ---------------------------------------------------------------------- |
---|
463 | ! 3. Dust entrainment by slope winds above the mountain |
---|
464 | ! ---------------------------------------------------------------------- |
---|
465 | !! rho on the layer |
---|
466 | rho(:,:)=pplay(:,:)/(r*pt(:,:)) |
---|
467 | !! rhobarz(l) = rho(l+1/2): rho on the levels |
---|
468 | rhobarz(:,1)=rho(:,1) |
---|
469 | DO l=2,nlayer |
---|
470 | rhobarz(:,l)=pplev(:,l)/(r*0.5*(pt(:,l)+pt(:,l-1))) |
---|
471 | ENDDO |
---|
472 | rhobarz(:,nlayer+1)=0. !top of the model |
---|
473 | !! Mass computation |
---|
474 | DO l=1,nlayer |
---|
475 | masse(:,l)=(pplev(:,l)-pplev(:,l+1))/g |
---|
476 | ENDDO |
---|
477 | |
---|
478 | DO ig=1,ngrid |
---|
479 | IF ( (mu0(ig) .gt. mu0lim) .and. (alpha_hmons(ig) .gt. 0.) ) THEN |
---|
480 | !! Total air mass within the PBL before entrainment (=> by PBL we mean between the surface and the layer where the vertical wind is maximum) |
---|
481 | masse_pbl(ig)=0. |
---|
482 | DO l=1,lwmax(ig)-1 |
---|
483 | masse_pbl(ig)=masse_pbl(ig)+(pplev(ig,l)-pplev(ig,l+1))/g |
---|
484 | ENDDO |
---|
485 | ENDIF ! ((mu0(ig) .gt. mu0lim)... |
---|
486 | ENDDO ! ig=1,ngrid |
---|
487 | ! ********************************************************************** |
---|
488 | ! 3.1. Entrainment |
---|
489 | ! ********************************************************************** |
---|
490 | DO ig=1,ngrid |
---|
491 | IF ( (mu0(ig) .gt. mu0lim) .and. (alpha_hmons(ig) .gt. 0.) .and. (masse_pbl(ig) .gt. 0.) ) THEN |
---|
492 | !! Transport of background dust + concentrated topdust above lwmax |
---|
493 | DO l=lwmax(ig),nlayer |
---|
494 | mr_dust_mass(ig,l) = zq_dust_mass(ig,l) |
---|
495 | mr_dust_number(ig,l) = zq_dust_number(ig,l) |
---|
496 | mr_topdust_mass(ig,l) = zq_dust_mass(ig,l) & |
---|
497 | +zq_topdust_mass(ig,l)/alpha_hmons(ig) |
---|
498 | mr_topdust_number(ig,l) = zq_dust_number(ig,l) & |
---|
499 | +zq_topdust_number(ig,l)/alpha_hmons(ig) |
---|
500 | ENDDO ! l=lwmax(ig),nlayer |
---|
501 | !! Entrainment flux from the PBL |
---|
502 | entr(ig) = alpha_hmons(ig)*wup(ig,lwmax(ig))*rhobarz(ig,lwmax(ig))/masse_pbl(ig) |
---|
503 | !! If entrains more than available masse_pbl (abs(entr) x mass_pbl x ptimestep > masse_pbl) |
---|
504 | if (abs(entr(ig)) .gt. 1./ptimestep) then |
---|
505 | entr(ig) = -(1./ptimestep) |
---|
506 | end if |
---|
507 | ! ********************************************************************** |
---|
508 | ! 3.2. Vertical transport: Van Leer |
---|
509 | ! ********************************************************************** |
---|
510 | !! Loop over the subtimesteps |
---|
511 | DO tsub=1,nsubtimestep(ig) |
---|
512 | !! qbar_pbl: total mixing ratio within the PBL |
---|
513 | masse_pbl_dust_mass(ig)=0. |
---|
514 | masse_pbl_dust_number(ig)=0. |
---|
515 | DO l=1,lwmax(ig)-1 |
---|
516 | masse_pbl_dust_mass(ig)=masse_pbl_dust_mass(ig)+zq_dust_mass(ig,l)*(pplev(ig,l)-pplev(ig,l+1))/g |
---|
517 | masse_pbl_dust_number(ig)=masse_pbl_dust_number(ig)+zq_dust_number(ig,l)*(pplev(ig,l)-pplev(ig,l+1))/g |
---|
518 | ENDDO |
---|
519 | qbar_mass_pbl(ig)=masse_pbl_dust_mass(ig)/masse_pbl(ig) |
---|
520 | qbar_number_pbl(ig)=masse_pbl_dust_mass(ig)/masse_pbl(ig) |
---|
521 | !! integration of the equation dq/dt = -(entr)q, with entr constant --> w0 < 0 when up so the equation is dq/dt = (entr)q |
---|
522 | dqm_mass_pbl(:)=0. |
---|
523 | dqm_number_pbl(:)=0. |
---|
524 | DO l=1,lwmax(ig)-1 ! this time we take the dust from the surface up to the level, where the velocity is maximum |
---|
525 | !! Total dust entrained from the PBL: |
---|
526 | dqm_mass_pbl(ig)=dqm_mass_pbl(ig)+(1.-exp(entr(ig)*subtimestep(ig)))*zq_dust_mass(ig,l)*(pplev(ig,l)-pplev(ig,l+1))/g |
---|
527 | dqm_number_pbl(ig)=dqm_number_pbl(ig)+(1.-exp(entr(ig)*subtimestep(ig)))*zq_dust_number(ig,l)*(pplev(ig,l)-pplev(ig,l+1))/g |
---|
528 | !! Background dust after entrainment (explicit: qdust(ig,l)=qdust(ig,l)-entr(ig)*qdust(ig,l)*ptimestep) |
---|
529 | zq_dust_mass(ig,l)=zq_dust_mass(ig,l)-(1.-exp(entr(ig)*subtimestep(ig)))*zq_dust_mass(ig,l) |
---|
530 | zq_dust_number(ig,l)=zq_dust_number(ig,l)-(1.-exp(entr(ig)*subtimestep(ig)))*zq_dust_number(ig,l) |
---|
531 | ENDDO |
---|
532 | !! Van Leer scheme (from lwmax to nlayer) |
---|
533 | wqmass(ig,:)=0. |
---|
534 | wqnumber(ig,:)=0. |
---|
535 | CALL van_leer(nlayer,mr_topdust_mass(ig,:),2.,lwmax(ig), & |
---|
536 | masse(ig,:),w(ig,:),wqmass(ig,:),masse_pbl(ig),dqm_mass_pbl(ig),alpha_hmons(ig),qbar_mass_pbl(ig)) |
---|
537 | CALL van_leer(nlayer,mr_topdust_number(ig,:),2.,lwmax(ig), & |
---|
538 | masse(ig,:),w(ig,:),wqnumber(ig,:),masse_pbl(ig),dqm_number_pbl(ig),alpha_hmons(ig),qbar_number_pbl(ig)) |
---|
539 | ENDDO !tsub=... |
---|
540 | ! ********************************************************************** |
---|
541 | ! 3.3. Re-calculation of the mixing ratios after vertical transport |
---|
542 | ! ********************************************************************** |
---|
543 | !! Redistribution from lwmax to nlayer |
---|
544 | DO l=lwmax(ig),nlayer |
---|
545 | !! General and "healthy" case |
---|
546 | IF (mr_topdust_mass(ig,l).ge.mr_dust_mass(ig,l)) THEN |
---|
547 | zq_dust_mass(ig,l) = mr_dust_mass(ig,l) |
---|
548 | zq_dust_number(ig,l) = mr_dust_number(ig,l) |
---|
549 | zq_topdust_mass(ig,l) = alpha_hmons(ig)*(mr_topdust_mass(ig,l)-mr_dust_mass(ig,l)) |
---|
550 | zq_topdust_number(ig,l) = alpha_hmons(ig)*(mr_topdust_number(ig,l)-mr_dust_number(ig,l)) |
---|
551 | !! Particular case: there is less than initial dust mixing ratio after the vertical transport |
---|
552 | ELSE |
---|
553 | zq_dust_mass(ig,l) = (1.-alpha_hmons(ig))*mr_dust_mass(ig,l)+alpha_hmons(ig)*mr_topdust_mass(ig,l) |
---|
554 | zq_dust_number(ig,l) = (1.-alpha_hmons(ig))*mr_dust_number(ig,l)+alpha_hmons(ig)*mr_topdust_number(ig,l) |
---|
555 | zq_topdust_mass(ig,l) = 0. |
---|
556 | zq_topdust_number(ig,l) = 0. |
---|
557 | ENDIF |
---|
558 | ENDDO ! DO l=lwmax(ig),nlayer |
---|
559 | ! ********************************************************************** |
---|
560 | ! 3.4. Calculation of the tendencies of the vertical transport from the surface up to the top of the atmosphere |
---|
561 | ! ********************************************************************** |
---|
562 | DO l=1,nlayer |
---|
563 | dqvl_topdust_mass(ig,l) = (zq_topdust_mass(ig,l)- & |
---|
564 | zq0(ig,l,igcm_topdust_mass)) /ptimestep |
---|
565 | dqvl_topdust_number(ig,l) = (zq_topdust_number(ig,l)- & |
---|
566 | zq0(ig,l,igcm_topdust_number)) /ptimestep |
---|
567 | dqvl_dust_mass(ig,l) = (zq_dust_mass(ig,l)-zq0(ig,l,igcm_dust_mass)) /ptimestep |
---|
568 | dqvl_dust_number(ig,l) = (zq_dust_number(ig,l)-zq0(ig,l,igcm_dust_number)) /ptimestep |
---|
569 | ENDDO |
---|
570 | |
---|
571 | ENDIF ! ((mu0(ig) .gt. mu0lim)... |
---|
572 | ENDDO ! ig=1,ngrid |
---|
573 | |
---|
574 | ! ---------------------------------------------------------------------- |
---|
575 | ! 4. Detrainment: topdust is converted to background dust |
---|
576 | ! ---------------------------------------------------------------------- |
---|
577 | |
---|
578 | !! ********************************************************************** |
---|
579 | !! 4.1 Compute the coefficient of detrainment |
---|
580 | !! ********************************************************************** |
---|
581 | DO ig=1,ngrid |
---|
582 | !DO l=lwmax(ig),nlayer-1 |
---|
583 | DO l=1,nlayer!-1 |
---|
584 | !! Detrainment during the day |
---|
585 | IF ( (mu0(ig) .gt. mu0lim) .and. (zq_topdust_mass(ig,l) .gt. zq_dust_mass(ig,l)*0.01)) THEN |
---|
586 | coefdetrain(ig,l)=1.*( rhobarz(ig,l+1)*abs(wfin(ig,l+1)) - rhobarz(ig,l)*abs(wfin(ig,l)) ) / masse(ig,l) |
---|
587 | !! Detrainment when abs(w(l)) > abs(w(l+1)), i.e. coefdetrain < 0 |
---|
588 | if ( coefdetrain(ig,l).lt.0.) then |
---|
589 | dqdet_topdust_mass(ig,l)=-(1.-exp(coefdetrain(ig,l)*ptimestep))*zq_topdust_mass(ig,l)/ptimestep |
---|
590 | dqdet_topdust_number(ig,l)=-(1.-exp(coefdetrain(ig,l)*ptimestep))*zq_topdust_number(ig,l)/ptimestep |
---|
591 | ! for diagnostics |
---|
592 | detrain_rate(ig,l)=(1.-exp(coefdetrain(ig,l)*ptimestep)) |
---|
593 | !! One cannot detrain more than available topdust |
---|
594 | if ( (abs(dqdet_topdust_mass(ig,l)).gt.zq_topdust_mass(ig,l)/ptimestep) .or. (abs(dqdet_topdust_number(ig,l)).gt.zq_topdust_number(ig,l)/ptimestep) ) then |
---|
595 | dqdet_topdust_mass(ig,l)=-zq_topdust_mass(ig,l)/ptimestep |
---|
596 | dqdet_topdust_number(ig,l)=-zq_topdust_number(ig,l)/ptimestep |
---|
597 | detrain_rate(ig,l)=1. |
---|
598 | endif |
---|
599 | !else!entrainment |
---|
600 | ! dqdet_topdust_mass(ig,l)=-(1.-exp(coefdetrain(ig,l)*ptimestep))*zq_dust_mass(ig,l)/ptimestep |
---|
601 | ! dqdet_topdust_number(ig,l)=-(1.-exp(coefdetrain(ig,l)*ptimestep))*zq_dust_number(ig,l)/ptimestep |
---|
602 | endif |
---|
603 | !! Full detrainment during the night imposed |
---|
604 | ELSE |
---|
605 | dqdet_topdust_mass(ig,l)=-zq_topdust_mass(ig,l)/ptimestep |
---|
606 | dqdet_topdust_number(ig,l)=-zq_topdust_number(ig,l)/ptimestep |
---|
607 | detrain_rate(ig,l)=1. |
---|
608 | ENDIF |
---|
609 | ENDDO ! DO l=1,nlayer |
---|
610 | ENDDO ! DO ig=1,ngrid |
---|
611 | |
---|
612 | ! ---------------------------------------------------------------------- |
---|
613 | ! 5. Final tendencies: vertical transport + detrainment |
---|
614 | ! ---------------------------------------------------------------------- |
---|
615 | |
---|
616 | DO ig=1,ngrid |
---|
617 | DO l=1,nlayer |
---|
618 | pdqtop(ig,l,igcm_topdust_mass)=dqdet_topdust_mass(ig,l) & |
---|
619 | +dqvl_topdust_mass(ig,l) |
---|
620 | pdqtop(ig,l,igcm_topdust_number)=dqdet_topdust_number(ig,l) & |
---|
621 | +dqvl_topdust_number(ig,l) |
---|
622 | pdqtop(ig,l,igcm_dust_mass)= -dqdet_topdust_mass(ig,l) & |
---|
623 | +dqvl_dust_mass(ig,l) |
---|
624 | pdqtop(ig,l,igcm_dust_number)= -dqdet_topdust_number(ig,l) & |
---|
625 | +dqvl_dust_number(ig,l) |
---|
626 | ENDDO ! DO l=1,nlayer |
---|
627 | ENDDO ! DO ig=1,ngrid |
---|
628 | |
---|
629 | |
---|
630 | ! ********************************************************************** |
---|
631 | ! WRITEDIAGFI |
---|
632 | ! ********************************************************************** |
---|
633 | IF (ngrid.eq.1) THEN |
---|
634 | CALL WRITEDIAGFI(ngrid,'wup_top', & |
---|
635 | 'wup_top','',1,wup(:,:)) |
---|
636 | CALL WRITEDIAGFI(ngrid,'wrad_top', & |
---|
637 | 'wrad_top','',1,wrad(:,:)) |
---|
638 | CALL WRITEDIAGFI(ngrid,'wfin_top', & |
---|
639 | 'wfin_top','',1,wfin(:,:)) |
---|
640 | CALL WRITEDIAGFI(ngrid,'wlwmax_top', & |
---|
641 | 'wlwmax_top','',0,wup(:,lwmax(1))) |
---|
642 | CALL WRITEDIAGFI(ngrid,'w0_top', & |
---|
643 | 'w0_top','',0,wup(:,lsummit(1)+1)) |
---|
644 | CALL WRITEDIAGFI(ngrid,'alpha_hmons', & |
---|
645 | 'alpha_hmons','',0,alpha_hmons) |
---|
646 | CALL WRITEDIAGFI(ngrid,'zt_top', & |
---|
647 | 'zt_top','',1,t_top(:,:)) |
---|
648 | CALL WRITEDIAGFI(ngrid,'dt_top', & |
---|
649 | 'dt_top','',0,dt_top(:)) |
---|
650 | CALL WRITEDIAGFI(ngrid,'envtemp', & |
---|
651 | 'envtemp','',1,zt(:,:)) |
---|
652 | CALL WRITEDIAGFI(ngrid,'t_env', & |
---|
653 | 't_env','',1,t_env(:,:)) |
---|
654 | CALL WRITEDIAGFI(ngrid,'newzt_top', & |
---|
655 | 'newzt_top','',1,newzt(:,:)) |
---|
656 | CALL WRITEDIAGFI(ngrid,'deltahr_top', & |
---|
657 | 'deltahr_top','',1,deltahr(:,:)) |
---|
658 | CALL WRITEDIAGFI(ngrid,'rholev', & |
---|
659 | 'rholev','',1,rho(:,:)) |
---|
660 | CALL WRITEDIAGFI(ngrid,'rhobarz', & |
---|
661 | 'rhobarz','',1,rhobarz(:,:)) |
---|
662 | CALL WRITEDIAGFI(ngrid,'zlsummit', & |
---|
663 | 'zlsummit','',0,pzlay(:,lsummit(1))) |
---|
664 | CALL WRITEDIAGFI(ngrid,'zlwmax', & |
---|
665 | 'zlwmax','',0,pzlay(:,lwmax(1))) |
---|
666 | CALL WRITEDIAGFI(ngrid,'pzlev_top', & |
---|
667 | 'pzlev_top','',1,pzlev(:,:)) |
---|
668 | CALL WRITEDIAGFI(ngrid,'coefdetrain', & |
---|
669 | 'coefdetrain','',1,coefdetrain(:,:)) |
---|
670 | CALL WRITEDIAGFI(ngrid,'zdtvert', & |
---|
671 | 'zdtvert','',1,zdtvert(:,:)) |
---|
672 | CALL WRITEDIAGFI(ngrid,'entr', & |
---|
673 | 'entr','',0,entr(:)) |
---|
674 | ELSE |
---|
675 | CALL WRITEDIAGFI(ngrid,'wup_top', & |
---|
676 | 'wup_top','',3,wup(:,:)) |
---|
677 | CALL WRITEDIAGFI(ngrid,'wup2_top', & |
---|
678 | 'wup2_top','',3,wup2(:,:)) |
---|
679 | CALL WRITEDIAGFI(ngrid,'wrad_top', & |
---|
680 | 'wrad_top','',3,wrad(:,:)) |
---|
681 | CALL WRITEDIAGFI(ngrid,'wfin_top', & |
---|
682 | 'wfin_top','',3,wfin(:,:)) |
---|
683 | CALL WRITEDIAGFI(ngrid,'alpha_hmons', & |
---|
684 | 'alpha_hmons','',2,alpha_hmons) |
---|
685 | CALL WRITEDIAGFI(ngrid,'zt_top', & |
---|
686 | 'zt_top','',3,t_top(:,:)) |
---|
687 | CALL WRITEDIAGFI(ngrid,'ztemp', & |
---|
688 | 'envtemp','',3,zt(:,:)) |
---|
689 | CALL WRITEDIAGFI(ngrid,'zt_env', & |
---|
690 | 't_env','',3,t_env(:,:)) |
---|
691 | CALL WRITEDIAGFI(ngrid,'zdtvert_top', & |
---|
692 | 'zdtvert_top','',3,zdtvert(:,:)) |
---|
693 | CALL WRITEDIAGFI(ngrid,'newzt_top', & |
---|
694 | 'newzt_top','',3,newzt(:,:)) |
---|
695 | CALL WRITEDIAGFI(ngrid,'deltahr_top', & |
---|
696 | 'deltahr_top','',3,deltahr(:,:)) |
---|
697 | CALL WRITEDIAGFI(ngrid,'rhobarz', & |
---|
698 | 'rhobarz','',3,rhobarz(:,:)) |
---|
699 | CALL WRITEDIAGFI(ngrid,'zlwmax', & |
---|
700 | 'zlwmax','',2,zlaywmax(:)) |
---|
701 | CALL WRITEDIAGFI(ngrid,'coefdetrain', & |
---|
702 | 'coefdetrain','',3,coefdetrain(:,:)) |
---|
703 | CALL WRITEDIAGFI(ngrid,'detrain_rate', & |
---|
704 | 'detrain_rate', & |
---|
705 | '',3,detrain_rate(:,:)) |
---|
706 | |
---|
707 | ENDIF |
---|
708 | |
---|
709 | END SUBROUTINE topmons |
---|
710 | |
---|
711 | !======================================================================= |
---|
712 | |
---|
713 | ! ********************************************************************** |
---|
714 | ! Subroutine to determine both temperature profiles above and near |
---|
715 | ! a sub-grid scale mountain |
---|
716 | !*********************************************************************** |
---|
717 | SUBROUTINE t_topmons(nlayer,summit,hsummit,hmons,zt,zlay,t_top,dt_top,t_env,lsummit,lmons) |
---|
718 | |
---|
719 | USE comcstfi_h, only: g,cpp |
---|
720 | |
---|
721 | IMPLICIT NONE |
---|
722 | |
---|
723 | !-------------------------------------------------------- |
---|
724 | ! Input Variables |
---|
725 | !-------------------------------------------------------- |
---|
726 | integer,intent(in) :: nlayer |
---|
727 | real,intent(in) :: summit ! Height of the mountain summit |
---|
728 | real,intent(in) :: hmons ! Height of the mountain slope |
---|
729 | real,intent(in) :: hsummit ! Height of the mountain summit from the GCM surface |
---|
730 | real,intent(in) :: zt(nlayer) ! large scale temperature profile |
---|
731 | real,intent(in) :: zlay(nlayer) ! height at the middle of each layer |
---|
732 | !-------------------------------------------------------- |
---|
733 | ! Output Variables |
---|
734 | !-------------------------------------------------------- |
---|
735 | real,intent(out) :: t_top(nlayer) ! temperature above the mountain |
---|
736 | real,intent(out) :: dt_top ! temperature difference between the slope and |
---|
737 | ! the environment |
---|
738 | real,intent(out) :: t_env(nlayer) ! temperature of the environment next to the mountain |
---|
739 | integer,intent(out) :: lsummit ! layer reached by the summit height in the GCM |
---|
740 | integer,intent(out) :: lmons ! layer of the real temperature to be considered near to the mountain top |
---|
741 | |
---|
742 | !-------------------------------------------------------- |
---|
743 | ! Local Variables |
---|
744 | !-------------------------------------------------------- |
---|
745 | integer l,ll |
---|
746 | real zlmons |
---|
747 | |
---|
748 | ! ********************************************************************** |
---|
749 | t_env(:)=zt(:) |
---|
750 | t_top(:)=zt(:) |
---|
751 | dt_top=0. |
---|
752 | lsummit=1 |
---|
753 | lmons=1 |
---|
754 | |
---|
755 | !! Layer where the dust is injected |
---|
756 | do while (hsummit .ge. zlay(lsummit)) |
---|
757 | !! highest layer reached by the mountain (we don't interpolate and define zlay(lsummit) as the top of the mountain) |
---|
758 | lsummit=lsummit+1 |
---|
759 | enddo |
---|
760 | !! temperature above the mountain is set to the surface temperature |
---|
761 | t_top(lsummit)=zt(1) |
---|
762 | do l=lsummit,nlayer-1 |
---|
763 | !! temperature above the mountain follows the adiabatic |
---|
764 | t_top(l+1)=t_top(l)-(zlay(l+1)-zlay(l))*g/cpp |
---|
765 | enddo |
---|
766 | |
---|
767 | !! Layer where to take the temperature above the mountain |
---|
768 | do while (hmons .ge. zlay(lmons)) |
---|
769 | lmons=lmons+1 |
---|
770 | enddo |
---|
771 | !! Real temperature of the environment near to the mountain is t_env(lsummit)=zt(lmons) |
---|
772 | !! (More simple and makes no real difference is to impose t_env(:)=zt(:)) |
---|
773 | if (lmons.gt.lsummit) then |
---|
774 | t_env(lsummit)=zt(lmons) |
---|
775 | zlmons=zlay(lmons) |
---|
776 | do l=0,nlayer-lsummit-1!2 |
---|
777 | zlmons=zlmons+(zlay(lsummit+l+1)-zlay(lsummit+l)) |
---|
778 | ll=lmons |
---|
779 | do while (zlmons .ge. zlay(ll) .and. ll.lt.nlayer ) |
---|
780 | ll=ll+1 |
---|
781 | enddo |
---|
782 | ll=ll-1 |
---|
783 | !! Interpolation above lsummit |
---|
784 | t_env(lsummit+l+1)= zt(ll) + ((zlmons-zlay(ll))*(zt(ll+1)-zt(ll))/(zlay(ll+1)-zlay(ll))) |
---|
785 | enddo |
---|
786 | t_env(nlayer)=zt(nlayer) |
---|
787 | endif ! (lmons.gt.lsummit) |
---|
788 | do l=1,nlayer |
---|
789 | !! temperature above the mountain follows the adiabatic up to reach the environment temperature |
---|
790 | if (t_top(l).le.t_env(l)) then |
---|
791 | t_top(l)=t_env(l) |
---|
792 | endif |
---|
793 | enddo |
---|
794 | !! Temperature delta created at the top of the mountain |
---|
795 | dt_top = t_top(lsummit)-t_env(lsummit) |
---|
796 | |
---|
797 | END SUBROUTINE t_topmons |
---|
798 | |
---|
799 | !======================================================================= |
---|
800 | ! ********************************************************************** |
---|
801 | ! Subroutine to determine the vertical transport with |
---|
802 | ! a Van Leer advection scheme (copied from the sedimentation scheme --> see vlz_fi.F) |
---|
803 | !*********************************************************************** |
---|
804 | SUBROUTINE van_leer(nlay,q,pente_max,lwmax,masse,w,wq,masse_pbl,dqm_pbl,alpha_hmons,qbar_pbl) |
---|
805 | |
---|
806 | IMPLICIT NONE |
---|
807 | |
---|
808 | !-------------------------------------------------------- |
---|
809 | ! Input/Output Variables |
---|
810 | !-------------------------------------------------------- |
---|
811 | INTEGER,INTENT(IN) :: nlay ! number of atmospheric layers |
---|
812 | REAL,INTENT(IN) :: masse(nlay) ! mass of the layer Dp/g |
---|
813 | REAL,INTENT(IN) :: pente_max != 2 conseillee |
---|
814 | INTEGER,INTENT(IN) :: lwmax ! layer of dust injection above the mountain, where the vertical velocity is maximum |
---|
815 | REAL,INTENT(INOUT) :: w(nlay) ! atmospheric mass "transferred" at each timestep (kg.m-2) |
---|
816 | REAL,INTENT(INOUT) :: wq(nlay+1) |
---|
817 | REAL,INTENT(INOUT) :: q(nlay) ! mixing ratio (kg/kg) |
---|
818 | REAL,INTENT(IN) :: masse_pbl |
---|
819 | REAL,INTENT(IN) :: dqm_pbl |
---|
820 | REAL,INTENT(IN) :: alpha_hmons |
---|
821 | REAL,INTENT(IN) :: qbar_pbl |
---|
822 | |
---|
823 | !-------------------------------------------------------- |
---|
824 | ! Local Variables |
---|
825 | !-------------------------------------------------------- |
---|
826 | |
---|
827 | INTEGER i,l,j,ii |
---|
828 | REAL dzq(nlay),dzqw(nlay),adzqw(nlay),dzqmax |
---|
829 | REAL sigw, Mtot, MQtot |
---|
830 | INTEGER m |
---|
831 | |
---|
832 | ! ********************************************************************** |
---|
833 | ! Mixing ratio vertical gradient at the levels |
---|
834 | ! ********************************************************************** |
---|
835 | do l=lwmax+1,nlay !l=2,nlay |
---|
836 | dzqw(l)=q(l-1)-q(l) |
---|
837 | adzqw(l)=abs(dzqw(l)) |
---|
838 | enddo |
---|
839 | |
---|
840 | do l=lwmax+1,nlay-1 !l=2,nlay-1 |
---|
841 | if(dzqw(l)*dzqw(l+1).gt.0.) then |
---|
842 | dzq(l)=0.5*(dzqw(l)+dzqw(l+1)) |
---|
843 | else |
---|
844 | dzq(l)=0. |
---|
845 | endif |
---|
846 | dzqmax=pente_max*min(adzqw(l),adzqw(l+1)) |
---|
847 | dzq(l)=sign(min(abs(dzq(l)),dzqmax),dzq(l)) |
---|
848 | enddo |
---|
849 | |
---|
850 | dzq(lwmax)=0. !! dzq(1)=0. => Attention: in the case of van leer only above the mountain it is very important to initialize the transport from lwmax !! |
---|
851 | dzq(nlay)=0. |
---|
852 | ! ********************************************************************** |
---|
853 | ! Vertical advection |
---|
854 | ! ********************************************************************** |
---|
855 | |
---|
856 | !! No flux at the model top: |
---|
857 | wq(nlay+1)=0. |
---|
858 | |
---|
859 | !! Mass flux injected at lwmax |
---|
860 | wq(lwmax) = -dqm_pbl/alpha_hmons ! dqm_pbl = alpha_hmons x wup(wmax) x rho(lwmax) x ptimestep x qbar_pbl |
---|
861 | ! /alpha_hmons because the variable considered is mr_topdust |
---|
862 | do l = lwmax,nlay-1 |
---|
863 | |
---|
864 | ! 1) Compute wq where w < 0 (up) (UPWARD TRANSPORT) |
---|
865 | ! =============================== |
---|
866 | |
---|
867 | if (w(l+1).le.0) then |
---|
868 | ! Regular scheme (transfered mass < 1 layer) |
---|
869 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
870 | if(-w(l+1).le.masse(l))then |
---|
871 | sigw=w(l+1)/masse(l) |
---|
872 | wq(l+1)=w(l+1)*(q(l)-0.5*(1.+sigw)*dzq(l)) |
---|
873 | !!------------------------------------------------------- |
---|
874 | ! The following part should not be needed in the |
---|
875 | ! case of an integration over subtimesteps |
---|
876 | !!------------------------------------------------------- |
---|
877 | !! Extended scheme (transfered mass > 1 layer) |
---|
878 | !! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
879 | !! else |
---|
880 | !! m = l-1 |
---|
881 | !! Mtot = masse(m+1) |
---|
882 | !! MQtot = masse(m+1)*q(m+1) |
---|
883 | !! if (m.lt.lwmax)goto 77!if (m.le.0)goto 77 |
---|
884 | !! do while(-w(l+1).gt.(Mtot+masse(m))) |
---|
885 | !! m=m-1 |
---|
886 | !! Mtot = Mtot + masse(m+1) |
---|
887 | !! MQtot = MQtot + masse(m+1)*q(m+1) |
---|
888 | !! if (m.lt.lwmax)goto 77!if (m.le.0)goto 77 |
---|
889 | !! end do |
---|
890 | !! 77 continue |
---|
891 | !! |
---|
892 | !! if (m.gt.lwmax) then!if (m.gt.0) then |
---|
893 | !! sigw=(w(l+1)+Mtot)/masse(m) |
---|
894 | !! wq(l+1)= -(MQtot + (-w(l+1)-Mtot)* & |
---|
895 | !! (q(m)-0.5*(1.+sigw)*dzq(m))) |
---|
896 | !! else if ((-w(l+1).gt.(Mtot))) then |
---|
897 | !! Mtot=Mtot+masse_pbl!*alpha_hmons |
---|
898 | !! MQtot=MQtot+masse_pbl*qbar_pbl!*alpha_hmons |
---|
899 | !! sigw=(w(l+1)+Mtot)/masse(m) |
---|
900 | !! wq(l+1)= -(MQtot + (-w(l+1)-Mtot)*qbar_pbl) |
---|
901 | !! else |
---|
902 | !! w(l+1) = -Mtot |
---|
903 | !! wq(l+1) = -MQtot |
---|
904 | !! end if |
---|
905 | !! |
---|
906 | endif |
---|
907 | |
---|
908 | ! 2) Compute wq where w > 0 (down) (DOWNWARD TRANSPORT) |
---|
909 | ! =============================== |
---|
910 | |
---|
911 | else if (w(l).gt.0.) then |
---|
912 | ! Regular scheme (transfered mass < 1 layer) |
---|
913 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
914 | if(w(l).le.masse(l))then |
---|
915 | sigw=w(l)/masse(l) |
---|
916 | wq(l)=w(l)*(q(l)+0.5*(1.-sigw)*dzq(l)) |
---|
917 | !!------------------------------------------------------- |
---|
918 | ! The following part should not be needed in the |
---|
919 | ! case of an integration over subtimesteps |
---|
920 | !!------------------------------------------------------- |
---|
921 | !! Extended scheme (transfered mass > 1 layer) |
---|
922 | !! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
923 | !! else |
---|
924 | !! m=l |
---|
925 | !! Mtot = masse(m) |
---|
926 | !! MQtot = masse(m)*q(m) |
---|
927 | !! if(m.ge.nlay)goto 88 |
---|
928 | !! do while(w(l).gt.(Mtot+masse(m+1))) |
---|
929 | !! m=m+1 |
---|
930 | !! Mtot = Mtot + masse(m) |
---|
931 | !! MQtot = MQtot + masse(m)*q(m) |
---|
932 | !! if(m.ge.nlay)goto 88 |
---|
933 | !! end do |
---|
934 | !! 88 continue |
---|
935 | !! if (m.lt.nlay) then |
---|
936 | !! sigw=(w(l)-Mtot)/masse(m+1) |
---|
937 | !! wq(l)=(MQtot + (w(l)-Mtot)* & |
---|
938 | !! (q(m+1)+0.5*(1.-sigw)*dzq(m+1)) ) |
---|
939 | !! else |
---|
940 | !! w(l) = Mtot |
---|
941 | !! wq(l) = MQtot |
---|
942 | !! end if |
---|
943 | end if ! (w(l).le.masse(l)) |
---|
944 | |
---|
945 | end if ! (w(l+1).le.0) |
---|
946 | |
---|
947 | enddo ! l = lwmax+1,nlay-1 |
---|
948 | |
---|
949 | ! ********************************************************************** |
---|
950 | ! Mixing ratio update after the vertical transport |
---|
951 | ! ********************************************************************** |
---|
952 | do l = lwmax,nlay-1 !l = 1,nlay-1 ! loop different than when w>0 |
---|
953 | |
---|
954 | if ( (wq(l+1)-wq(l)) .lt. -(masse(l)*q(l)) .and. (abs(wq(l+1)).gt. 0.) .and. (abs(wq(l)).gt. 0.) ) then |
---|
955 | wq(l+1) = wq(l)-masse(l)*q(l) |
---|
956 | end if |
---|
957 | |
---|
958 | q(l)=q(l) + (wq(l+1)-wq(l))/masse(l) |
---|
959 | |
---|
960 | enddo |
---|
961 | |
---|
962 | |
---|
963 | end subroutine van_leer |
---|
964 | |
---|
965 | !******************************************************************************** |
---|
966 | ! initialization module variables |
---|
967 | subroutine ini_topmons_mod(ngrid,nlayer) |
---|
968 | |
---|
969 | implicit none |
---|
970 | |
---|
971 | integer, intent(in) :: ngrid |
---|
972 | integer, intent(in) :: nlayer |
---|
973 | |
---|
974 | allocate(alpha_hmons(ngrid)) |
---|
975 | |
---|
976 | end subroutine ini_topmons_mod |
---|
977 | |
---|
978 | subroutine end_topmons_mod |
---|
979 | |
---|
980 | implicit none |
---|
981 | |
---|
982 | if (allocated(alpha_hmons)) deallocate(alpha_hmons) |
---|
983 | |
---|
984 | end subroutine end_topmons_mod |
---|
985 | |
---|
986 | END MODULE topmons_mod |
---|