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

Last change on this file since 2091 was 2091, checked in by mvals, 6 years ago

Mars GCM:

  • follow-up of the last change in rocketduststorm_mod.F90: corrections in vector definitions

MV

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