| [3083] | 1 | subroutine calc_condlh(ngrid,nlayer,pt,Lcx) |
|---|
| 2 | |
|---|
| 3 | use tracer_h |
|---|
| 4 | |
|---|
| 5 | implicit none |
|---|
| 6 | |
|---|
| 7 | !==================================================================== |
|---|
| 8 | ! |
|---|
| 9 | ! Purpose |
|---|
| 10 | ! ------- |
|---|
| 11 | ! Computes latent heat of condensation [J.kg-1] for ice tracers |
|---|
| 12 | ! |
|---|
| 13 | ! INPUT |
|---|
| 14 | ! ----- |
|---|
| 15 | ! ngrid = Number of grid points in the chunk [-] |
|---|
| 16 | ! nlayer = Number of vertical layers [-] |
|---|
| 17 | ! pt = Temperature [K] |
|---|
| 18 | ! |
|---|
| 19 | ! OUTPUT |
|---|
| 20 | ! ------ |
|---|
| 21 | ! Lcx = Condensation latente heat [J.kg-1] |
|---|
| 22 | ! |
|---|
| [3090] | 23 | ! Author |
|---|
| 24 | ! ------ |
|---|
| [3083] | 25 | ! B. de Batz de Trenquelléon (07/2023) |
|---|
| 26 | ! |
|---|
| 27 | ! !!! WARNING !!! |
|---|
| 28 | ! --------------- |
|---|
| 29 | ! We suppose a given order of tracers ! |
|---|
| [3090] | 30 | ! Tracers come from microphysics (ices_indx) : |
|---|
| [3083] | 31 | ! CH4 --> 1 |
|---|
| 32 | ! C2H2 --> 2 |
|---|
| 33 | ! C2H6 --> 3 |
|---|
| [3090] | 34 | ! HCN --> 4 |
|---|
| [3083] | 35 | !==================================================================== |
|---|
| 36 | |
|---|
| 37 | |
|---|
| 38 | !------------------------------------ |
|---|
| 39 | ! 0. DECLARATIONS |
|---|
| 40 | !------------------------------------ |
|---|
| 41 | |
|---|
| 42 | ! Inputs : |
|---|
| 43 | !--------- |
|---|
| 44 | integer, intent(in) :: ngrid |
|---|
| 45 | integer, intent(in) :: nlayer |
|---|
| 46 | real, intent(in) :: pt(ngrid,nlayer) |
|---|
| 47 | |
|---|
| 48 | ! Outputs : |
|---|
| 49 | !---------- |
|---|
| 50 | real, intent(out) :: Lcx(ngrid,nlayer,size(ices_indx)) |
|---|
| 51 | |
|---|
| 52 | ! Parameters : |
|---|
| 53 | !------------- |
|---|
| 54 | real, parameter :: mmolCH4 = 16.e-3 ! Molar mass of CH4 [kg.mol-1] |
|---|
| 55 | real, parameter :: mmolC2H2 = 26.e-3 ! Molar mass of C2H2 [kg.mol-1] |
|---|
| 56 | real, parameter :: mmolC2H6 = 30.e-3 ! Molar mass of C2H6 [kg.mol-1] |
|---|
| 57 | real, parameter :: mmolHCN = 27.e-3 ! Molar mass of HCN [kg.mol-1] |
|---|
| 58 | |
|---|
| 59 | real, parameter :: TcCH4 = 190.56 ! Critical point of CH4 [K] |
|---|
| 60 | real, parameter :: TcC2H2 = 308.35 ! Critical point of C2H2 [K] |
|---|
| 61 | real, parameter :: TcC2H6 = 305.35 ! Critical point of C2H6 [K] |
|---|
| 62 | real, parameter :: TcHCN = 456.70 ! Critical point of HCN [K] |
|---|
| 63 | |
|---|
| 64 | ! Local variables : |
|---|
| 65 | !------------------ |
|---|
| 66 | integer :: iq, ig, l |
|---|
| 67 | real*8 :: ftx ! Variables for latente heat [-] |
|---|
| 68 | |
|---|
| 69 | |
|---|
| 70 | !------------------------------------ |
|---|
| 71 | ! 1. INITIALISATION |
|---|
| 72 | !------------------------------------ |
|---|
| 73 | |
|---|
| 74 | do iq = 1, size(ices_indx) ! For each species |
|---|
| 75 | |
|---|
| 76 | ! CH4 : |
|---|
| 77 | !------ |
|---|
| 78 | if(trim(nameOfTracer(gazs_indx(iq))) .eq. "CH4") then |
|---|
| 79 | |
|---|
| 80 | do ig = 1, ngrid ! Main loop on the horizontal grid |
|---|
| 81 | do l = 1, nlayer ! Loop on vertical layers |
|---|
| 82 | ftx = (1. - pt(ig,l)) / TcCH4 |
|---|
| 83 | if (ftx .le. 1.e-3) ftx = 1.e-3 |
|---|
| 84 | Lcx(ig,l,iq) = 8.314 * TcCH4 * (7.08 * ftx**0.354 + 10.95 * 1.1e-2 * ftx**0.456) / mmolCH4 |
|---|
| 85 | enddo ! End of loop on vertical layers |
|---|
| 86 | enddo ! End of main loop on the horizontal grid |
|---|
| 87 | |
|---|
| 88 | |
|---|
| 89 | ! C2H2 : |
|---|
| 90 | !------- |
|---|
| 91 | else if(trim(nameOfTracer(gazs_indx(iq))) .eq. "C2H2") then |
|---|
| 92 | |
|---|
| 93 | do ig = 1, ngrid ! Main loop on the horizontal grid |
|---|
| 94 | do l = 1, nlayer ! Loop on vertical layers |
|---|
| 95 | ftx = (1. - pt(ig,l)) / TcC2H2 |
|---|
| 96 | if (ftx .le. 1.e-3) ftx = 1.e-3 |
|---|
| 97 | Lcx(ig,l,iq) = 8.314 * TcC2H2 * (7.08 * ftx**0.354 + 10.95 * 1.9e-1 * ftx**0.456) / mmolC2H2 |
|---|
| 98 | enddo ! End of loop on vertical layers |
|---|
| 99 | enddo ! End of main loop on the horizontal grid |
|---|
| 100 | |
|---|
| 101 | |
|---|
| 102 | ! C2H6 : |
|---|
| 103 | !------- |
|---|
| 104 | else if(trim(nameOfTracer(gazs_indx(iq))) .eq. "C2H6") then |
|---|
| 105 | |
|---|
| 106 | do ig = 1, ngrid ! Main loop on the horizontal grid |
|---|
| 107 | do l = 1, nlayer ! Loop on vertical layers |
|---|
| 108 | ftx = (1. - pt(ig,l)) / TcC2H6 |
|---|
| 109 | if (ftx .le. 1.e-3) ftx = 1.e-3 |
|---|
| 110 | Lcx(ig,l,iq) = 8.314 * TcC2H6 * (7.08 * ftx**0.354 + 10.95 * 9.9e-2 * ftx**0.456) / mmolC2H6 |
|---|
| 111 | enddo ! End of loop on vertical layers |
|---|
| 112 | enddo ! End of main loop on the horizontal grid |
|---|
| 113 | |
|---|
| 114 | |
|---|
| 115 | ! HCN : |
|---|
| 116 | !------ |
|---|
| 117 | else if(trim(nameOfTracer(gazs_indx(iq))) .eq. "HCN") then |
|---|
| 118 | |
|---|
| 119 | do ig = 1, ngrid ! Main loop on the horizontal grid |
|---|
| 120 | do l = 1, nlayer ! Loop on vertical layers |
|---|
| 121 | ftx = (1. - pt(ig,l)) / TcHCN |
|---|
| 122 | if (ftx .le. 1.e-3) ftx = 1.e-3 |
|---|
| 123 | Lcx(ig,l,iq) = 8.314 * TcHCN * (7.08 * ftx**0.354 + 10.95 * 3.88e-1 * ftx**0.456) / mmolHCN |
|---|
| 124 | enddo ! End of loop on vertical layers |
|---|
| 125 | enddo ! End of main loop on the horizontal grid |
|---|
| 126 | |
|---|
| 127 | |
|---|
| 128 | else |
|---|
| 129 | write(*,*) 'Species ',trim(nameOfTracer(gazs_indx(iq))),' not considered... Check traceur.def' |
|---|
| 130 | Lcx(ig,l,iq) = 0.0D0 |
|---|
| 131 | endif |
|---|
| 132 | |
|---|
| 133 | enddo ! End of species |
|---|
| 134 | |
|---|
| 135 | return |
|---|
| 136 | |
|---|
| 137 | end subroutine calc_condlh |
|---|