Changeset 3963 for trunk/LMDZ.TITAN
- Timestamp:
- Nov 18, 2025, 1:41:14 PM (4 days ago)
- File:
-
- 1 edited
-
trunk/LMDZ.TITAN/libf/phytitan/calc_condlh.F90 (modified) (7 diffs)
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.TITAN/libf/phytitan/calc_condlh.F90
r3497 r3963 24 24 ! ------ 25 25 ! B. de Batz de Trenquelléon (07/2023) 26 !27 ! !!! WARNING !!!28 ! ---------------29 ! We suppose a given order of tracers !30 ! Tracers come from microphysics (ices_indx) :31 ! CH4 --> 132 ! C2H2 --> 233 ! C2H6 --> 334 ! HCN --> 435 ! AC6H6--> 536 26 !==================================================================== 37 27 … … 56 46 real, parameter :: mmolC2H2 = 26.e-3 ! Molar mass of C2H2 [kg.mol-1] 57 47 real, parameter :: mmolC2H6 = 30.e-3 ! Molar mass of C2H6 [kg.mol-1] 48 real, parameter :: mmolC6H6 = 78.e-3 ! Molar mass of AC6H6 [kg.mol-1] 58 49 real, parameter :: mmolHCN = 27.e-3 ! Molar mass of HCN [kg.mol-1] 59 real, parameter :: mmol AC6H6 = 78.e-3 ! Molar mass of AC6H6[kg.mol-1]50 real, parameter :: mmolHC3N = 51.e-3 ! Molar mass of HC3N [kg.mol-1] 60 51 61 52 real, parameter :: TcCH4 = 190.56 ! Critical point of CH4 [K] 62 53 real, parameter :: TcC2H2 = 308.35 ! Critical point of C2H2 [K] 63 54 real, parameter :: TcC2H6 = 305.35 ! Critical point of C2H6 [K] 55 real, parameter :: TcC6H6 = 562.10 ! Critical point of AC6H6 [K] 64 56 real, parameter :: TcHCN = 456.70 ! Critical point of HCN [K] 65 real, parameter :: Tc AC6H6 = 562.10 ! Critical point of AC6H6[K]57 real, parameter :: TcHC3N = 527.00 ! Critical point of HC3N [K] 66 58 67 59 ! Local variables : … … 83 75 do ig = 1, ngrid ! Main loop on the horizontal grid 84 76 do l = 1, nlayer ! Loop on vertical layers 85 ftx = (1. - pt(ig,l) ) / TcCH477 ftx = (1. - pt(ig,l) / TcCH4) 86 78 if (ftx .le. 1.e-3) ftx = 1.e-3 87 79 Lcx(ig,l,iq) = 8.314 * TcCH4 * (7.08 * ftx**0.354 + 10.95 * 1.1e-2 * ftx**0.456) / mmolCH4 … … 96 88 do ig = 1, ngrid ! Main loop on the horizontal grid 97 89 do l = 1, nlayer ! Loop on vertical layers 98 ftx = (1. - pt(ig,l) ) / TcC2H290 ftx = (1. - pt(ig,l) / TcC2H2) 99 91 if (ftx .le. 1.e-3) ftx = 1.e-3 100 92 Lcx(ig,l,iq) = 8.314 * TcC2H2 * (7.08 * ftx**0.354 + 10.95 * 1.9e-1 * ftx**0.456) / mmolC2H2 … … 109 101 do ig = 1, ngrid ! Main loop on the horizontal grid 110 102 do l = 1, nlayer ! Loop on vertical layers 111 ftx = (1. - pt(ig,l) ) / TcC2H6103 ftx = (1. - pt(ig,l) / TcC2H6) 112 104 if (ftx .le. 1.e-3) ftx = 1.e-3 113 105 Lcx(ig,l,iq) = 8.314 * TcC2H6 * (7.08 * ftx**0.354 + 10.95 * 9.9e-2 * ftx**0.456) / mmolC2H6 … … 115 107 enddo ! End of main loop on the horizontal grid 116 108 117 ! AC6H6 :118 !------- -109 ! C6H6 : 110 !------- 119 111 else if(trim(nameOfTracer(gazs_indx(iq))) .eq. "AC6H6") then 120 112 121 113 do ig = 1, ngrid ! Main loop on the horizontal grid 122 114 do l = 1, nlayer ! Loop on vertical layers 123 ftx = (1. - pt(ig,l) ) / TcAC6H6115 ftx = (1. - pt(ig,l) / TcC6H6) 124 116 if (ftx .le. 1.e-3) ftx = 1.e-3 125 Lcx(ig,l,iq) = 8.314 * Tc AC6H6 * (7.08 * ftx**0.354 + 10.95 * 9.9e-2 * ftx**0.456) / mmolAC6H6117 Lcx(ig,l,iq) = 8.314 * TcC6H6 * (7.08 * ftx**0.354 + 10.95 * 9.9e-2 * ftx**0.456) / mmolC6H6 126 118 enddo ! End of loop on vertical layers 127 119 enddo ! End of main loop on the horizontal grid … … 133 125 do ig = 1, ngrid ! Main loop on the horizontal grid 134 126 do l = 1, nlayer ! Loop on vertical layers 135 ftx = (1. - pt(ig,l) ) / TcHCN127 ftx = (1. - pt(ig,l) / TcHCN) 136 128 if (ftx .le. 1.e-3) ftx = 1.e-3 137 129 Lcx(ig,l,iq) = 8.314 * TcHCN * (7.08 * ftx**0.354 + 10.95 * 3.88e-1 * ftx**0.456) / mmolHCN 130 enddo ! End of loop on vertical layers 131 enddo ! End of main loop on the horizontal grid 132 133 ! HC3N : 134 !------- 135 else if(trim(nameOfTracer(gazs_indx(iq))) .eq. "HC3N") then 136 137 do ig = 1, ngrid ! Main loop on the horizontal grid 138 do l = 1, nlayer ! Loop on vertical layers 139 ftx = (1. - pt(ig,l) / TcHC3N) 140 if (ftx .le. 1.e-3) ftx = 1.e-3 141 Lcx(ig,l,iq) = 8.314 * TcHC3N * (7.08 * ftx**0.354 + 10.95 * 3.88e-1 * ftx**0.456) / mmolHC3N 138 142 enddo ! End of loop on vertical layers 139 143 enddo ! End of main loop on the horizontal grid
Note: See TracChangeset
for help on using the changeset viewer.
