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

Last change on this file since 2683 was 2622, checked in by slebonnois, 3 years ago

SL: VENUS update (i) bug correction (2 bugs, phytrac and physiq), affected meam molec mass computations... (ii) updates for VCD 2.0 (iii) aeropacity: for latitudinal variations of the cloud distribution

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