source: trunk/LMDZ.MARS/libf/phymars/simpleclouds.F @ 1246

Last change on this file since 1246 was 1246, checked in by aslmd, 11 years ago

LMDZ.MARS. Made number of scatterers a free dimension not in need to be prescribe at compiling time. Instead it must be set in callphys.def. See README for further information about this commit.

File size: 6.5 KB
Line 
1      subroutine simpleclouds(ngrid,nlay,ptimestep,
2     &             pplay,pzlay,pt,pdt,
3     &             pq,pdq,pdqcloud,pdtcloud,
4     &             nq,tau,rice)
5      USE updaterad
6      use tracer_mod, only: igcm_h2o_vap, igcm_h2o_ice
7      USE comcstfi_h
8      use dimradmars_mod, only: naerkind
9      implicit none
10c------------------------------------------------------------------
11c  This routine is used to form clouds when a parcel of the GCM is
12c    saturated. It is a simplified scheme, and there is almost no
13c    microphysics involved. When the air is saturated, water-ice
14c    clouds form on a fraction of the dust particles, specified by
15c    the constant called "ccn_factor". There is no supersaturation,
16c    and no nucleation rates computed. A more accurate scheme can
17c    be found in the routine called "improvedclouds.F".
18
19c  Modif de zq si saturation dans l'atmosphere
20c  si zq(ig,l)> zqsat(ig,l) ->    zq(ig,l)=zqsat(ig,l)
21c  Le test est effectue de bas en haut. L'eau condensee
22c    (si saturation) est remise dans la couche en dessous.
23c  L'eau condensee dans la couche du bas est deposee a la surface
24
25c  Authors: Franck Montmessin (water ice scheme)
26c           Francois Forget (changed nuclei density & outputs)
27c           Ehouarn Millour (sept.2008, tracers are now handled
28c                                   by name and not fixed index)
29c           J.-B. Madeleine (developed a single routine called
30c                            simpleclouds.F, and corrected calculations
31c                            of the typical CCN profile, Oct. 2011)
32c------------------------------------------------------------------
33!#include "dimensions.h"
34!#include "dimphys.h"
35#include "callkeys.h"
36!#include "tracer.h"
37!#include "comgeomfi.h"
38!#include "dimradmars.h"
39
40c------------------------------------------------------------------
41c     Arguments:
42c     ---------
43c     Inputs:
44      INTEGER ngrid,nlay
45      integer nq                 ! nombre de traceurs
46      REAL ptimestep             ! pas de temps physique (s)
47      REAL pplay(ngrid,nlay)     ! pression au milieu des couches (Pa)
48      REAL pzlay(ngrid,nlay)     ! altitude at the middle of the layers
49      REAL pt(ngrid,nlay)        ! temperature at the middle of the
50                                 !   layers (K)
51      REAL pdt(ngrid,nlay)       ! tendance temperature des autres
52                                 !   param.
53      real pq(ngrid,nlay,nq)     ! traceur (kg/kg)
54      real pdq(ngrid,nlay,nq)    ! tendance avant condensation
55                                 !   (kg/kg.s-1)
56      REAL tau(ngrid,naerkind)   ! Column dust optical depth at each point
57
58c     Output:
59      REAL rice(ngrid,nlay)      ! Ice mass mean radius (m)
60                                 ! (r_c in montmessin_2004)
61      real pdqcloud(ngrid,nlay,nq) ! tendance de la condensation
62                                   !   H2O(kg/kg.s-1)
63      REAL pdtcloud(ngrid,nlay)    ! tendance temperature due
64                                   !   a la chaleur latente
65
66c------------------------------------------------------------------
67c     Local variables:
68
69      LOGICAL,SAVE :: firstcall = .true.
70           
71      REAL rhocloud(ngrid,nlay)  ! Cloud density (kg.m-3)
72
73      INTEGER ig,l
74
75      REAL zq(ngrid,nlay,nq)    ! local value of tracers
76      REAL zq0(ngrid,nlay,nq)   ! local initial value of tracers
77      REAL zt(ngrid,nlay)       ! local value of temperature
78      REAL zqsat(ngrid,nlay)    ! saturation
79      REAL*8 dzq                      ! masse de glace echangee (kg/kg)
80      REAL lw                         !Latent heat of sublimation (J.kg-1)
81      REAL,PARAMETER :: To=273.15     ! reference temperature, T=273.15 K
82      real rdusttyp(ngrid,nlay) ! Typical dust geom. mean radius (m)
83      REAL ccntyp(ngrid,nlay)
84                                      ! Typical dust number density (#/kg)
85c     CCN reduction factor
86c      REAL, PARAMETER :: ccn_factor = 4.5  !! comme TESTS_JB // 1. avant
87     
88
89c-----------------------------------------------------------------------
90c    1. initialisation
91c    -----------------
92
93c    On "update" la valeur de q(nq) (water vapor) et temperature.
94c    On effectue qqes calculs preliminaires sur les couches :
95
96      do l=1,nlay
97        do ig=1,ngrid
98          zq(ig,l,igcm_h2o_vap)=
99     &      pq(ig,l,igcm_h2o_vap)+pdq(ig,l,igcm_h2o_vap)*ptimestep
100          zq(ig,l,igcm_h2o_vap)=max(zq(ig,l,igcm_h2o_vap),1.E-30) ! FF 12/2004
101          zq0(ig,l,igcm_h2o_vap)=zq(ig,l,igcm_h2o_vap)
102          zt(ig,l)=pt(ig,l)+ pdt(ig,l)*ptimestep
103
104          zq(ig,l,igcm_h2o_ice)=
105     &      pq(ig,l,igcm_h2o_ice)+pdq(ig,l,igcm_h2o_ice)*ptimestep
106          zq(ig,l,igcm_h2o_ice)=max(zq(ig,l,igcm_h2o_ice),0.) ! FF 12/2004
107          zq0(ig,l,igcm_h2o_ice)=zq(ig,l,igcm_h2o_ice)
108        enddo
109      enddo
110
111
112      pdqcloud(1:ngrid,1:nlay,1:nq)=0
113      pdtcloud(1:ngrid,1:nlay)=0
114
115c     ----------------------------------------------
116c
117c
118c     Rapport de melange a saturation dans la couche l : -------
119c     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
120
121      call watersat(ngrid*nlay,zt,pplay,zqsat)
122
123c     taux de condensation (kg/kg/s-1) dans les differentes couches
124c     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
125
126      do l=1,nlay
127        do ig=1,ngrid
128
129          if (zq(ig,l,igcm_h2o_vap).ge.zqsat(ig,l))then  !  Condensation
130            dzq=zq(ig,l,igcm_h2o_vap)-zqsat(ig,l)               
131          elseif(zq(ig,l,igcm_h2o_vap).lt.zqsat(ig,l))then  ! Sublimation
132            dzq=-min(zqsat(ig,l)-zq(ig,l,igcm_h2o_vap),
133     &               zq(ig,l,igcm_h2o_ice))
134          endif
135
136c         Water Mass change
137c         ~~~~~~~~~~~~~~~~~
138          zq(ig,l,igcm_h2o_ice)=zq(ig,l,igcm_h2o_ice)+dzq
139          zq(ig,l,igcm_h2o_vap)=zq(ig,l,igcm_h2o_vap)-dzq
140         
141
142        enddo ! of do ig=1,ngrid
143      enddo ! of do l=1,nlay
144
145c     Tendance finale
146c     ~~~~~~~~~~~~~~~
147      do l=1, nlay
148        do ig=1,ngrid
149          pdqcloud(ig,l,igcm_h2o_vap)=(zq(ig,l,igcm_h2o_vap)
150     &                            -zq0(ig,l,igcm_h2o_vap))/ptimestep
151          pdqcloud(ig,l,igcm_h2o_ice) =
152     &      (zq(ig,l,igcm_h2o_ice) - zq0(ig,l,igcm_h2o_ice))/ptimestep
153          lw=(2834.3-0.28*(zt(ig,l)-To)-0.004*(zt(ig,l)-To)**2)*1.e+3
154          pdtcloud(ig,l)=-pdqcloud(ig,l,igcm_h2o_vap)*lw/cpp
155        end do
156      end do
157
158c     ice crystal radius
159      do l=1, nlay
160        do ig=1,ngrid
161          call updaterice_typ(zq(ig,l,igcm_h2o_ice),
162     &       tau(ig,1),pzlay(ig,l),rice(ig,l))
163        end do
164      end do
165
166c------------------------------------------------------------------
167      return
168      end
Note: See TracBrowser for help on using the repository browser.