subroutine calc_condlh(ngrid,nlayer,pt,Lcx) use tracer_h implicit none !==================================================================== ! ! Purpose ! ------- ! Computes latent heat of condensation [J.kg-1] for ice tracers ! ! INPUT ! ----- ! ngrid = Number of grid points in the chunk [-] ! nlayer = Number of vertical layers [-] ! pt = Temperature [K] ! ! OUTPUT ! ------ ! Lcx = Condensation latente heat [J.kg-1] ! ! Author ! ------ ! B. de Batz de Trenquelléon (07/2023) ! ! !!! WARNING !!! ! --------------- ! We suppose a given order of tracers ! ! Tracers come from microphysics (ices_indx) : ! CH4 --> 1 ! C2H2 --> 2 ! C2H6 --> 3 ! HCN --> 4 ! AC6H6--> 5 !==================================================================== !------------------------------------ ! 0. DECLARATIONS !------------------------------------ ! Inputs : !--------- integer, intent(in) :: ngrid integer, intent(in) :: nlayer real, intent(in) :: pt(ngrid,nlayer) ! Outputs : !---------- real, intent(out) :: Lcx(ngrid,nlayer,size(ices_indx)) ! Parameters : !------------- real, parameter :: mmolCH4 = 16.e-3 ! Molar mass of CH4 [kg.mol-1] real, parameter :: mmolC2H2 = 26.e-3 ! Molar mass of C2H2 [kg.mol-1] real, parameter :: mmolC2H6 = 30.e-3 ! Molar mass of C2H6 [kg.mol-1] real, parameter :: mmolHCN = 27.e-3 ! Molar mass of HCN [kg.mol-1] real, parameter :: mmolAC6H6 = 78.e-3 ! Molar mass of AC6H6 [kg.mol-1] real, parameter :: TcCH4 = 190.56 ! Critical point of CH4 [K] real, parameter :: TcC2H2 = 308.35 ! Critical point of C2H2 [K] real, parameter :: TcC2H6 = 305.35 ! Critical point of C2H6 [K] real, parameter :: TcHCN = 456.70 ! Critical point of HCN [K] real, parameter :: TcAC6H6 = 562.10 ! Critical point of AC6H6 [K] ! Local variables : !------------------ integer :: iq, ig, l real*8 :: ftx ! Variables for latente heat [-] !------------------------------------ ! 1. INITIALISATION !------------------------------------ do iq = 1, size(ices_indx) ! For each species ! CH4 : !------ if(trim(nameOfTracer(gazs_indx(iq))) .eq. "CH4") then do ig = 1, ngrid ! Main loop on the horizontal grid do l = 1, nlayer ! Loop on vertical layers ftx = (1. - pt(ig,l)) / TcCH4 if (ftx .le. 1.e-3) ftx = 1.e-3 Lcx(ig,l,iq) = 8.314 * TcCH4 * (7.08 * ftx**0.354 + 10.95 * 1.1e-2 * ftx**0.456) / mmolCH4 enddo ! End of loop on vertical layers enddo ! End of main loop on the horizontal grid ! C2H2 : !------- else if(trim(nameOfTracer(gazs_indx(iq))) .eq. "C2H2") then do ig = 1, ngrid ! Main loop on the horizontal grid do l = 1, nlayer ! Loop on vertical layers ftx = (1. - pt(ig,l)) / TcC2H2 if (ftx .le. 1.e-3) ftx = 1.e-3 Lcx(ig,l,iq) = 8.314 * TcC2H2 * (7.08 * ftx**0.354 + 10.95 * 1.9e-1 * ftx**0.456) / mmolC2H2 enddo ! End of loop on vertical layers enddo ! End of main loop on the horizontal grid ! C2H6 : !------- else if(trim(nameOfTracer(gazs_indx(iq))) .eq. "C2H6") then do ig = 1, ngrid ! Main loop on the horizontal grid do l = 1, nlayer ! Loop on vertical layers ftx = (1. - pt(ig,l)) / TcC2H6 if (ftx .le. 1.e-3) ftx = 1.e-3 Lcx(ig,l,iq) = 8.314 * TcC2H6 * (7.08 * ftx**0.354 + 10.95 * 9.9e-2 * ftx**0.456) / mmolC2H6 enddo ! End of loop on vertical layers enddo ! End of main loop on the horizontal grid ! AC6H6 : !-------- else if(trim(nameOfTracer(gazs_indx(iq))) .eq. "AC6H6") then do ig = 1, ngrid ! Main loop on the horizontal grid do l = 1, nlayer ! Loop on vertical layers ftx = (1. - pt(ig,l)) / TcAC6H6 if (ftx .le. 1.e-3) ftx = 1.e-3 Lcx(ig,l,iq) = 8.314 * TcAC6H6 * (7.08 * ftx**0.354 + 10.95 * 9.9e-2 * ftx**0.456) / mmolAC6H6 enddo ! End of loop on vertical layers enddo ! End of main loop on the horizontal grid ! HCN : !------ else if(trim(nameOfTracer(gazs_indx(iq))) .eq. "HCN") then do ig = 1, ngrid ! Main loop on the horizontal grid do l = 1, nlayer ! Loop on vertical layers ftx = (1. - pt(ig,l)) / TcHCN if (ftx .le. 1.e-3) ftx = 1.e-3 Lcx(ig,l,iq) = 8.314 * TcHCN * (7.08 * ftx**0.354 + 10.95 * 3.88e-1 * ftx**0.456) / mmolHCN enddo ! End of loop on vertical layers enddo ! End of main loop on the horizontal grid else write(*,*) 'Species ',trim(nameOfTracer(gazs_indx(iq))),' not considered in calc_condlh.F90... Check traceur.def' Lcx(ig,l,iq) = 0.0D0 endif enddo ! End of species return end subroutine calc_condlh