Ignore:
Timestamp:
Feb 23, 2024, 10:22:20 AM (9 months ago)
Author:
gmilcareck
Message:

Add the possibility to use a fixed vertical molar fraction profile for the
collision-induced absorption like the correlated-k.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phystd/dyn1d/kcm1d.F90

    r3105 r3233  
    88  use comsaison_h, only: mu0, fract, dist_star
    99  use planete_mod
    10   use callkeys_mod, only: pceil, tstrat, tracer, global1d
     10  use callkeys_mod, only: pceil, tstrat, tracer, global1d, &
     11                          varspec, varspec_data, nvarlayer
    1112  use inifis_mod, only: inifis
    1213  use aerosol_mod, only: iniaerosol
     
    2021  use geometry_mod, only: init_geometry
    2122  use dimphy, only : init_dimphy
     23  use gases_h, only: ngasmx
    2224
    2325  implicit none
     
    8789  real int_dtauv(1,llm,L_NSPECTV)
    8890  real Eatmtot
     91 
     92  ! For fixed variable molar mass
     93  real, dimension(:),allocatable,save :: p_var
     94  real, dimension(:),allocatable,save :: mu_var
     95  real, dimension(:,:),allocatable,save :: frac_var
    8996
    9097  integer ierr
     
    116123  ALLOCATE(fract(1))
    117124  ALLOCATE(glat(1))
     125
     126! Initialize fixed variable mu
     127!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~   
     128
     129  if(varspec) then
     130    IF (.NOT.ALLOCATED(p_var))  ALLOCATE(p_var(nvarlayer))
     131    IF (.NOT.ALLOCATED(mu_var))  ALLOCATE(mu_var(nvarlayer))
     132    IF (.NOT.ALLOCATED(frac_var))  ALLOCATE(frac_var(nvarlayer,ngasmx))
     133    p_var = 0.0
     134    mu_var = 0.0
     135    frac_var = 0.0
     136    dt_file=TRIM(varspec_data)
     137    open(34,file=dt_file,form='formatted',status='old',iostat=ios)
     138    if (ios.ne.0) then        ! file not found
     139      write(*,*) 'Error from varspec'
     140      write(*,*) 'Data file ',trim(varspec_data),' not found.'
     141      write(*,*) 'Check that the data is in your run repository.'
     142      call abort_physic !a verifier
     143    else
     144      do k=1,nvarlayer
     145        read(34,*) p_var(k), mu_var(k),frac_var(k,1:ngasmx)
     146        !The order of columns in frac_var must correspond to order of molecules gases.def
     147        !The format of your file must be:
     148        ! pressure(k) molar_mass(k), molar_fraction_of_gas_1(k), molar_fraction_of_gas_2(k),..., molar_fraction_of_gas_ngasmx(k)
     149      enddo
     150    endif
     151    close(34)
     152  endif
    118153
    119154  !  Load parameters from "run.def"
     
    386421          int_dtaui,int_dtauv,                            &
    387422          tau_col,cloudfrac,totcloudfrac,                 &
    388           .false.,firstcall,lastcall)
     423          .false.,p_var,frac_var,firstcall,lastcall)
    389424
    390425     !write(*,*) 'BASE 3'
     
    438473       int_dtaui,int_dtauv,                                    &
    439474       tau_col,cloudfrac,totcloudfrac,                         &
    440        .false.,firstcall,lastcall)
     475       .false.,p_var,frac_var,firstcall,lastcall)
    441476
    442477
Note: See TracChangeset for help on using the changeset viewer.