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

Last change on this file since 3652 was 3650, checked in by afalco, 4 months ago

Pluto: fix bug in aerosol_mod in parallel. Setting deprecated aerohaze=true gives an error.
AF

File size: 4.0 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!$OMP THREADPRIVATE(iaero_haze,i_haze)
13
14!==================================================================
15
16contains
17
18  !==================================================================
19   subroutine haze_prof(ngrid,nlayer,zzlay,pplay,pt,reffrad,profmmr)
20!==================================================================
21!     Purpose
22!     -------
23!     Get fixed haze properties
24!     profile of haze (from txt file) and fixed radius profile
25!
26!==================================================================
27      use radinc_h, only: naerkind
28      use datafile_mod
29      use tracer_h
30      use comcstfi_mod, only: r, pi
31      use mod_phys_lmdz_para, only : is_master
32
33!-----------------------------------------------------------------------
34!     Arguments
35      Implicit none
36
37      integer,intent(in) :: ngrid
38      integer,intent(in) :: nlayer
39      real,intent(in) :: zzlay(ngrid,nlayer)
40      real,intent(in) :: pplay(ngrid,nlayer)
41      real,intent(in) :: pt(ngrid,nlayer)
42      real, intent(in) :: reffrad(ngrid,nlayer,naerkind)      ! haze particles radii (m)
43
44      real, intent(out) :: profmmr(ngrid,nlayer)              ! mmr haze kg/kg
45
46!     Local variables
47      integer :: iaer,l,ig,ifine
48
49      LOGICAL,SAVE :: firstcall=.true.
50!$OMP THREADPRIVATE(firstcall)
51
52      !!read altitudes and haze mmrs
53      integer Nfine
54      !parameter(Nfine=21)
55      parameter(Nfine=701)
56      character(len=100) :: file_path
57      character(len=100) :: file_name
58    !   real,save :: levdat(Nfine),densdat(Nfine)
59      real,save,allocatable :: levdat(:)
60      real,save,allocatable :: densdat(:)
61
62!---------------- INPUT ------------------------------------------------
63
64      !! Read data
65      IF (firstcall) then
66        firstcall=.false.
67
68        if (hazemmr_file/='None')then
69            file_name = hazemmr_file
70            print*, 'Read Haze MMR file: ',hazemmr_file
71        else if(hazedens_file/='None')then
72            file_name = hazedens_file
73            print*, 'Read Haze density: ',hazedens_file
74        else
75            STOP "No filename given for haze profile. Either set hazemmr_file or hazedens_file"
76        endif
77!$OMP MASTER
78        if(.not.allocated(levdat)) then
79            allocate(levdat(Nfine))
80        endif
81        if(.not.allocated(densdat)) then
82            allocate(densdat(Nfine))
83        endif
84
85
86        file_path=trim(datadir)//'/haze_prop/'//file_name
87        open(224,file=file_path,form='formatted')
88        do ifine=1,Nfine
89        read(224,*) levdat(ifine), densdat(ifine)
90        enddo
91        close(224)
92        print*, 'Read Haze profile: ',file_path
93!$OMP END MASTER
94!$OMP BARRIER
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.