Ignore:
Timestamp:
Nov 23, 2022, 4:09:32 PM (2 years ago)
Author:
abierjon
Message:

Mars GCM:
The program util/aeroptical.F90 can now compute column-integrated optical depths
of the aerosols, in addition to the opacity profiles. The user has to specify it
('yes' or 'no') in aeroptical.def
Detailed changes :

  • update of the aeroptical.def file to ask the user if the column optical depth should be computed (non-retrocompatible change)
  • add in the init2 subroutine a computation of the layers' height delta_z, with adaptable method depending on the availability of some variables in the input file and on wether or not the input file has been zrecasted before. Preliminary validation shows that these different methods yield a +/-3% error on the output column optical depth tau_[aer]. The delta_z variable is also written in the output file.
  • add a log file for warnings in the interpolation subroutine in aeropt_mod.F90 + add some comments in the code
  • add the ouput "zzlev" (= interlayer altitudes) in libf/phymars/physiq_mod.F so that it can be directly used by aeroptical.F90 to compute delta_z

AB

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/util/aeropt_mod.F90

    r2817 r2830  
    1212 real,dimension(:,:),save,allocatable :: Qext    ! Extinction efficiency coefficient [/]
    1313 real,dimension(:,:),save,allocatable :: omeg    ! Single scattering albedo [/]
     14 
     15 integer,save :: aeroptlogfileID = 100           ! Log file's unit ID
    1416
    1517 CONTAINS
     
    4648  integer :: read_ok                         ! to know if the line is well read
    4749  integer :: iwvl,isize                      ! Wavelength and Particle size loop indices
     50 
     51  !==============================================================================
     52  ! aeropt_mod module variables used here:
     53  !==============================================================================
     54  ! nwvl,nsize,wvl,radiusdyn,Qext,omeg
    4855
    4956  !==========================================================================
     
    228235  real,dimension(2) :: reffQext         ! Qext after reff interpolation
    229236  real,dimension(2) :: reffomeg         ! omeg after reff interpolation
     237 
     238  !==============================================================================
     239  ! aeropt_mod module variables used here:
     240  !==============================================================================
     241  ! nwvl,nsize,wvl,radiusdyn,Qext,omeg,aeroptlogfileID
    230242
    231243  !==========================================================================
     
    262274  ENDDO
    263275
    264   if ((radiusdyn(nsize).lt.reff_val).or.(radiusdyn(1).gt.reff_val)) then
    265     write(*,*) "WARNING: the effective radius (",reff_val,"m) is out of the bounds"
    266     write(*,*) "of the optical properties' tables"
    267     write(*,*) "(reff_inf=",radiusdyn(1),"m ; reff_sup=",radiusdyn(nsize),"m)"
    268     write(*,*) "A missing value will be written for opacity."
    269 
    270     ! NB: this should especially handle cases when reff=0
     276  if (reff_val.eq.missval) then
     277    ! if the effective radius is a missing value,
     278    ! put a missing value for the optical properties
     279    Qext_val = missval
     280    omeg_val = missval
     281!    write(*,*) "reff missing value detected"
     282   
     283  else if ((radiusdyn(nsize).lt.reff_val).or.(radiusdyn(1).gt.reff_val)) then
     284    write(aeroptlogfileID,*) "WARNING: the effective radius (",reff_val,"m) is out of the bounds"
     285    write(aeroptlogfileID,*) "of the optical properties' tables"
     286    write(aeroptlogfileID,*) "(reff_inf=",radiusdyn(1),"m ; reff_sup=",radiusdyn(nsize),"m)"
     287    write(aeroptlogfileID,*) "A missing value will be written for optical properties."
     288
    271289    Qext_val = missval
    272290    omeg_val = missval
     
    333351
    334352  implicit none
     353 
     354  !==============================================================================
     355  ! aeropt_mod module variables used here:
     356  !==============================================================================
     357  ! wvl,radiusdyn,Qext,omeg
    335358
    336359  if (allocated(wvl)) deallocate(wvl)
     
    342365 
    343366!*******************************************************************************
     367 SUBROUTINE create_logfile(filename)
     368  !==============================================================================
     369  ! Purpose:
     370  ! Deallocate module variables
     371  !==============================================================================
     372
     373  implicit none
     374 
     375  !==============================================================================
     376  ! Arguments:
     377  !==============================================================================
     378  character(len=128),intent(in) :: filename ! name of the aeropt log file
     379  !==============================================================================
     380  ! aeropt_mod module variables used here:
     381  !==============================================================================
     382  ! aeroptlogfileID
     383 
     384  ! Create output file
     385  write(*,*) "Creating aeropt log file: ",trim(filename)
     386  write(*,*) "It will overwrite any existing file with the same name."
     387  OPEN(UNIT=aeroptlogfileID,FILE=trim(filename),FORM='formatted',STATUS='replace')
     388  ! NB: STATUS='replace' will overwrite any existing file with the same name
     389 
     390 END SUBROUTINE create_logfile
     391 
     392!*******************************************************************************
    344393
    345394END MODULE aeropt_mod
Note: See TracChangeset for help on using the changeset viewer.