source: trunk/LMDZ.MARS/libf/aeronomars/surfacearea.F @ 3567

Last change on this file since 3567 was 3466, checked in by emillour, 2 months ago

Mars PCM:
More tidying in aeronomars:

  • remove unused "inv.F" and remove "dtridgl.F" which is not used here and is a duplicate of the "dtridgl" routine in phymars/swr_toon.F
  • turn aeronomars routines to modules, for those which aren't in modules yet.

EM

File size: 5.5 KB
Line 
1      module surfacearea_mod
2     
3      implicit none
4     
5      contains
6     
7      subroutine surfacearea(ngrid, nlay, naerkind, ptimestep,
8     $                       pplay, pzlay,
9     $                       pt, pq, pdq, nq, rdust, rice, tau,
10     $                       tauscaling,
11     $                       surfdust, surfice)
12
13      use tracer_mod, only: nuice_sed, igcm_dust_number,
14     &                      igcm_ccn_number, varian, ccn_factor
15      use conc_mod, only: rnew
16      use comcstfi_h, only: pi
17      use wstats_mod, only: wstats
18      implicit none
19
20!==========================================================================
21!     calculation of the ice and dust surface area (m2/m3)
22!     available for heterogeneous reactions
23!
24!     Franck Lefevre
25!     version 1.2 april 2012
26!==========================================================================
27
28      include "callkeys.h"
29
30! input
31
32      integer,intent(in) :: ngrid, nlay, naerkind
33      integer,intent(in) :: nq               ! number of tracers
34      real,intent(in) :: ptimestep           ! physics time step (s)
35      real,intent(in) :: pplay(ngrid,nlay)   ! pressure at mid-layers (Pa)
36      real,intent(in) :: pzlay(ngrid,nlay)   ! altitude at mid-layers (m)
37      real,intent(in) :: pt(ngrid,nlay)      ! temperature at mid-layers (K)
38      real,intent(in) :: pq(ngrid,nlay,nq)   ! tracers (kg/kg)
39      real,intent(in) :: pdq(ngrid,nlay,nq)  ! physical tendency (kg/kg.s-1)
40      real,intent(in) :: rdust(ngrid,nlay)   ! dust geometric mean radius (m)
41      real,intent(in) :: rice(ngrid,nlay)    ! ice mass mean radius (m)
42      real,intent(in) :: tau(ngrid,naerkind) ! column dust optical depth at each point
43      real,intent(in) :: tauscaling(ngrid)   ! conversion factor for dust amount
44
45! output
46
47      real,intent(out) :: surfdust(ngrid,nlay) ! dust surface area (m2/m3)
48      real,intent(out) :: surfice(ngrid,nlay)  ! water-ice surface area (m2/m3)
49
50! local
51
52      integer    :: l, ig
53      real       :: rho                     ! density (kg/m3)
54      real       :: dustnd, icend           ! uodated dust and ice number densities (kg/kg)
55      real, save :: factor_ice, factor_dust ! multiplying factor to compute total surface area
56                                            ! from the mass-mean radius
57      real       :: sigma_ice, sigma_dust   ! variance of the ice and dust distributions
58      real       :: ccntyp                  ! typical dust number density (#/kg)
59                                            ! (microphys = false)
60      real       :: rdusttyp                ! typical dust radius (m)
61                                            ! (microphys = false)
62
63      logical, save :: firstcall = .true.
64
65!$OMP THREADPRIVATE(factor_ice, factor_dust,firstcall)
66
67!==========================================================================
68
69      if (firstcall) then ! compute the multiplying factors
70         sigma_dust  = varian
71         sigma_ice   = sqrt(log(nuice_sed + 1.))
72         factor_dust = exp(0.5*(log(sigma_dust))**2)
73         factor_ice  = exp(0.5*(log(sigma_ice))**2)
74         write(*,*) 'surfacearea : factor_dust = ', factor_dust
75         write(*,*) 'surfacearea : factor_ice  = ', factor_ice
76         firstcall = .false.
77      end if
78
79      if (microphys) then ! improvedclouds
80         do l = 1,nlay
81            do ig = 1,ngrid
82!              atmospheric density
83               rho = pplay(ig,l)/(rnew(ig,l)*pt(ig,l))
84!              updated dust number density
85               dustnd = pq(ig,l,igcm_dust_number)
86     $                + pdq(ig,l,igcm_dust_number)*ptimestep
87!              updated ice number density
88               icend  = pq(ig,l,igcm_ccn_number)
89     $                + pdq(ig,l,igcm_ccn_number)*ptimestep
90!              dust surface area
91               surfdust(ig,l) = factor_dust*dustnd*rho*tauscaling(ig)
92     $                          *4.*pi*rdust(ig,l)**2
93!              ice surface area
94               surfice(ig,l)  = factor_ice*icend*rho*tauscaling(ig)
95     $                          *4.*pi*rice(ig,l)**2
96            end do
97         end do
98      else               ! simpleclouds
99         do l = 1,nlay
100            do ig = 1,ngrid
101!              atmospheric density
102               rho = pplay(ig,l)/(rnew(ig,l)*pt(ig,l))
103!              typical dust radius
104               rdusttyp = max(.8e-6*exp(-pzlay(ig,l)/18000.),1.e-9)
105!              typical dust number density
106               ccntyp = 1.3e+8*max(tau(ig,1),0.001)/0.1
107     $                  *exp(-pzlay(ig,l)/10000.)
108               ccntyp = ccntyp/ccn_factor
109               if (rice(ig,l) .gt. rdust(ig,l)) then
110                  surfdust(ig,l) = factor_dust*ccntyp*(ccn_factor - 1.)
111     $                             *rho*4.*pi*rdusttyp**2
112                  surfice(ig,l)  = factor_ice*ccntyp*4.*pi*rice(ig,l)**2
113               else
114                  surfdust(ig,l) = factor_dust*ccntyp*ccn_factor
115     $                             *rho*4.*pi*rdusttyp**2
116                  surfice(ig,l)  = 0.
117               end if
118            end do
119         end do
120      end if         ! of microphys
121
122! write diagnostics in micron2/cm3
123     
124      call wstats(ngrid,"surfdust", "Dust surface area",
125     $            "micron2 cm-3",3,surfdust*1.e6)
126      call wstats(ngrid,"surfice", "Ice cloud surface area",
127     $            "micron2 cm-3",3,surfice*1.e6)
128
129      call writediagfi(ngrid,"surfdust", "Dust surface area",
130     $            "micron2 cm-3",3,surfdust*1.e6)
131      call writediagfi(ngrid,"surfice", "Ice cloud surface area",
132     $            "micron2 cm-3",3,surfice*1.e6)
133
134      end subroutine surfacearea
135
136      end module surfacearea_mod
Note: See TracBrowser for help on using the repository browser.