Ignore:
Timestamp:
Nov 7, 2024, 11:13:46 AM (2 weeks ago)
Author:
debatzbr
Message:

Add AC6H6 condensation in the microphysics

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.TITAN/libf/phytitan/calc_condlh.F90

    r3090 r3497  
    3333!   C2H6 --> 3
    3434!   HCN  --> 4
     35!   AC6H6--> 5
    3536!====================================================================   
    3637
     
    5657real, parameter :: mmolC2H6 = 30.e-3 ! Molar mass of C2H6 [kg.mol-1]
    5758real, parameter :: mmolHCN  = 27.e-3 ! Molar mass of HCN [kg.mol-1]
     59real, parameter :: mmolAC6H6 = 78.e-3 ! Molar mass of AC6H6 [kg.mol-1]
    5860
    5961real, parameter :: TcCH4  = 190.56 ! Critical point of CH4  [K]
     
    6163real, parameter :: TcC2H6 = 305.35 ! Critical point of C2H6 [K]
    6264real, parameter :: TcHCN  = 456.70 ! Critical point of HCN [K]
     65real, parameter :: TcAC6H6 = 562.10 ! Critical point of AC6H6 [K]
    6366
    6467! Local variables :
     
    111114            enddo ! End of loop on vertical layers
    112115        enddo ! End of main loop on the horizontal grid
    113 
     116   
     117    ! AC6H6 :
     118    !--------
     119    else if(trim(nameOfTracer(gazs_indx(iq))) .eq. "AC6H6") then
     120       
     121        do ig = 1, ngrid ! Main loop on the horizontal grid
     122            do l = 1, nlayer ! Loop on vertical layers
     123                ftx = (1. - pt(ig,l)) / TcAC6H6
     124                if (ftx .le. 1.e-3) ftx = 1.e-3
     125                Lcx(ig,l,iq) = 8.314 * TcAC6H6 * (7.08 * ftx**0.354 + 10.95 * 9.9e-2 * ftx**0.456) / mmolAC6H6
     126            enddo ! End of loop on vertical layers
     127        enddo ! End of main loop on the horizontal grid
    114128
    115129    ! HCN :
     
    127141
    128142    else
    129         write(*,*) 'Species ',trim(nameOfTracer(gazs_indx(iq))),' not considered... Check traceur.def'
     143        write(*,*) 'Species ',trim(nameOfTracer(gazs_indx(iq))),' not considered in calc_condlh.F90... Check traceur.def'
    130144        Lcx(ig,l,iq) = 0.0D0
    131145    endif
Note: See TracChangeset for help on using the changeset viewer.