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

Last change on this file since 2627 was 2616, checked in by romain.vande, 3 years ago

LMDZ_MARS RV : Open_MP;
Put all the "save" variables as "!$OMP THREADPRIVATE" in phymars.
The code can now be tested, see README for more info

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