Ignore:
Timestamp:
Jan 31, 2024, 4:36:51 PM (10 months ago)
Author:
afalco
Message:

Pluto PCM:
Imported condense n2 from pluto.old.
Aerosol data from Pluto.old not yet working.
AF

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.PLUTO/libf/phypluto/aerosol_mod.F90

    r3193 r3195  
    88!                 corresponding aerosol was not activated in callphys.def
    99!                 -- otherwise a value is set via iniaerosol
    10       integer, save, protected :: iaero_n2 = 0
    11       integer, save, protected :: iaero_dust = 0
    12       integer, save, protected :: iaero_h2so4 = 0
     10      integer, save :: iaero_haze = 0
     11      integer, save :: i_haze = 0
    1312      logical, save, protected :: noaero = .false.
    14 !$OMP THREADPRIVATE(iaero_n2,iaero_h2o,iaero_dust,iaero_h2so4,noaero)
     13!$OMP THREADPRIVATE(iaero_haze,i_haze,noaero)
    1514
    1615! two-layer simple aerosol model
    1716      integer, save, protected :: iaero_back2lay = 0
    18  ! NH3 cloud
    19       integer, save, protected :: iaero_nh3 = 0
    2017! N-layer aerosol model (replaces the 2-layer and hard-coded clouds)
    2118      integer,dimension(:), allocatable, save, protected :: iaero_nlay
    22 ! Auroral aerosols
    23       integer, save, protected :: iaero_aurora = 0
    24 !$OMP THREADPRIVATE(iaero_back2lay,iaero_nh3,iaero_nlay,iaero_aurora)
     19!$OMP THREADPRIVATE(iaero_back2lay,iaero_nlay)
    2520
    2621! Generic aerosols
     
    3328contains
    3429
    35   SUBROUTINE iniaerosol
     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
    3643
    37   use mod_phys_lmdz_para, only : is_master
    38   use radinc_h, only: naerkind
    39   use tracer_h, only: n_rgcs, nqtot, is_rgcs
    40   use callkeys_mod, only: aeron2, dusttau, &
    41                           aeroback2lay, aeronh3, nlayaero, aeronlay, &
    42                           aeroaurora, aerogeneric
     44!-----------------------------------------------------------------------
     45!     Arguments
     46      Implicit none
    4347
    44   IMPLICIT NONE
    45 !=======================================================================
    46 !   subject:
    47 !   --------
    48 !   Initialization related to aerosols
    49 !   (N2 aerosols, dust, water, chemical species, ice...)   
    50 !
    51 !   author: Laura Kerber, S. Guerlet
    52 !   ------
    53 !       
    54 !=======================================================================
     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)
    5554
    56   integer :: i, ia, iq
     55      real, intent(out) :: profmmr(ngrid,nlayer)              ! mmr haze kg/kg
    5756
    58   ! Special case, dyn. allocation for n-layer depending on callphys.def
    59   IF(.NOT.ALLOCATED(iaero_nlay)) ALLOCATE(iaero_nlay(nlayaero))
    60   iaero_nlay(:) = 0
    61   ! Do the same for iaero_generic and i_rgcs_ice
    62   IF (.not. allocated(iaero_generic)) allocate(iaero_generic(aerogeneric))
    63   if (.not. allocated(i_rgcs_ice)) allocate(i_rgcs_ice(aerogeneric))
     57!     Local variables
     58      integer :: iaer,l,ig,ifine
    6459
    65   ! Init of i_rgcs_ice
    66   i_rgcs_ice(:) =0
    67   ia = 1
    68   do iq=1,nqtot
    69     if (is_rgcs(iq) .eq. 1) then
    70         i_rgcs_ice(ia)=iq
    71         ia = ia+1
    72      endif
    73   enddo
     60      LOGICAL firstcall
     61      SAVE firstcall
     62      DATA firstcall/.true./
    7463
    75   iaero_generic(:)=0
    76   ia=0
    77   if (aeron2) then
    78      ia=ia+1
    79      iaero_n2=ia
    80   endif
    81   if (is_master) write(*,*) '--- N2 aerosol = ', iaero_n2
     64      !!read altitudes and haze mmrs
     65      integer Nfine
     66      !parameter(Nfine=21)
     67      parameter(Nfine=701)
     68      character(len=100) :: file_path
     69      real,save :: levdat(Nfine),densdat(Nfine)
    8270
    83   if (dusttau.gt.0) then
    84      ia=ia+1
    85      iaero_dust=ia
    86   endif
    87   if (is_master) write(*,*) '--- Dust aerosol = ', iaero_dust
     71!---------------- INPUT ------------------------------------------------
    8872
    89      
    90   if (aeroback2lay) then
    91      ia=ia+1
    92      iaero_back2lay=ia
    93   endif
    94   if (is_master) write(*,*) '--- Two-layer aerosol = ', iaero_back2lay
     73      !! Read data
     74      IF (firstcall) then
     75        firstcall=.false.
     76        file_path=trim(datadir)//'/haze_prop/hazemmr.txt'
     77        open(224,file=file_path,form='formatted')
     78        do ifine=1,Nfine
     79           read(224,*) levdat(ifine), densdat(ifine)
     80        enddo
     81        close(224)
     82        print*, 'TB22 read Haze MMR profile'
     83      ENDIF
    9584
    96   if (aeronh3) then
    97      ia=ia+1
    98      iaero_nh3=ia
    99   endif
    100   if (is_master) write(*,*) '--- NH3 Cloud = ', iaero_nh3
     85      !! Interpolate on the model vertical grid
     86      do ig=1,ngrid
     87        CALL interp_line(levdat,densdat,Nfine,zzlay(ig,:)/1000.,profmmr(ig,:),nlayer)
     88      enddo
    10189
    102   if (aeronlay) then
    103      do i=1,nlayaero
    104        ia=ia+1
    105        iaero_nlay(i)=ia
    106      enddo
    107   endif
    108   if (is_master) write(*,*) '--- N-layer aerosol = ', iaero_nlay
     90      !! Get profile Mass mixing ratio from number density:    part.cm-3 --> m-3 --> m3 m-3
     91      !                                --> kg m-3 --> kg/kg
     92      do iaer=1,naerkind
     93            if(iaer.eq.iaero_haze.and.1.eq.2) then !TB22 activate/deactivate mmr or part density
     94              !print*, 'Haze profile is fixed'
     95              do ig=1,ngrid
     96                 do l=1,nlayer
     97                    !from number density in cm-3
     98                    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)))
     99                 enddo
     100              enddo
     101            endif
     102      enddo
     103   end subroutine haze_prof
     104!==================================================================
    109105
    110   if (aeroaurora) then
    111      ia=ia+1
    112      iaero_aurora=ia
    113   endif
    114   if (is_master) write(*,*) '--- Auroral aerosols = ', iaero_aurora
    115 
    116   if (aerogeneric .ne. 0) then
    117      do i=1,aerogeneric
    118         ia = ia+1
    119         iaero_generic(i) = ia
    120      enddo
    121   endif
    122      
    123   if (is_master) then
    124     write(*,*)'--- Radiative Generic Condensable Species = ',iaero_generic
    125 
    126     write(*,*) '=== Number of aerosols= ', ia
    127   endif ! of is_master
    128 
    129 ! For the zero aerosol case, we currently make a dummy n2 aerosol which is zero everywhere.
    130 ! (See aeropacity.F90 for how this works). A better solution would be to turn off the
    131 ! aerosol machinery in the no aerosol case, but this would be complicated. LK
    132 
    133   if (ia.eq.0) then  !For the zero aerosol case.
    134      ia = 0
    135      noaero = .true.
    136   endif
    137 
    138   if (.not.noaero .and. ia.ne.naerkind) then
    139     if (is_master) then
    140       print*, 'Aerosols counted not equal to naerkind'
    141       print*, 'set correct value for nearkind in callphys.def'
    142       print*, 'which should be ',ia
    143       print*, 'according to current options in callphys.def'
    144       print*, 'or change/correct incompatible options there'
    145       print*, 'Abort in iniaerosol'
    146     endif
    147     call abort_physic("iniaerosl",'wrong number of aerosols',1)
    148   endif ! of if (ia.ne.naerkind)
    149 
    150   END SUBROUTINE iniaerosol
    151106
    152107end module aerosol_mod
Note: See TracChangeset for help on using the changeset viewer.