source: trunk/LMDZ.MARS/libf/phymars/rocketduststorm_mod.F90 @ 2800

Last change on this file since 2800 was 2685, checked in by emillour, 2 years ago

Mars GCM:
Add possibility to output either upward or downward SW flux at the surface
and top of atmosphere from physiq. Required adding some output arguments
to callradite.
EM

File size: 31.4 KB
Line 
1      MODULE rocketduststorm_mod
2
3      IMPLICIT NONE
4
5      REAL, SAVE, ALLOCATABLE :: dustliftday(:) ! dust lifting rate (s-1)
6
7!$OMP THREADPRIVATE(dustliftday)
8     
9      CONTAINS
10
11!=======================================================================
12! ROCKET DUST STORM - vertical transport and detrainment
13!=======================================================================
14! calculation of the vertical flux
15! call of van_leer : Van Leer transport scheme of the dust tracers
16! detrainement of stormdust in the background dust
17! -----------------------------------------------------------------------
18! Authors: M. Vals; C. Wang; F. Forget; T. Bertrand
19! Institution: Laboratoire de Meteorologie Dynamique (LMD) Paris, France
20! -----------------------------------------------------------------------
21
22      SUBROUTINE rocketduststorm(ngrid,nlayer,nq,ptime,ptimestep,      &
23                                 pq,pdqfi,pt,pdtfi,pplev,pplay,pzlev,  &
24                                 pzlay,pdtsw,pdtlw,                    &
25!             input for radiative transfer
26                                 clearatm,icount,zday,zls,             &
27                                 tsurf,igout,totstormfract,            &
28                                 tauscaling,dust_rad_adjust,           &
29                                 IRtoVIScoef,                          &
30!             input sub-grid scale cloud
31                                 clearsky,totcloudfrac,                &
32!             input sub-grid scale topography
33                                 nohmons,                              &
34!             output
35                                 pdqrds,wrad,dsodust,dsords,dsotop,    &
36                                 tau_pref_scenario,tau_pref_gcm)
37
38      USE tracer_mod, only: igcm_stormdust_mass,igcm_stormdust_number &
39                            ,igcm_dust_mass,igcm_dust_number          &
40                            ,rho_dust
41      USE comcstfi_h, only: r,g,cpp,rcp
42      USE dimradmars_mod, only: albedo,naerkind
43      USE comsaison_h, only: dist_sol,mu0,fract
44      USE surfdat_h, only: emis,co2ice,zmea, zstd, zsig, hmons
45      USE callradite_mod, only: callradite
46      IMPLICIT NONE
47
48      include "callkeys.h"
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) :: pdqfi(ngrid,nlayer,nq)! tendancy field pq
62      REAL, INTENT(IN) :: pt(ngrid,nlayer)      ! temperature at mid-layer (K)
63      REAL, INTENT(IN) :: pdtfi(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      LOGICAL, INTENT(IN) :: clearatm
75      INTEGER, INTENT(INOUT) :: icount
76      REAL, INTENT(IN) :: zday
77      REAL, INTENT(IN) :: zls
78      REAL, INTENT(IN) :: tsurf(ngrid)
79      INTEGER, INTENT(IN) :: igout
80      REAL, INTENT(IN) :: totstormfract(ngrid)
81      REAL, INTENT(INOUT) :: tauscaling(ngrid)
82      REAL,INTENT(INOUT) :: dust_rad_adjust(ngrid)                             
83      REAL,INTENT(INOUT) :: IRtoVIScoef(ngrid) ! NB: not modified by this call to callradite,
84                                               ! the OUT is just here because callradite needs it
85     
86!     subgrid scale water ice clouds
87      logical, intent(in) :: clearsky
88      real, intent(in) :: totcloudfrac(ngrid)
89
90!     subgrid scale topography
91      LOGICAL, INTENT(IN) :: nohmons
92 
93!--------------------------------------------------------
94! Output Variables
95!--------------------------------------------------------
96
97      REAL, INTENT(OUT) :: pdqrds(ngrid,nlayer,nq) ! tendancy field for dust when detraining
98      REAL, INTENT(OUT) :: wrad(ngrid,nlayer+1)   ! vertical speed within the rocket dust storm
99      REAL, INTENT(OUT) :: dsodust(ngrid,nlayer) ! density scaled opacity of env. dust
100      REAL, INTENT(OUT) :: dsords(ngrid,nlayer) ! density scaled opacity of storm dust
101      REAL, INTENT(OUT) :: dsotop(ngrid,nlayer) ! density scaled opacity of topmons dust
102      REAL,INTENT(OUT) :: tau_pref_scenario(ngrid) ! prescribed dust column
103                                               ! visible opacity at odpref
104      REAL,INTENT(OUT) :: tau_pref_gcm(ngrid) ! dust column visible opacity at
105                                              ! odpref in the GCM
106!--------------------------------------------------------
107! Local variables
108!--------------------------------------------------------
109      INTEGER l,ig,iq,ll
110!     local variables from callradite.F       
111      REAL zdtlw1(ngrid,nlayer)    ! (K/s) storm
112      REAL zdtsw1(ngrid,nlayer)    ! (K/s) storm
113      REAL zt(ngrid,nlayer)        ! actual temperature at mid-layer (K)
114      REAL zdtvert(ngrid,nlayer)   ! dT/dz , lapse rate
115      REAL ztlev(ngrid,nlayer)     ! temperature at intermediate levels l+1/2 without last level
116
117      REAL zdtlw1_lev(nlayer),zdtsw1_lev(nlayer) ! rad. heating rate at intermediate levels l+1/2 for stormdust
118      REAL zdtlw_lev(nlayer),zdtsw_lev(nlayer)   ! rad. heating rate at intermediate levels l+1/2 for background dust
119
120      REAL zq_stormdust_mass(ngrid,nlayer) ! intermediate tracer stormdust mass
121      REAL zq_stormdust_number(ngrid,nlayer) ! intermediate tracer stormdust number
122      REAL zq_dust_mass(ngrid,nlayer)           ! intermediate tracer dust mass
123      REAL zq_dust_number(ngrid,nlayer)         ! intermediate tracer dust number
124
125      REAL mr_stormdust_mass(ngrid,nlayer) ! intermediate mixing ratio to calculate van leer transport with the "real" concentration (stormdust mass)
126      REAL mr_stormdust_number(ngrid,nlayer) ! intermediate mixing ratio to calculate van leer transport with the "real" concentration (stormdust number)
127      REAL mr_dust_mass(ngrid,nlayer) ! intermediate mixing ratio to calculate van leer transport with the "real" concentration (dust mass)
128      REAL mr_dust_number(ngrid,nlayer) ! intermediate mixing ratio to calculate van leer transport with the "real" concentration (sdust number)
129                   
130      REAL dqvl_stormdust_mass(ngrid,nlayer)    ! tendancy of vertical transport (stormdust mass)
131      REAL dqvl_stormdust_number(ngrid,nlayer)  ! tendancy of vertical transport (stormdust number)
132      REAL dqvl_dust_mass(ngrid,nlayer)    ! tendancy of vertical transport (dust mass)
133      REAL dqvl_dust_number(ngrid,nlayer)  ! tendancy of vertical transport (dust number)
134      REAL dqdet_stormdust_mass(ngrid,nlayer)   ! tendancy of detrainement (stormdust mass)
135      REAL dqdet_stormdust_number(ngrid,nlayer) ! tendancy of detrainement (stormdust number)
136
137      REAL masse_col(nlayer)     ! mass of atmosphere (kg/m2)
138      REAL zq(ngrid,nlayer,nq)   ! updated tracers
139     
140      REAL w(ngrid,nlayer)          ! air mass flux (calculated with the vertical wind velocity profile) used as input in Van Leer (kgair/m2)
141      REAL wqmass(ngrid,nlayer+1)   ! tracer (dust_mass) mass flux in Van Leer (kg/m2)
142      REAL wqnumber(ngrid,nlayer+1) ! tracer (dust_number) mass flux in Van Leer (kg/m2)
143
144      LOGICAL storm(ngrid)    ! true when there is a dust storm (if the opacity is high): trigger the rocket dust storm scheme
145      REAL detrain(ngrid,nlayer)  ! coefficient for detrainment : % of stormdust detrained
146      INTEGER scheme(ngrid)   ! triggered scheme
147           
148      REAL,PARAMETER:: coefmin =0.025 ! 0<coefmin<1 Minimum fraction of stormdust detrained 
149      REAL,PARAMETER:: wmin =0.25 ! stormdust detrainment if wrad < wmin 
150      REAL,PARAMETER:: wmax =10.   ! maximum vertical velocity of the rocket dust storms (m/s)
151
152!     subtimestep
153      INTEGER tsub
154      INTEGER nsubtimestep    !number of subtimestep when calling van_leer
155      REAL subtimestep        !ptimestep/nsubtimestep
156      REAL dtmax              !considered time needed for dust to cross one layer
157      REAL,PARAMETER:: secu=3.!3.      !coefficient on wspeed to avoid dust crossing many layers during subtimestep
158
159!     diagnostics
160      REAL lapserate(ngrid,nlayer)
161      REAL deltahr(ngrid,nlayer+1)
162   
163      LOGICAL,SAVE :: firstcall=.true.
164     
165!$OMP THREADPRIVATE(firstcall)
166
167!     variables for the radiative transfer
168      REAL  fluxsurf_lw1(ngrid)
169      REAL  fluxsurf_dn_sw1(ngrid,2),fluxsurf_up_sw1(ngrid,2)
170      REAL  fluxtop_lw1(ngrid)
171      REAL  fluxtop_dn_sw1(ngrid,2),fluxtop_up_sw1(ngrid,2)
172      REAL  tau(ngrid,naerkind)
173      REAL  aerosol(ngrid,nlayer,naerkind)
174      REAL  taucloudtes(ngrid)
175      REAL  rdust(ngrid,nlayer)
176      REAL  rstormdust(ngrid,nlayer)
177      REAL  rtopdust(ngrid,nlayer)
178      REAL  rice(ngrid,nlayer)
179      REAL  nuice(ngrid,nlayer)
180      DOUBLE PRECISION  riceco2(ngrid,nlayer)
181      REAL  nuiceco2(ngrid,nlayer)
182
183      ! **********************************************************************
184      ! **********************************************************************
185      ! Rocket dust storm parametrization to reproduce the detached dust layers
186      ! during the dust storm season:
187      !     The radiative warming due to the presence of storm dust is
188      !     balanced by the adiabatic cooling. The tracer "storm dust" 
189      !     is transported by the upward/downward flow.
190      ! **********************************************************************
191      ! **********************************************************************     
192      !!    1. Radiative transfer in storm dust
193      !!    2. Compute vertical velocity for storm dust
194      !!      case 1 storm = false: nothing to do
195      !!      case 2 rocket dust storm (storm=true)
196      !!    3. Vertical transport (Van Leer)
197      !!    4. Detrainment
198      ! **********************************************************************
199      ! **********************************************************************
200
201     
202      ! **********************************************************************
203      ! Initializations
204      ! **********************************************************************
205      storm(:)=.false.
206      pdqrds(:,:,:) = 0.
207      mr_dust_mass(:,:)=0.
208      mr_dust_number(:,:)=0.
209      mr_stormdust_mass(:,:)=0.
210      mr_stormdust_number(:,:)=0.
211      dqvl_dust_mass(:,:)=0.
212      dqvl_dust_number(:,:)=0.
213      dqvl_stormdust_mass(:,:)=0.
214      dqvl_stormdust_number(:,:)=0.
215      dqdet_stormdust_mass(:,:)=0.
216      dqdet_stormdust_number(:,:)=0.
217      wrad(:,:)=0.
218      w(:,:)=0.
219      wqmass(:,:)=0.
220      wqnumber(:,:)=0.
221      zdtvert(:,:)=0.
222      lapserate(:,:)=0.
223      deltahr(:,:)=0.
224      scheme(:)=0
225      detrain(:,:)=1.
226
227      !! no update for the stormdust tracer and temperature fields
228      !! because previous callradite was for background dust
229      zq(1:ngrid,1:nlayer,1:nq)=pq(1:ngrid,1:nlayer,1:nq)
230      zt(1:ngrid,1:nlayer)=pt(1:ngrid,1:nlayer)
231
232
233      zq_dust_mass(1:ngrid,1:nlayer)=zq(1:ngrid,1:nlayer,igcm_dust_mass)
234      zq_dust_number(1:ngrid,1:nlayer)=zq(1:ngrid,1:nlayer,igcm_dust_number)
235      zq_stormdust_mass(1:ngrid,1:nlayer)=zq(1:ngrid,1:nlayer,igcm_stormdust_mass)
236      zq_stormdust_number(1:ngrid,1:nlayer)=zq(1:ngrid,1:nlayer,igcm_stormdust_number)
237
238      ! *********************************************************************
239      ! 0. Check if there is a storm
240      ! *********************************************************************
241      DO ig=1,ngrid
242        storm(ig)=.false.
243        DO l=1,nlayer
244          IF (zq(ig,l,igcm_stormdust_mass) &
245          .gt. zq(ig,l,igcm_dust_mass)*(1.E-4)) THEN
246            storm(ig)=.true.
247            EXIT
248          ENDIF
249        ENDDO     
250      ENDDO
251     
252      ! *********************************************************************
253      ! 1. Call the second radiative transfer for stormdust, obtain the extra heating
254      ! *********************************************************************
255      CALL callradite(icount,ngrid,nlayer,nq,zday,zls,pq,albedo,          &
256                 emis,mu0,pplev,pplay,pt,tsurf,fract,dist_sol,igout,      &
257                 zdtlw1,zdtsw1,fluxsurf_lw1,fluxsurf_dn_sw1,fluxsurf_up_sw1, &
258                 fluxtop_lw1,fluxtop_dn_sw1,fluxtop_up_sw1, &
259                 tau_pref_scenario,tau_pref_gcm, &
260                 tau,aerosol,dsodust,tauscaling,dust_rad_adjust,IRtoVIScoef, &
261                 taucloudtes,rdust,rice,nuice,riceco2,nuiceco2,co2ice,rstormdust,rtopdust, &
262                 totstormfract,clearatm,dsords,dsotop,nohmons,&
263                 clearsky,totcloudfrac)
264
265      ! **********************************************************************
266      ! 2. Compute vertical velocity for storm dust
267      ! **********************************************************************
268        !! **********************************************************************
269        !! 2.1 Nothing to do when no storm
270        !!             no storm
271        DO ig=1,ngrid       
272          IF (.NOT.(storm(ig))) then
273            scheme(ig)=1
274            cycle
275          ENDIF ! IF (storm(ig))
276        ENDDO ! DO ig=1,ngrid                 
277       
278        !! **********************************************************************
279        !! 2.2 Calculation of the extra heating : computing heating rates
280        !! gradient at boundaries of each layer, start from surface
281        DO ig=1,ngrid
282          IF (storm(ig)) THEN
283
284            scheme(ig)=2
285     
286            !! computing heating rates gradient at boundraies of each layer
287            !! start from surface
288            zdtlw1_lev(1)=0.
289            zdtsw1_lev(1)=0.
290            zdtlw_lev(1)=0.
291            zdtsw_lev(1)=0.
292            ztlev(ig,1)=zt(ig,1)
293
294            DO l=1,nlayer-1
295              !! Calculation for the dust storm fraction
296              zdtlw1_lev(l+1)=(zdtlw1(ig,l)*(pzlay(ig,l+1)-pzlev(ig,l+1))+ &
297                           zdtlw1(ig,l+1)*(pzlev(ig,l+1)-pzlay(ig,l)))  /  &
298                              (pzlay(ig,l+1)-pzlay(ig,l))
299           
300              zdtsw1_lev(l+1)=(zdtsw1(ig,l)*(pzlay(ig,l+1)-pzlev(ig,l+1))+ &
301                           zdtsw1(ig,l+1)*(pzlev(ig,l+1)-pzlay(ig,l)))  /  &
302                              (pzlay(ig,l+1)-pzlay(ig,l))
303              !! Calculation for the background dust fraction
304              zdtlw_lev(l+1)=(pdtlw(ig,l)*(pzlay(ig,l+1)-pzlev(ig,l+1))+   &
305                           pdtlw(ig,l+1)*(pzlev(ig,l+1)-pzlay(ig,l)))  /   &
306                              (pzlay(ig,l+1)-pzlay(ig,l))
307           
308              zdtsw_lev(l+1)=(pdtsw(ig,l)*(pzlay(ig,l+1)-pzlev(ig,l+1))+   &
309                           pdtsw(ig,l+1)*(pzlev(ig,l+1)-pzlay(ig,l)))  /   &
310                              (pzlay(ig,l+1)-pzlay(ig,l))
311           
312              ztlev(ig,l+1)=(zt(ig,l)*(pzlay(ig,l+1)-pzlev(ig,l+1))+          &
313                           zt(ig,l+1)*(pzlev(ig,l+1)-pzlay(ig,l)))  /      &
314                              (pzlay(ig,l+1)-pzlay(ig,l))
315            ENDDO ! DO l=1,nlayer-1
316
317            !! This is the env. lapse rate
318            zdtvert(ig,1)=0.
319            DO l=1,nlayer-1
320              zdtvert(ig,l+1)=(ztlev(ig,l+1)-ztlev(ig,l))/(pzlay(ig,l+1)-pzlay(ig,l))
321              lapserate(ig,l+1)=zdtvert(ig,l+1)
322            ENDDO
323
324            !! **********************************************************************
325            !! 2.3 Calculation of the vertical velocity : extra heating
326            !! balanced by adiabatic cooling
327           
328            DO l=1,nlayer
329              deltahr(ig,l)=(zdtlw1_lev(l)+zdtsw1_lev(l))  &
330                                          -(zdtlw_lev(l)+zdtsw_lev(l))
331              wrad(ig,l)=-deltahr(ig,l)/(g/cpp+   &
332                                         max(zdtvert(ig,l),-0.99*g/cpp))       
333              !! Limit vertical wind in case of lapse rate close to adiabatic
334              wrad(ig,l)=max(wrad(ig,l),-wmax)
335              wrad(ig,l)=min(wrad(ig,l),wmax)
336            ENDDO
337
338          ENDIF ! IF (storm(ig))
339        ENDDO ! DO ig=1,ngrid
340
341      ! **********************************************************************
342      ! 3. Vertical transport
343      ! **********************************************************************
344        !! **********************************************************************
345        !! 3.1 Transport of background dust + storm dust (concentrated)
346        DO ig=1,ngrid
347          IF (storm(ig)) THEN
348            DO l=1,nlayer
349              mr_dust_mass(ig,l) = zq_dust_mass(ig,l)
350              mr_dust_number(ig,l) = zq_dust_number(ig,l)
351              mr_stormdust_mass(ig,l) = zq_dust_mass(ig,l)+ &
352                                      zq_stormdust_mass(ig,l)/totstormfract(ig)
353              mr_stormdust_number(ig,l) = zq_dust_number(ig,l)+ &
354                                      zq_stormdust_number(ig,l)/totstormfract(ig)
355            ENDDO ! DO l=1,nlayer
356          ENDIF ! IF (storm(ig))
357        ENDDO ! DO ig=1,ngrid
358
359        DO ig=1,ngrid
360          IF (storm(ig)) THEN
361        !! **********************************************************************
362        !! 3.2 Compute the subtimestep to conserve the mass in the Van Leer transport
363            dtmax=ptimestep
364            DO l=2,nlayer
365              IF (wrad(ig,l).lt.0.) THEN
366                 dtmax=min(dtmax,(pzlev(ig,l)-pzlev(ig,l-1))/  &
367                                   (secu*abs(wrad(ig,l))))
368              ELSE IF (wrad(ig,l).gt.0.) then
369                 dtmax=min(dtmax,(pzlev(ig,l+1)-pzlev(ig,l))/  &
370                                   (secu*abs(wrad(ig,l))))
371              ENDIF
372            ENDDO
373            nsubtimestep= int(ptimestep/dtmax)
374            subtimestep=ptimestep/float(nsubtimestep)
375            !! Mass flux generated by wup in kg/m2
376            DO l=1,nlayer
377               w(ig,l)=wrad(ig,l)*pplev(ig,l)/(r*ztlev(ig,l)) &
378                        *subtimestep
379            ENDDO ! l=1,nlayer
380
381        !! **********************************************************************
382        !! 3.3 Vertical transport by a Van Leer scheme
383            !! Mass of atmosphere in the layer
384            DO l=1,nlayer
385               masse_col(l)=(pplev(ig,l)-pplev(ig,l+1))/g
386            ENDDO
387            !! Mass flux in kg/m2 if you are not using the subtimestep
388            !DO l=1,nlayer
389            !   w(ig,l)=wrad(ig,l)*(pplev(ig,l)/(r*ztlev(ig,l)))*ptimestep
390            !ENDDO
391            !! Loop over the subtimestep
392            DO tsub=1,nsubtimestep
393              !! Van Leer scheme
394              wqmass(ig,:)=0.
395              wqnumber(ig,:)=0.
396              CALL van_leer(nlayer,mr_stormdust_mass(ig,:),2.,   &
397                    masse_col,w(ig,:),wqmass(ig,:))
398              CALL van_leer(nlayer,mr_stormdust_number(ig,:),2.,  &
399                    masse_col,w(ig,:),wqnumber(ig,:))
400            ENDDO !tsub=...           
401             
402          ENDIF ! IF storm(ig)
403        ENDDO ! DO ig=1,ngrid 
404
405        !! **********************************************************************
406        !! 3.4 Re-calculation of the mixing ratios after vertical transport
407        DO ig=1,ngrid
408         IF (storm(ig)) THEN
409           DO l=1,nlayer
410           
411             !! General and "healthy" case
412             IF (mr_stormdust_mass(ig,l).ge.mr_dust_mass(ig,l)) THEN
413               zq_dust_mass(ig,l) = mr_dust_mass(ig,l)
414               zq_dust_number(ig,l) = mr_dust_number(ig,l)
415               zq_stormdust_mass(ig,l) = totstormfract(ig)*(mr_stormdust_mass(ig,l)-mr_dust_mass(ig,l))
416               zq_stormdust_number(ig,l) = totstormfract(ig)*(mr_stormdust_number(ig,l)-mr_dust_number(ig,l))
417             !! Particular case: there is less than initial dust mixing ratio after the vertical transport
418             ELSE
419               zq_dust_mass(ig,l) = (1.-totstormfract(ig))*mr_dust_mass(ig,l)+totstormfract(ig)*mr_stormdust_mass(ig,l)
420               zq_dust_number(ig,l) = (1.-totstormfract(ig))*mr_dust_number(ig,l)+totstormfract(ig)*mr_stormdust_number(ig,l)
421               zq_stormdust_mass(ig,l) = 0.
422               zq_stormdust_number(ig,l) = 0.
423             ENDIF
424             
425           ENDDO ! DO l=1,nlayer           
426         ENDIF ! IF storm(ig)
427        ENDDO ! DO ig=1,ngrid
428
429        !! **********************************************************************
430        !! 3.5 Calculation of the tendencies of the vertical transport
431        DO ig=1,ngrid
432         IF (storm(ig)) THEN
433           DO l=1,nlayer
434             dqvl_stormdust_mass(ig,l) = (zq_stormdust_mass(ig,l)- &
435                                   zq(ig,l,igcm_stormdust_mass)) /ptimestep
436             dqvl_stormdust_number(ig,l) = (zq_stormdust_number(ig,l)- &
437                                     zq(ig,l,igcm_stormdust_number)) /ptimestep
438             dqvl_dust_mass(ig,l) = (zq_dust_mass(ig,l)-zq(ig,l,igcm_dust_mass)) /ptimestep
439             dqvl_dust_number(ig,l) = (zq_dust_number(ig,l)-zq(ig,l,igcm_dust_number)) /ptimestep
440           ENDDO
441         ENDIF ! IF storm(ig)
442        ENDDO ! DO ig=1,ngrid           
443
444      ! **********************************************************************
445      ! 4. Detrainment: stormdust is converted to background dust
446      ! **********************************************************************
447        !! **********************************************************************
448        !! 4.1 Compute the coefficient of detrainmen
449        DO ig=1,ngrid
450          DO l=1,nlayer
451            IF ((max(abs(wrad(ig,l)),abs(wrad(ig,l+1))) .lt.  &
452                          wmin) .or. (zq_dust_mass(ig,l) .gt.  &
453                                     10000.*zq_stormdust_mass(ig,l))) THEN
454               detrain(ig,l)=1.
455            ELSE IF (max(abs(wrad(ig,l)),abs(wrad(ig,l+1)))   &
456                                                       .le. wmax) THEN
457               detrain(ig,l)=coeff_detrainment*                 &
458                             (((1-coefmin)/(wmin-wmax)**2)*     &
459                             (max(abs(wrad(ig,l)),abs(wrad(ig,l+1)))-wmax)**2 &
460                              +coefmin)
461            ELSE IF (max(abs(wrad(ig,l)),abs(wrad(ig,l+1))).gt. wmax ) THEN
462               detrain(ig,l)=coefmin
463            ELSE
464               detrain(ig,l)=coefmin
465            ENDIF
466          ENDDO ! DO l=1,nlayer
467        ENDDO ! DO ig=1,ngrid
468       
469        !! **********************************************************************
470        !! 4.2 Calculation of the tendencies of the detrainment
471        DO ig=1,ngrid
472          DO l=1,nlayer
473            dqdet_stormdust_mass(ig,l)=-detrain(ig,l)*zq_stormdust_mass(ig,l) &
474                                                        /ptimestep
475            dqdet_stormdust_number(ig,l)=-detrain(ig,l)*zq_stormdust_number(ig,l) &
476                                                        /ptimestep
477          ENDDO ! DO l=1,nlayer
478        ENDDO ! DO ig=1,ngrid
479           
480      ! **********************************************************************
481      ! 5. Final tendencies: vertical transport + detrainment
482      ! **********************************************************************
483        DO ig=1,ngrid
484          DO l=1,nlayer     
485          pdqrds(ig,l,igcm_stormdust_mass)=dqdet_stormdust_mass(ig,l) &
486                                                 +dqvl_stormdust_mass(ig,l)
487          pdqrds(ig,l,igcm_stormdust_number)=dqdet_stormdust_number(ig,l) &
488                                                 +dqvl_stormdust_number(ig,l)
489          pdqrds(ig,l,igcm_dust_mass)= -dqdet_stormdust_mass(ig,l) &
490                                       +dqvl_dust_mass(ig,l)
491          pdqrds(ig,l,igcm_dust_number)= -dqdet_stormdust_number(ig,l) &
492                                       +dqvl_dust_number(ig,l)
493          ENDDO ! DO l=1,nlayer
494        ENDDO ! DO ig=1,ngrid
495
496!      ! **********************************************************************
497!      ! 6. To prevent from negative values
498!      ! **********************************************************************
499!        DO ig=1,ngrid
500!          DO l=1,nlayer
501!            IF ((pq(ig,l,igcm_stormdust_mass) &
502!               +pdqrds(ig,l,igcm_stormdust_mass)*ptimestep .le. 0.) .or. &
503!              (pq(ig,l,igcm_stormdust_number) &
504!               +pdqrds(ig,l,igcm_stormdust_number)*ptimestep .le. 0.)) THEN
505!               pdqrds(ig,l,igcm_stormdust_mass)=-pq(ig,l,igcm_stormdust_mass)/ptimestep
506!               pdqrds(ig,l,igcm_stormdust_number)=-pq(ig,l,igcm_stormdust_number)/ptimestep
507!            ENDIF
508!          ENDDO ! nlayer
509!        ENDDO ! DO ig=1,ngrid
510!
511!        DO ig=1,ngrid
512!          DO l=1,nlayer           
513!            IF ((pq(ig,l,igcm_dust_mass) &
514!               +pdqrds(ig,l,igcm_dust_mass)*ptimestep .le. 0.) .or. &
515!              (pq(ig,l,igcm_dust_number) &
516!               +pdqrds(ig,l,igcm_dust_number)*ptimestep .le. 0.)) THEN
517!               pdqrds(ig,l,igcm_dust_mass)=-pq(ig,l,igcm_dust_mass)/ptimestep
518!               pdqrds(ig,l,igcm_dust_number)=-pq(ig,l,igcm_dust_number)/ptimestep
519!            ENDIF
520!          ENDDO ! nlayer
521!        ENDDO ! DO ig=1,ngrid
522       
523!=======================================================================
524! WRITEDIAGFI
525
526        call WRITEDIAGFI(ngrid,'lapserate','lapse rate in the storm', &
527           &                        'k/m',3,lapserate)
528        call WRITEDIAGFI(ngrid,'deltahr','extra heating rates', &
529           &                        'K/s',3,deltahr)
530        call writediagfi(ngrid,'scheme','which scheme',&
531                                                   ' ',2,real(scheme))
532
533        END SUBROUTINE rocketduststorm
534
535!======================================================================= 
536! **********************************************************************
537! Subroutine to determine the vertical transport with
538! a Van Leer advection scheme (copied from the sedimentation scheme --> see vlz_fi.F)
539!***********************************************************************
540      SUBROUTINE van_leer(nlay,q,pente_max,masse,w,wq)
541
542      IMPLICIT NONE
543
544!--------------------------------------------------------
545! Input/Output Variables
546!--------------------------------------------------------
547      INTEGER,INTENT(IN) :: nlay       ! number of atmospheric layers
548      REAL,INTENT(IN) ::  masse(nlay)  ! mass of the layer Dp/g
549      REAL,INTENT(IN) :: pente_max     != 2 conseillee
550      REAL,INTENT(INOUT) :: q(nlay)    ! mixing ratio (kg/kg)
551      REAL,INTENT(INOUT) :: w(nlay)    ! atmospheric mass "transferred" at each timestep (kg.m-2)
552      REAL,INTENT(INOUT) :: wq(nlay+1)
553
554!--------------------------------------------------------
555! Local Variables
556!--------------------------------------------------------
557
558      INTEGER i,l,j,ii
559      REAL dzq(nlay),dzqw(nlay),adzqw(nlay),dzqmax
560      REAL newmasse
561      REAL sigw, Mtot, MQtot
562      INTEGER m
563
564! **********************************************************************
565!  Mixing ratio vertical gradient at the levels
566! **********************************************************************
567      do l=2,nlay
568            dzqw(l)=q(l-1)-q(l)
569            adzqw(l)=abs(dzqw(l))
570      enddo
571
572      do l=2,nlay-1
573            if(dzqw(l)*dzqw(l+1).gt.0.) then
574                dzq(l)=0.5*(dzqw(l)+dzqw(l+1))
575            else
576                dzq(l)=0.
577            endif
578            dzqmax=pente_max*min(adzqw(l),adzqw(l+1))
579            dzq(l)=sign(min(abs(dzq(l)),dzqmax),dzq(l))
580      enddo
581
582      dzq(1)=0.
583      dzq(nlay)=0.
584
585! **********************************************************************
586!  Vertical advection
587! **********************************************************************
588
589       !! No flux at the model top:
590       wq(nlay+1)=0.
591
592       !! Surface flux up:
593       if(w(1).lt.0.) wq(1)=0. ! warning : not always valid
594
595       do l = 1,nlay-1
596
597!      1) Compute wq where w < 0 (up) (UPWARD TRANSPORT)
598!      ===============================
599
600         if(w(l+1).le.0)then
601!         Regular scheme (transfered mass < 1 layer)
602!         ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
603          if(-w(l+1).le.masse(l))then
604            sigw=w(l+1)/masse(l)
605            wq(l+1)=w(l+1)*(q(l)-0.5*(1.+sigw)*dzq(l))
606!!-------------------------------------------------------
607!          The following part should not be needed in the
608!          case of an integration over subtimesteps
609!!-------------------------------------------------------
610!         Extended scheme (transfered mass > 1 layer)
611!         ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
612          else
613             m = l-1
614             Mtot = masse(m+1)
615             MQtot = masse(m+1)*q(m+1)
616             if (m.le.0)goto 77
617             do while(-w(l+1).gt.(Mtot+masse(m)))
618!             do while(-w(l+1).gt.Mtot)
619                m=m-1
620                Mtot = Mtot + masse(m+1)
621                MQtot = MQtot + masse(m+1)*q(m+1)
622                if (m.le.0)goto 77
623             end do
624 77          continue
625
626             if (m.gt.0) then
627                sigw=(w(l+1)+Mtot)/masse(m)
628                wq(l+1)= -(MQtot + (-w(l+1)-Mtot)*         &
629                (q(m)-0.5*(1.+sigw)*dzq(m))  )
630             else
631                w(l+1) = -Mtot
632                wq(l+1) = -MQtot
633             end if
634          endif ! (-w(l+1).le.masse(l))
635     
636!      2) Compute wq where w > 0 (down) (DOWNWARD TRANSPORT)     
637!      ===============================
638
639         else if(w(l).gt.0.)then
640
641!         Regular scheme (transfered mass < 1 layer)
642!         ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
643          if(w(l).le.masse(l))then
644            sigw=w(l)/masse(l)
645            wq(l)=w(l)*(q(l)+0.5*(1.-sigw)*dzq(l))           
646!!-------------------------------------------------------
647!          The following part should not be needed in the
648!          case of an integration over subtimesteps
649!!-------------------------------------------------------
650!         Extended scheme (transfered mass > 1 layer)
651!         ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
652          else
653            m=l
654            Mtot = masse(m)
655            MQtot = masse(m)*q(m)
656            if(m.ge.nlay)goto 88
657            do while(w(l).gt.(Mtot+masse(m+1)))
658                m=m+1
659                Mtot = Mtot + masse(m)
660                MQtot = MQtot + masse(m)*q(m)
661                if(m.ge.nlay)goto 88
662            end do
663 88         continue
664            if (m.lt.nlay) then
665                sigw=(w(l)-Mtot)/masse(m+1)
666                wq(l)=(MQtot + (w(l)-Mtot)* &
667                          (q(m+1)+0.5*(1.-sigw)*dzq(m+1)) )
668            else
669                w(l) = Mtot
670                wq(l) = MQtot
671            end if
672          end if
673
674         end if ! w<0 (up)
675
676       enddo ! l = 1,nlay-1
677       
678       do l = 1,nlay
679
680!         it cannot entrain more than available mass !
681          if ( (wq(l+1)-wq(l)) .lt. -(masse(l)*q(l)) ) then
682            wq(l+1) = wq(l)-masse(l)*q(l)
683          end if
684
685          q(l)=q(l) +  (wq(l+1)-wq(l))/masse(l)
686
687       enddo
688       
689      END SUBROUTINE van_leer
690
691!=======================================================================
692! Initialization of the module variables
693
694       subroutine ini_rocketduststorm_mod(ngrid)
695       
696       implicit none
697       
698       integer, intent(in) :: ngrid
699       
700       allocate(dustliftday(ngrid))
701       
702       end subroutine ini_rocketduststorm_mod
703       
704       subroutine end_rocketduststorm_mod
705       
706       implicit none
707       
708       if (allocated(dustliftday)) deallocate(dustliftday)
709
710       end subroutine end_rocketduststorm_mod       
711     
712      END MODULE rocketduststorm_mod
Note: See TracBrowser for help on using the repository browser.