Changeset 4047 for trunk


Ignore:
Timestamp:
Feb 5, 2026, 5:50:37 PM (3 weeks ago)
Author:
aslmd
Message:

Titan PCM:

  • add possibility to release latent heat to the atmosphere due to gases condensation.
  • correct sign error on temperature tendency due to latent heat release in phytitan/cond_muphy.F90
  • add latent heat release flag in callphys.def ("latent_heat" : .true. or .false.).

EMo

Location:
trunk/LMDZ.TITAN
Files:
5 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.TITAN/deftank/callphys.def.allmuparams

    r3657 r4047  
    104104# If yes, number of ices ? (must be compatible with traceur.def AND microphysical model)
    105105nices      = 1
     106# Apply condensation heating rate?
     107latent_heat = .true.
    106108# Use new optics for clouds ?
    107109opt4clouds = .true.
     
    137139rho_aer                = 600.
    138140# Enable/disable Haze production process
    139 haze_production        = T
     141haze_production        = .true.
    140142# Enable/disable Haze coagulation process
    141 haze_coagulation       = T
     143haze_coagulation       = .true.
    142144# Coagulation interactions, a combination of:
    143145#    0 - no interactions (same as haze_coagulation == F)
     
    148150haze_coag_interactions = 7
    149151# Enable/disable Haze sedimentation process
    150 haze_sedimentation     = T
     152haze_sedimentation     = .true.
    151153# Disable Fiadero correction for sedimentation process
    152 no_fiadero             = T
     154no_fiadero             = .true.
    153155# Fiadero correction minimum ratio threshold
    154156fiadero_min_ratio      = 0.1
     
    156158fiadero_max_ratio      = 10.
    157159# Force settling velocity to M0
    158 wsed_m0                = T
     160wsed_m0                = .true.
    159161# Force settling velocity to M3
    160 wsed_m3                = F
     162wsed_m3                = .false.
    161163# Enable/disable clouds sedimentation process
    162164# (automatically set to F if clouds microphysics is not enabled)
    163 clouds_sedimentation   = T
     165clouds_sedimentation   = .true.
    164166# Enable/disable clouds nucleation and condensation processes
    165167# (automatically set to F if clouds microphysics is not enabled)
    166 clouds_nuc_cond        = T
     168clouds_nuc_cond        = .true.
    167169# Condensible species configuration file
    168170# (not needed if clouds microphysics is not enabled)
    169171specie_cfg             = mp2m_species.cfg
    170172# Enable/disable spherical mode transfert probability
    171 transfert_probability = T
     173transfert_probability = .true.
    172174# Spherical mode transfert probability look-up tables file
    173175# (optional if 'transfert_probability' is False)
     
    175177# Electric charging coagulation correction
    176178# If set to .false. then no correction is assumed
    177 electric_charging     = T
     179electric_charging     = .true.
    178180# File for the electric charging correction factor.
    179181# (optional if 'electric_charging' is False)
  • trunk/LMDZ.TITAN/libf/phytitan/callkeys_mod.F90

    r3657 r4047  
    2121      logical,save :: callchim, callmufi, callclouds
    2222!$OMP THREADPRIVATE(callchim,callmufi,callclouds)
     23      logical,save :: latent_heat
     24!$OMP THREADPRIVATE(latent_heat)
    2325      logical,save :: global1d
    2426!$OMP THREADPRIVATE(global1d)
  • trunk/LMDZ.TITAN/libf/phytitan/cond_muphy.F90

    r3090 r4047  
    6969! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    7070do iq = 1, size(ices_indx)
    71     dtlc(:,:) = dtlc(:,:) + (dqmuficond(:,:,iq) * Lc(:,:,iq))
     71    dtlc(:,:) = dtlc(:,:) - (dqmuficond(:,:,iq) * Lc(:,:,iq))
    7272enddo
    7373
     
    7575! Condensation heating rate :
    7676! ~~~~~~~~~~~~~~~~~~~~~~~~~~~
    77 ! If ice formation  : dqmuficond > 0 --> dtlc > 0
    78 ! Else vaporisation : dqmuficond < 0 --> dtlc < 0
     77! If ice formation  : dqmuficond < 0 --> dtlc > 0
     78! Else vaporisation : dqmuficond > 0 --> dtlc < 0
    7979dtlc(:,:) = dtlc(:,:) / cpp ! [K.s-1]
    8080
  • trunk/LMDZ.TITAN/libf/phytitan/inifis_mod.F90

    r3657 r4047  
    490490     endif
    491491
     492     write(*,*) "Apply condensation heating rate?"
     493     latent_heat=.false. ! default value
     494     call getin_p("latent_heat",latent_heat)
     495     write(*,*)" latent_heat = ",latent_heat
     496
    492497     write(*,*) "Disable the coupling of microphysics within rad. transf. ?"
    493498     write(*,*) "If disabled we will assume a planetwide vert. profile of extinction ..."
  • trunk/LMDZ.TITAN/libf/phytitan/physiq_mod.F90

    r4046 r4047  
    196196      real,intent(in) :: pv(ngrid,nlayer)      ! Meridional wind component (m/s).
    197197      real,intent(in) :: pt(ngrid,nlayer)      ! Temperature (K).
    198       real,intent(in) :: pq(ngrid,nlayer,nq)   ! Tracers (kg/kg_of_air).
     198      real,intent(in) :: pq(ngrid,nlayer,nq)   ! Tracers (X/kg_of_air). (X is kg or nothing (number of particles))
    199199      real,intent(in) :: flxw(ngrid,nlayer)    ! Vertical mass flux (ks/s) at lower boundary of layer
    200200
     
    431431
    432432        ! Allocate saved arrays (except for 1D model, where this has already been done)
     433        ! For mesoscale it is done in the interface
    433434#ifndef MESOSCALE
    434435        if (ngrid>1) call phys_state_var_init(nq)
     
    11771178               enddo
    11781179               call cond_muphy(ngrid,nlayer,pt,dmuficond,zdtlc)
    1179                !pdt(:,:) = pdt(:,:) + zdtlc(:,:)
     1180               if (latent_heat) then
     1181                  pdt(:,:) = pdt(:,:) + zdtlc(:,:)
     1182               endif
    11801183            endif
    11811184         endif ! callmufi
Note: See TracChangeset for help on using the changeset viewer.