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

Last change on this file since 1112 was 1047, checked in by emillour, 11 years ago

Mars GCM:

  • IMPORTANT CHANGE: Removed all reference/use of ngridmx (dimphys.h) in routines (necessary prerequisite to using parallel dynamics); in most cases this just means adding 'ngrid' as routine argument, and making local saved variables allocatable (and allocated at first call). In the process, had to convert many *.h files to equivalent modules: yomaer.h => yomaer_h.F90 , surfdat.h => surfdat_h.F90 , comsaison.h => comsaison_h.F90 , yomlw.h => yomlw_h.F90 , comdiurn.h => comdiurn_h.F90 , dimradmars.h => dimradmars_mod.F90 , comgeomfi.h => comgeomfi_h.F90, comsoil.h => comsoil_h.F90 , slope.h => slope_mod.F90
  • Also updated EOF routines, everything is now in eofdump_mod.F90
  • Removed unused routine lectfux.F (in dyn3d)

EM

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