source: trunk/LMDZ.PLUTO/libf/phypluto/aerosol_mod.F90 @ 3506

Last change on this file since 3506 was 3353, checked in by afalco, 6 months ago

Pluto PCM:
Added zrecast & old sedim ;
Choose haze file ;
AF

File size: 4.2 KB
Line 
1!==================================================================
2module aerosol_mod
3implicit none
4
5!==================================================================
6
7!  aerosol indexes: these are initialized to be 0 if the
8!                 corresponding aerosol was not activated in callphys.def
9!                 -- otherwise a value is set via iniaerosol
10      integer, save :: iaero_haze = 0
11      integer, save :: i_haze = 0
12      logical, save, protected :: noaero = .false.
13!$OMP THREADPRIVATE(iaero_haze,i_haze,noaero)
14
15! two-layer simple aerosol model
16      integer, save, protected :: iaero_back2lay = 0
17! N-layer aerosol model (replaces the 2-layer and hard-coded clouds)
18      integer,dimension(:), allocatable, save, protected :: iaero_nlay
19!$OMP THREADPRIVATE(iaero_back2lay,iaero_nlay)
20
21! Generic aerosols
22      integer, dimension(:), allocatable, save, protected :: iaero_generic
23      integer, dimension(:), allocatable, save, protected :: i_rgcs_ice
24!$OMP THREADPRIVATE(iaero_generic,i_rgcs_ice)
25
26!==================================================================
27
28contains
29
30  !==================================================================
31   subroutine haze_prof(ngrid,nlayer,zzlay,pplay,pt,reffrad,profmmr)
32!==================================================================
33!     Purpose
34!     -------
35!     Get fixed haze properties
36!     profile of haze (from txt file) and fixed radius profile
37!
38!==================================================================
39      use radinc_h, only: naerkind
40      use datafile_mod
41      use tracer_h
42      use comcstfi_mod, only: r, pi
43
44!-----------------------------------------------------------------------
45!     Arguments
46      Implicit none
47
48      integer,intent(in) :: ngrid
49      integer,intent(in) :: nlayer
50      real,intent(in) :: zzlay(ngrid,nlayer)
51      real,intent(in) :: pplay(ngrid,nlayer)
52      real,intent(in) :: pt(ngrid,nlayer)
53      real, intent(in) :: reffrad(ngrid,nlayer,naerkind)      ! haze particles radii (m)
54
55      real, intent(out) :: profmmr(ngrid,nlayer)              ! mmr haze kg/kg
56
57!     Local variables
58      integer :: iaer,l,ig,ifine
59
60      LOGICAL firstcall
61      SAVE firstcall
62      DATA firstcall/.true./
63
64      !!read altitudes and haze mmrs
65      integer Nfine
66      !parameter(Nfine=21)
67      parameter(Nfine=701)
68      character(len=100) :: file_path
69      character(len=100) :: file_name
70      real,save :: levdat(Nfine),densdat(Nfine)
71
72!---------------- INPUT ------------------------------------------------
73
74      !! Read data
75      IF (firstcall) then
76        firstcall=.false.
77
78        if (hazemmr_file/='None')then
79            file_name = hazemmr_file
80            print*, 'Read Haze MMR file: ',hazemmr_file
81        else if(hazedens_file/='None')then
82            file_name = hazedens_file
83            print*, 'Read Haze density: ',hazedens_file
84        else
85            STOP "No filename given for haze profile. Either set hazemmr_file or hazedens_file"
86        endif
87
88        file_path=trim(datadir)//'/haze_prop/'//file_name
89        open(224,file=file_path,form='formatted')
90        do ifine=1,Nfine
91           read(224,*) levdat(ifine), densdat(ifine)
92        enddo
93        close(224)
94        print*, 'Read Haze profile: ',file_path
95      ENDIF
96
97      !! Interpolate on the model vertical grid
98      do ig=1,ngrid
99        CALL interp_line(levdat,densdat,Nfine,zzlay(ig,:)/1000.,profmmr(ig,:),nlayer)
100      enddo
101
102      !! Get profile Mass mixing ratio from number density:    part.cm-3 --> m-3 --> m3 m-3
103      !                                --> kg m-3 --> kg/kg
104      do iaer=1,naerkind
105            if(iaer.eq.iaero_haze.and.hazedens_file/='None') then !AF24 activate/deactivate mmr or part density
106              !print*, 'Haze profile is fixed'
107              do ig=1,ngrid
108                 do l=1,nlayer
109                    !from number density in cm-3
110                    profmmr(ig,l)=profmmr(ig,l)*1.e6*4./3.*pi*reffrad(ig,l,iaer)**3*rho_q(i_haze)/(pplay(ig,l)/(r*pt(ig,l)))
111                ! print*, profmmr(ig,l)
112                 enddo
113              enddo
114            endif
115      enddo
116   end subroutine haze_prof
117!==================================================================
118
119
120end module aerosol_mod
121!==================================================================
Note: See TracBrowser for help on using the repository browser.