source: trunk/LMDZ.VENUS/libf/phyvenus/phytrac_chimie.F @ 2560

Last change on this file since 2560 was 2523, checked in by slebonnois, 4 years ago

SL: Venus GCM outputs for VCD + bug correction in moldiff_red

  • Property svn:executable set to *
File size: 11.5 KB
Line 
1      subroutine phytrac_chimie (
2     $                    debutphy,
3     $                    gmtime,
4     $                    nqmax,
5     $                    nlon,
6     $                    lat,
7     $                    lon,
8     $                    nlev,
9     $                    pdtphys,
10     $                    temp,
11     $                    pplay,
12     $                    trac,
13     $                    d_tr_chem,
14     $                    iter)
15
16      use chemparam_mod
17      use conc, only: mmean,rnew
18
19      implicit none
20     
21#include "clesphys.h"
22#include "YOMCST.h"
23
24!===================================================================
25!     input
26!===================================================================
27
28      integer :: nlon, nlev                     ! number of gridpoints and levels
29      integer :: nqmax                          ! number of tracers
30
31      real :: gmtime                            ! day fraction
32      real :: pdtphys                           ! phytrac_chimie timestep (s)
33      real, dimension(nlon,nlev) :: temp        ! temperature (k)
34      real, dimension(nlon,nlev) :: pplay       ! pressure (pa)
35      real, dimension(nlon,nlev,nqmax) :: trac  ! tracer mass mixing ratio
36
37      logical :: debutphy                       ! first call flag
38
39!===================================================================
40!     output
41!===================================================================
42
43      real, dimension(nlon,nlev,nqmax) :: d_tr_chem  ! chemical tendency for each tracer
44      integer, dimension(nlon,nlev) :: iter          ! chemical iterations
45
46!===================================================================
47!     local
48!===================================================================
49
50      real :: sza_local   ! solar zenith angle (deg)
51      real :: lon_sun
52
53      integer :: i, iq
54      integer :: ilon, ilev
55
56      real  lat(nlon), lat_local(nlon)
57      real  lon(nlon), lon_local(nlon)
58
59      real, dimension(nlon,nlev) :: mrtwv, mrtsa ! total water and total sulfuric acid
60      real, dimension(nlon,nlev) :: mrwv, mrsa   ! gas-phase water and gas-phase sulfuric acid
61      real, dimension(nlon,nlev) :: trac_sum
62      real, dimension(nlon,nlev,nqmax) :: ztrac  ! local tracer mixing ratio
63     
64!===================================================================
65!     first call : initialisations
66!===================================================================
67
68      if (debutphy) then
69     
70!!--- Adjustment of Helium amount
71!        if (i_he/=0) then
72!           trac(:,:,i_he)=trac(:,:,i_he)*20.
73!        endif
74!!---
75
76!-------------------------------------------------------------------
77!        case of tracers re-initialisation with chemistry
78!-------------------------------------------------------------------
79         if (reinit_trac .and. ok_chem) then
80
81!!! in this reinitialisation, trac is VOLUME mixing ratio
82            print*, "Tracers are re-initialised"
83            trac(:,:,:) = 1.0e-30
84
85            if ((i_ocs /= 0) .and. (i_co /= 0) .and. (i_hcl /= 0)
86     $           .and. (i_so2 /= 0) .and. (i_h2o /= 0) .and. (i_n2/ = 0)
87     $           .and. (i_co2 /= 0)) then
88
89               trac(:,1:22,i_ocs) = 3.e-6
90               trac(:,1:22,i_co)  = 25.e-6
91               trac(:,:,i_hcl)    = 0.4e-6
92               trac(:,1:22,i_so2) = 7.e-6
93               trac(:,1:22,i_h2o) = 30.e-6
94               trac(:,:,i_n2)     = 0.35e-1
95   
96!          ensure that sum of mixing ratios = 1
97
98               trac_sum(:,:) = 0.
99
100               do iq = 1,nqmax - nmicro
101                  if (iq /= i_co2) then
102                     trac_sum(:,:) = trac_sum(:,:) + trac(:,:,iq)
103                  end if
104               end do
105
106!          initialise co2
107
108               trac(:,:,i_co2) = 1. - trac_sum(:,:)
109
110            else
111               write(*,*) "at least one tracer is missing : stop"
112               stop
113            end if
114       
115
116! update mmean
117            mmean(:,:) = 0.
118            do iq = 1,nqmax - nmicro
119               mmean(:,:) = mmean(:,:)+trac(:,:,iq)*m_tr(iq)
120               rnew(:,:) = 8.314/mmean(:,:)*1.e3     ! J/kg K
121            enddo
122
123!           convert volume to mass mixing ratio
124            do iq = 1,nqmax - nmicro
125               trac(:,:,iq) = trac(:,:,iq)*m_tr(iq)/mmean(:,:)
126            end do
127   
128         end if  ! reinit_trac
129
130!-------------------------------------------------------------------
131!        case of detailed microphysics without chemistry
132!-------------------------------------------------------------------
133         if (.not. ok_chem .and. ok_cloud .and. cl_scheme == 2) then
134
135!           convert mass to volume mixing ratio
136
137            do iq = 1,nqmax - nmicro
138               ztrac(:,:,iq) = trac(:,:,iq)*mmean(:,:)/m_tr(iq)
139            end do
140         
141!           initialise microphysics
142 
143            call vapors4muphy_ini(nlon,nlev,ztrac)
144
145!           convert volume to mass mixing ratio
146
147            do iq = 1,nqmax - nmicro
148               trac(:,:,iq) = ztrac(:,:,iq)*m_tr(iq)/mmean(:,:)
149            end do
150   
151         end if
152
153      end if  ! debutphy
154
155!===================================================================
156!     convert mass to volume mixing ratio : gas phase
157!===================================================================
158
159      do iq = 1,nqmax - nmicro
160         ztrac(:,:,iq) = max(trac(:,:,iq)*mmean(:,:)/m_tr(iq), 1.e-30)
161      end do
162
163!===================================================================
164!     microphysics: simplified scheme (phd aurelien stolzenbach)
165!===================================================================
166
167      if (ok_cloud .and. cl_scheme == 1) then
168
169!     convert mass to volume mixing ratio : liquid phase
170
171         ztrac(:,:,i_h2so4liq) = max(trac(:,:,i_h2so4liq)
172     $                             *mmean(:,:)/m_tr(i_h2so4liq), 1.e-30)
173         ztrac(:,:,i_h2oliq) = max(trac(:,:,i_h2oliq)
174     $                             *mmean(:,:)/m_tr(i_h2oliq), 1.e-30)
175             
176!     total water and sulfuric acid (gas + liquid)
177
178         mrtwv(:,:) = ztrac(:,:,i_h2o) + ztrac(:,:,i_h2oliq)
179         mrtsa(:,:) = ztrac(:,:,i_h2so4) + ztrac(:,:,i_h2so4liq)
180
181!     all water and sulfuric acid is put in the gas-phase
182
183         mrwv(:,:) = mrtwv(:,:)
184         mrsa(:,:) = mrtsa(:,:)
185
186!     call microphysics
187
188         call new_cloud_venus(nlev, nlon, temp, pplay,
189     $                        mrtwv, mrtsa, mrwv, mrsa)
190
191!     update water vapour and sulfuric acid
192
193         ztrac(:,:,i_h2o) = mrwv(:,:)
194         ztrac(:,:,i_h2oliq) = mrtwv(:,:) - ztrac(:,:,i_h2o)
195     
196         ztrac(:,:,i_h2so4) = mrsa(:,:)
197         ztrac(:,:,i_h2so4liq) = mrtsa(:,:) - ztrac(:,:,i_h2so4)
198
199      end if  ! simplified scheme
200
201!===================================================================
202!     microphysics: detailed scheme (phd sabrina guilbon)
203!     !!! to be confirmed whether mad_muphy expects mmr or vmr for h2o and h2so4
204!===================================================================
205
206      if (ok_cloud .and. cl_scheme == 2) then
207
208         do iq = nqmax-nmicro+1,nqmax
209            ztrac(:,:,iq) = trac(:,:,iq)
210         end do
211
212         do ilon = 1,nlon       
213            do ilev = 1, nlev
214               if (temp(ilon,ilev) < 500.) then
215                  call mad_muphy(pdtphys,                               ! timestep
216     $                           temp(ilon,ilev),pplay(ilon,ilev),      ! temperature and pressure
217     $                           ztrac(ilon,ilev,i_h2o),
218     $                           ztrac(ilon,ilev,i_h2so4),     
219     $                           ztrac(ilon,ilev,i_m0_aer),
220     $                           ztrac(ilon,ilev,i_m3_aer),   
221     $                           ztrac(ilon,ilev,i_m0_mode1drop),
222     $                           ztrac(ilon,ilev,i_m0_mode1ccn),
223     $                           ztrac(ilon,ilev,i_m3_mode1sa),
224     $                           ztrac(ilon,ilev,i_m3_mode1w),   
225     $                           ztrac(ilon,ilev,i_m3_mode1ccn),   
226     $                           ztrac(ilon,ilev,i_m0_mode2drop),
227     $                           ztrac(ilon,ilev,i_m0_mode2ccn),
228     $                           ztrac(ilon,ilev,i_m3_mode2sa),
229     $                           ztrac(ilon,ilev,i_m3_mode2w),
230     $                           ztrac(ilon,ilev,i_m3_mode2ccn))
231               else
232                  ztrac(ilon,ilev,i_m0_aer)       = 0.
233                  ztrac(ilon,ilev,i_m3_aer)       = 0.
234                  ztrac(ilon,ilev,i_m0_mode1drop) = 0.
235                  ztrac(ilon,ilev,i_m0_mode1ccn)  = 0.
236                  ztrac(ilon,ilev,i_m3_mode1sa)   = 0.
237                  ztrac(ilon,ilev,i_m3_mode1w)    = 0.
238                  ztrac(ilon,ilev,i_m3_mode1ccn)  = 0.
239                  ztrac(ilon,ilev,i_m0_mode2drop) = 0.
240                  ztrac(ilon,ilev,i_m0_mode2ccn)  = 0.
241                  ztrac(ilon,ilev,i_m3_mode2sa)   = 0.
242                  ztrac(ilon,ilev,i_m3_mode2w)    = 0.
243                  ztrac(ilon,ilev,i_m3_mode2ccn)  = 0.
244               end if
245            end do
246         end do
247
248      end if  ! detailed scheme
249           
250!===================================================================
251!     photochemistry
252!===================================================================
253
254      if (ok_chem) then
255
256         lon_sun = (0.5 - gmtime)*2.*rpi
257         lon_local(:) = lon(:)*rpi/180.
258         lat_local(:) = lat(:)*rpi/180.
259
260         do ilon = 1,nlon
261
262!           solar zenith angle
263
264            sza_local = acos(cos(lat_local(ilon))*cos(lon_local(ilon))
265     $                 *cos(lon_sun) + cos(lat_local(ilon))
266     $                 *sin(lon_local(ilon))*sin(lon_sun))*180./rpi
267     
268            call photochemistry_venus(nlev, nlon, pdtphys,
269     $                                pplay(ilon,:)/100.,
270     $                                temp(ilon,:),
271     $                                ztrac(ilon,:,:),
272     $                                mmean(ilon,:),
273     $                                sza_local,
274     $                                lon(ilon), lat(ilon),
275     $                                nqmax, iter(ilon,:))
276
277         end do
278
279      end if  ! ok_chem
280
281!===================================================================
282!     compute tendencies
283!===================================================================
284
285! update mmean
286            mmean(:,:) = 0.
287            do iq = 1,nqmax - nmicro
288               mmean(:,:) = mmean(:,:)+trac(:,:,iq)*m_tr(iq)
289               rnew(:,:) = 8.314/mmean(:,:)*1.e3     ! J/kg K
290            enddo
291
292!     gas phase
293
294      do iq = 1,nqmax - nmicro
295         ztrac(:,:,iq) = max(ztrac(:,:,iq)*m_tr(iq)/mmean(:,:),
296     $                       1.e-30)
297         d_tr_chem(:,:,iq) = (ztrac(:,:,iq) - trac(:,:,iq))/pdtphys
298      end do
299
300!     liquid phase or moments
301
302      if (ok_cloud .and. cl_scheme == 1) then
303         ztrac(:,:,i_h2so4liq) = max(ztrac(:,:,i_h2so4liq)
304     $                               *m_tr(i_h2so4liq)/mmean(:,:),
305     $                               1.e-30)
306         ztrac(:,:,i_h2oliq)   = max(ztrac(:,:,i_h2oliq)
307     $                               *m_tr(i_h2oliq)/mmean(:,:),
308     $                               1.e-30)
309         d_tr_chem(:,:,i_h2so4liq) = (ztrac(:,:,i_h2so4liq)
310     $                              - trac(:,:,i_h2so4liq))/pdtphys
311         d_tr_chem(:,:,i_h2oliq) = (ztrac(:,:,i_h2oliq)
312     $                            - trac(:,:,i_h2oliq))/pdtphys
313      end if
314
315
316      if (ok_cloud .and. cl_scheme == 2) then
317         do iq = nqmax-nmicro+1,nqmax
318            d_tr_chem(:,:,iq) = (ztrac(:,:,iq) - trac(:,:,iq))/pdtphys
319         end do
320      end if
321
322      end subroutine phytrac_chimie
Note: See TracBrowser for help on using the repository browser.