Ignore:
Timestamp:
Sep 5, 2025, 2:46:19 PM (4 months ago)
Author:
jbclement
Message:

Mars PCM 1D:
Prescription of the atmospheric water ice profile in 1D. The boolean 'ctrl_h2oice' activates the option, the file "profile_ref_h2ovap" defines the prescribed profile and the value 'relaxtime_h2oice', if positive, activates the relaxation and gives the constant of time.
JBC

File:
1 edited

Legend:

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

    r3838 r3912  
    66use phyredem,            only: physdem0, physdem1
    77use watersat_mod,        only: watersat
    8 use tracer_mod,          only: igcm_h2o_vap, noms
     8use tracer_mod,          only: igcm_h2o_vap, igcm_h2o_ice, noms
    99use comcstfi_h,          only: pi, g, rcp, cpp
    1010use time_phylmdz_mod,    only: daysec
     
    8383logical :: startfiles_1D, therestart1D, therestartfi, there
    8484
    85 ! Prescription of atmospheric water profile
    86 logical                         :: prescribed_h2ovap
    87 real                            :: h2ovap_relax_time
    88 real, dimension(nlayer)         :: q_prescribed_h2o_vap
     85! Prescription of atmospheric water/ice profiles
     86logical                         :: ctrl_h2ovap, ctrl_h2oice
     87real                            :: relaxtime_h2ovap, relaxtime_h2oice
     88real, dimension(nlayer)         :: qref_h2ovap, qref_h2oice
    8989real, dimension(:), allocatable :: zqsat
    9090!=======================================================================
     
    142142call init_testphys1d('start1D.txt','startfi.nc',therestart1D,therestartfi,ngrid,nlayer,odpref,nq,q, &
    143143                     time,psurf,u,v,temp,ndt,ptif,pks,dttestphys,zqsat,dq,dqdyn,day0,day,gru,grv,w, &
    144                      play,plev,latitude,longitude,cell_area,prescribed_h2ovap,h2ovap_relax_time,q_prescribed_h2o_vap)
     144                     play,plev,latitude,longitude,cell_area,                                        &
     145                     ctrl_h2ovap,relaxtime_h2ovap,qref_h2ovap,ctrl_h2oice,relaxtime_h2oice,qref_h2oice)
    145146
    146147! Write a "startfi" file
     
    208209    ! Water tracer increment: specific for 1D
    209210    ! ---------------------------------------
    210     ! The physics computes the tendency on q(1,:,igcm_h2o_vap), here we can mimic an effect of the "3D dynamics"
    211     ! either by forcing a profile or by relaxing towards a prescribed profile
    212     if (water .and. prescribed_h2ovap) then ! If atmospheric water profile is prescribed
    213         if (h2ovap_relax_time < 0.) then ! Forcing
    214             ! For some tests: unless it reaches saturation (maximal value)
    215             !call watersat(nlayer,temp,play,zqsat)
    216             !dq(1,:,igcm_h2o_vap) = (min(zqsat(:),q_prescribed_h2o_vap(:)) - q(1,:,igcm_h2o_vap))/dttestphys
    217             dq(1,:,igcm_h2o_vap) = (q_prescribed_h2o_vap(:) - q(1,:,igcm_h2o_vap))/dttestphys
    218         else ! Relaxation
    219             dq(1,:,igcm_h2o_vap) = dq(1,:,igcm_h2o_vap) - (q(1,:,igcm_h2o_vap) - q_prescribed_h2o_vap(:))/h2ovap_relax_time
     211    ! The physics computes the tendency on q for igcm_h2o_vap/igcm_h2o_ice, here we can mimic an effect
     212    ! of the "3D dynamics" either by forcing a profile or by relaxing towards a prescribed profile
     213    if (water) then
     214        if (ctrl_h2ovap) then ! If atmospheric water vapour profile is prescribed
     215            if (relaxtime_h2ovap < 0.) then ! Forcing
     216                ! For some tests: unless it reaches saturation (maximal value)
     217                !call watersat(nlayer,temp,play,zqsat)
     218                !dq(1,:,igcm_h2o_vap) = (min(zqsat(:),qref_h2ovap(:)) - q(1,:,igcm_h2o_vap))/dttestphys
     219                dq(1,:,igcm_h2o_vap) = (qref_h2ovap(:) - q(1,:,igcm_h2o_vap))/dttestphys
     220            else ! Relaxation
     221                dq(1,:,igcm_h2o_vap) = dq(1,:,igcm_h2o_vap) - (q(1,:,igcm_h2o_vap) - qref_h2ovap(:))/relaxtime_h2ovap
     222            endif
     223        endif
     224        if (ctrl_h2oice) then ! If atmospheric water ice profile is prescribed
     225            if (relaxtime_h2oice < 0.) then ! Forcing
     226                dq(1,:,igcm_h2o_ice) = (qref_h2oice(:) - q(1,:,igcm_h2o_ice))/dttestphys
     227            else ! Relaxation
     228                dq(1,:,igcm_h2o_ice) = dq(1,:,igcm_h2o_ice) - (q(1,:,igcm_h2o_ice) - qref_h2oice(:))/relaxtime_h2oice
     229            endif
    220230        endif
    221231     endif
Note: See TracChangeset for help on using the changeset viewer.