source: trunk/LMDZ.MARS/libf/phymars/topmons_mod.F90 @ 2415

Last change on this file since 2415 was 2415, checked in by emillour, 4 years ago

Mars GCM:
Code cleanup: remove "tauref" and replace it with two distinct variables,
"tau_pref_gcm" and "tau_pref_scenario" which are repectively the visible
dust opacity column (at 610 Pa) in the GCM and prescribed by a dust scenario.
EM

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