Ignore:
Timestamp:
Jul 15, 2025, 12:13:09 PM (5 weeks ago)
Author:
jbclement
Message:

Mars PCM 1D:
In the 1D model, if no input profile found, the 'dust_mass' and 'dust_number' tracers are initialized according to "Conrath dust" in line with what it's done for CO2 and HDO, instead of being set to 0 which caused bugs.
JBC

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/libf/phymars/dyn1d/read_profile_mod.F90

    r3734 r3848  
    1919!          - q(co2) = 0.95
    2020!          - q(hdo) = q(h2o)*2*155.76e-6*5
     21!          - q(dust_mass) and q(dust_number) initialized with Conrath dust
    2122!
    2223!     - igcm is not available at this part of the code
     
    3637!   4. Traitment for minor isotopologue if the major isotopologue input profile does not exist
    3738!=====================================================================================================================!
    38   subroutine read_profile(nb_tracer, nb_layer, qsurf, q)
    39 
    40   use infotrac, only: tname
    41   use tracer_mod, only: igcm_co2, igcm_h2o_vap, igcm_h2o_ice, &
    42                         igcm_dust_number, igcm_dust_mass,     &
    43                         igcm_ccn_number, igcm_ccn_mass
     39  subroutine read_profile(nb_tracer, nb_layer, qsurf, q, play)
     40
     41  use infotrac,       only: tname
     42  use tracer_mod,     only: igcm_co2, igcm_h2o_vap, igcm_h2o_ice, &
     43                            igcm_dust_number, igcm_dust_mass,     &
     44                            igcm_ccn_number, igcm_ccn_mass
     45  use aeropacity_mod, only: topdustref
     46  use dust_param_mod, only: odpref
     47  use comcstfi_h,     only: pi
     48
    4449  implicit none
    4550!---------------------------------------------------------------------------------------------------------------------!
     
    5055  integer, intent(in) :: &
    5156     nb_layer, & ! number of layer
    52      nb_tracer  ! number of traceur read from traceur.def
     57     nb_tracer   ! number of traceur read from traceur.def
     58  real, dimension(nb_layer), intent(in) :: play ! Pressure at the middle of the layers (Pa)
    5359!---------------------------------------------------------------------------------------------------------------------!
    5460!  Output arguments:
     
    6167!--------
    6268  integer :: &
    63      iq,                 & ! loop over nb_tracer
    64      ilayer,             & ! loop over nb_layer
    65      ierr,               & ! open file iostat
    66      indice_h2o_vap = 0, & ! indice of h2o_vap tracer
    67      indice_h2o_ice = 0, & ! indice of h2o_ice tracer
    68      indice_hdo_vap = 0, & ! indice of hdo_vap tracer
    69      indice_hdo_ice = 0    ! indice of hdo_ice tracer
     69     iq,                     & ! loop over nb_tracer
     70     ilayer,                 & ! loop over nb_layer
     71     ierr,                   & ! open file iostat
     72     indice_h2o_vap     = 0, & ! indice of h2o_vap tracer
     73     indice_h2o_ice     = 0, & ! indice of h2o_ice tracer
     74     indice_hdo_vap     = 0, & ! indice of hdo_vap tracer
     75     indice_hdo_ice     = 0, & ! indice of hdo_ice tracer
     76     indice_dust_mass   = 0, & ! indice of hdo_vap tracer
     77     indice_dust_number = 0    ! indice of hdo_ice tracer
    7078
    7179  character(len=80), dimension(nb_tracer) :: &
     
    7381
    7482  logical :: &
    75      hdo_vap = .false., & ! used to compute hdo_vap profile if its input profile is missing (= .true)
    76      hdo_ice = .false.    ! used to compute hdo_ice profile if its input profile is missing (= .true)
     83     hdo_vap = .false., & ! used to compute hdo_vap profile if its input profile is missing (= .true.)
     84     hdo_ice = .false., & ! used to compute hdo_ice profile if its input profile is missing (= .true.)
     85     dust = .false.       ! used to compute dust mass & number profiles if their input profiles are missing (= .true.)
     86
     87  real :: zp, varian, r3n_q, ref_r0, r0_lift, mass2number_lift ! Dust variables
     88  real, parameter ::    &
     89     nueff_lift = 0.5,  &
     90     reff_lift = 3.e-6, & ! Effective radius of lifted dust (m)
     91     rho_dust = 2500      ! Mars dust density (kg.m-3)
    7792!=====================================================================================================================!
    7893!=== BEGIN                                                                                                            !
     
    96111    else if (trim(name_tracer(iq)) == 'hdo_ice') then
    97112      indice_hdo_ice = iq
     113    else if (trim(name_tracer(iq)) == 'dust_mass') then
     114      indice_dust_mass = iq
     115    else if (trim(name_tracer(iq)) == 'dust_number') then
     116      indice_dust_number = iq
    98117    end if
    99118  end do
     
    126145!---------------------------------------------------------------------------------------------------------------------!
    127146    else
    128       write(*,*)"File profile_"//trim(name_tracer(iq))//" not found!"
    129       write(*,*)"q(traceur) and qsurf(traceur) are set to 0! Except:"
    130       write(*,*)"    - q(co2) = 0.95"
    131       write(*,*)"    - q(hdo) = q(h2o)*2*155.76e-6*5"
     147      write(*,*)"Warning: file profile_"//trim(name_tracer(iq))//" not found!"
    132148!---------------------------------------------------------------------------------------------------------------------!
    133149! 3.1.2.a. Some cases require a special initialization
     
    135151      select case(trim(name_tracer(iq)))
    136152        case("co2")
    137           q(1:nb_layer, iq) = 0.95
     153          q(1:nb_layer,iq) = 0.95
     154          write(*,*)"Initilization: q(co2) = 0.95"
    138155
    139156        case("hdo_vap")
    140157          hdo_vap = .true.
     158          write(*,*)"Initialization: q(hdo) = q(h2o)*2*155.76e-6*5"
    141159
    142160        case("hdo_ice")
    143161          hdo_ice = .true.
     162          write(*,*)"Initialization: q(hdo) = q(h2o)*2*155.76e-6*5"
     163
     164        case("dust_mass","dust_number")
     165          dust = .true.
     166          write(*,*)"q(dust_mass) and q(dust_number) initialized with Conrath dust"
     167
     168        case default
     169          write(*,*)"q("//trim(name_tracer(iq))//") and qsurf("//trim(name_tracer(iq))//") are set to 0!"
    144170
    145171       end select
     
    163189    end do
    164190  end if
     191
     192  if (dust) then
     193    varian = sqrt(log(1. + nueff_lift))
     194    r3n_q = exp(-4.5*varian**2)*(3./4.)/(pi*rho_dust)
     195    ref_r0 = exp(2.5*varian**2)
     196    r0_lift = reff_lift/ref_r0
     197    mass2number_lift = r3n_q/r0_lift**3
     198    do ilayer = 1,nb_layer
     199      zp = (odpref/play(ilayer))**(70./topdustref)
     200      q(ilayer,indice_dust_mass) = max(exp(0.007*(1. - max(zp,1.))),1.e-3)
     201      q(ilayer,indice_dust_number) = q(ilayer,indice_dust_mass)*mass2number_lift
     202    enddo
     203  endif
    165204!=====================================================================================================================!
    166205!=== END                                                                                                              !
Note: See TracChangeset for help on using the changeset viewer.