[92fab65] | 1 | ! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt |
---|
| 2 | ! SPDX-License-Identifier: GPL-3.0-or-later |
---|
[332fbbd] | 3 | |
---|
[8a65cb0] | 4 | subroutine cbl(wp,zp,ust,wst,h,rhoa,rhograd,sigmaw,dsigmawdz,tlw,ptot,Q,phi,ath,bth,ol,flagrein) |
---|
| 5 | ! i i i i i i i i i i o o o o o i o |
---|
| 6 | |
---|
| 7 | use par_mod, only:pi |
---|
| 8 | use com_mod, only:ldirect |
---|
| 9 | |
---|
| 10 | implicit none |
---|
| 11 | !======================================================================================================================================================= |
---|
| 12 | !=============== CBL skewed vertical profiles and formulation of LHH 1996 with profile of w3 from LHB 2000 ======== |
---|
| 13 | !=============== LHH formulation has been modified to account for variable density profiles and backward in time or forward in time simulations ======== |
---|
| 14 | !=============== see Cassiani et al. BLM 2014 doi for explanations and references ======== |
---|
| 15 | !======================================================================================================================================================= |
---|
| 16 | |
---|
| 17 | real :: usurad2,usurad2p,C0,costluar4,eps |
---|
| 18 | parameter (usurad2=0.7071067812,usurad2p=0.3989422804,C0=3,costluar4=0.66667,eps=0.000001) |
---|
| 19 | |
---|
| 20 | |
---|
| 21 | integer :: flagrein |
---|
| 22 | real :: wp,zp,ust,wst,h,dens,ddens,sigmaw,dsigmawdz,tlw,rhoa,rhograd |
---|
| 23 | real ::fluarw,fluarw2 |
---|
| 24 | real ::w3,w2 |
---|
| 25 | real ::dw3,dw2 |
---|
| 26 | real ::wb,wa |
---|
| 27 | real ::deltawa,deltawb |
---|
| 28 | real ::wold,wold2,wold_z |
---|
| 29 | real ::pa,pb,alfa |
---|
| 30 | real ::Phi,Q,ptot |
---|
| 31 | real :: timedir |
---|
| 32 | real ::cuberoot |
---|
| 33 | real ::z0,ol,transition |
---|
| 34 | |
---|
| 35 | |
---|
| 36 | real :: & |
---|
| 37 | erf, & |
---|
| 38 | aperfa, & |
---|
| 39 | aperfb, & |
---|
| 40 | ath, & |
---|
| 41 | bth |
---|
| 42 | |
---|
| 43 | real :: & |
---|
| 44 | pow, & |
---|
| 45 | z, & |
---|
| 46 | skew, & |
---|
| 47 | skew2, & |
---|
| 48 | radw2, & |
---|
| 49 | rluarw, & |
---|
| 50 | xluarw, & |
---|
| 51 | aluarw, & |
---|
| 52 | bluarw, & |
---|
| 53 | sigmawa, & |
---|
| 54 | sigmawb, & |
---|
| 55 | dskew, & |
---|
| 56 | dradw2, & |
---|
| 57 | dfluarw, & |
---|
| 58 | drluarw, & |
---|
| 59 | dxluarw, & |
---|
| 60 | daluarw, & |
---|
| 61 | dbluarw, & |
---|
| 62 | dsigmawa, & |
---|
| 63 | dsigmawb, & |
---|
| 64 | dwa, & |
---|
| 65 | dwb, & |
---|
| 66 | sigmawa2, & |
---|
| 67 | sigmawb2 |
---|
| 68 | |
---|
| 69 | |
---|
| 70 | dens=rhoa !put to 1 for test constant density simulation |
---|
| 71 | ddens=rhograd !put to 0 for test constant density simulation |
---|
| 72 | |
---|
| 73 | |
---|
| 74 | timedir=ldirect !ldirect contains direction of time forward (1) or backward(-1) |
---|
| 75 | !========================= assegnazione z =========================================================== |
---|
| 76 | z=(zp/h) |
---|
| 77 | |
---|
| 78 | !================== stability transition function see Cassiani et al(2015) BLM ====================== |
---|
| 79 | transition=1. |
---|
| 80 | !if (ol.lt.-50) transition=((sin(((ol+100.)/100.)*pi))-1.)/2. |
---|
| 81 | if (-h/ol.lt.15) transition=((sin((((-h/ol)+10.)/10.)*pi)))/2.+0.5 |
---|
| 82 | |
---|
| 83 | !========================= momento secondo ========================================================== |
---|
| 84 | |
---|
| 85 | w2=(sigmaw*sigmaw) |
---|
| 86 | dw2=(2.*sigmaw*dsigmawdz) |
---|
| 87 | |
---|
| 88 | |
---|
| 89 | !=================== dissipazione ======================================= |
---|
| 90 | |
---|
| 91 | alfa=2.*w2/(C0*tlw) |
---|
| 92 | |
---|
| 93 | !======================================================================== |
---|
| 94 | |
---|
| 95 | wold=timedir*wp |
---|
| 96 | |
---|
| 97 | !========================================================================= |
---|
| 98 | !------------------------------ momento terzo ============================ |
---|
| 99 | |
---|
| 100 | w3=((1.2*z*((1.-z)**(3./2.)))+eps)*(wst**3)*transition |
---|
| 101 | dw3=(1.2*(((1.-z)**(3./2.))+z*1.5*((1.-z)**(1./2.))*(-1.)))*(wst**3)*(1./h)*transition |
---|
| 102 | |
---|
| 103 | !===========================================================================0 |
---|
| 104 | |
---|
| 105 | skew=w3/(w2**1.5) |
---|
| 106 | skew2=skew*skew |
---|
| 107 | dskew=(dw3*w2**(1.5)-w3*1.5*w2**0.5*dw2)/w2**3 |
---|
| 108 | radw2=w2**0.5 |
---|
| 109 | dradw2=0.5*w2**(-0.5)*dw2 |
---|
| 110 | !costluar4=0.66667 ! costante da LHH |
---|
[23547f3] | 111 | fluarw=costluar4*(cuberoot(skew)) !skew**(1./3.) |
---|
[8a65cb0] | 112 | fluarw2=fluarw*fluarw |
---|
| 113 | |
---|
| 114 | if (skew.ne.0) then |
---|
| 115 | |
---|
| 116 | dfluarw=costluar4*(1./3.)*cuberoot(skew**(-2.))*dskew |
---|
| 117 | |
---|
| 118 | rluarw=(1.+fluarw2)**3.*skew2/((3.+fluarw2)**2.*fluarw2) !-> r |
---|
| 119 | xluarw=(1.+fluarw2)**1.5*skew/((3.+fluarw2)*fluarw) !----> r^1/2 |
---|
| 120 | |
---|
| 121 | drluarw=( ((3.*(1.+fluarw2)**2*(2.*fluarw*dfluarw)*skew2)+ & |
---|
| 122 | (1.+fluarw2)**3*2.*skew*dskew) *(3.+fluarw2)**2.*fluarw2 - & |
---|
| 123 | (1.+fluarw2)**3*skew2* & |
---|
| 124 | ( (2.*(3.+fluarw2)*(2.*fluarw*dfluarw)*fluarw2) + & |
---|
| 125 | (3.+fluarw2)**2*2.*fluarw*dfluarw) )/ & |
---|
| 126 | (((3.+fluarw2)**2.*fluarw2)**2) |
---|
| 127 | |
---|
| 128 | dxluarw=( ((1.5*(1.+fluarw2)**0.5*(2.*fluarw*dfluarw)*skew)+ & |
---|
| 129 | (1.+fluarw2)**1.5*dskew) *(3.+fluarw2)*fluarw - & |
---|
| 130 | (1.+fluarw2)**1.5*skew* & |
---|
| 131 | (3.*dfluarw+3*fluarw2*dfluarw) )/ & |
---|
| 132 | (((3.+fluarw2)*fluarw)**2) |
---|
| 133 | |
---|
| 134 | else |
---|
| 135 | dfluarw=0. |
---|
| 136 | rluarw=0. |
---|
| 137 | drluarw=0. |
---|
| 138 | xluarw=0. |
---|
| 139 | dxluarw=0. |
---|
| 140 | end if |
---|
| 141 | |
---|
| 142 | |
---|
| 143 | |
---|
| 144 | aluarw=0.5*(1.-xluarw/(4.+rluarw)**0.5) |
---|
| 145 | bluarw=1.-aluarw |
---|
| 146 | |
---|
| 147 | daluarw=-0.5*( (dxluarw*(4.+rluarw)**0.5) - & |
---|
| 148 | (0.5*xluarw*(4.+rluarw)**(-0.5)*drluarw) ) & |
---|
| 149 | /(4.+rluarw) |
---|
| 150 | dbluarw=-daluarw |
---|
| 151 | |
---|
| 152 | sigmawa=radw2*(bluarw/(aluarw*(1.+fluarw2)))**0.5 |
---|
| 153 | sigmawb=radw2*(aluarw/(bluarw*(1.+fluarw2)))**0.5 |
---|
| 154 | |
---|
| 155 | dsigmawa=dradw2*(bluarw/(aluarw*(1.+fluarw2)))**0.5+ & |
---|
| 156 | radw2*( & |
---|
| 157 | (0.5*(bluarw/(aluarw*(1.+fluarw2)))**(-0.5)) * & |
---|
| 158 | ( & |
---|
| 159 | (dbluarw*(aluarw*(1.+fluarw2))- & |
---|
| 160 | bluarw*(daluarw*(1.+fluarw2)+aluarw*2.*fluarw*dfluarw)) & |
---|
| 161 | /((aluarw*(1.+fluarw2))**2) & |
---|
| 162 | ) & |
---|
| 163 | ) |
---|
| 164 | dsigmawb=dradw2*(aluarw/(bluarw*(1.+fluarw2)))**0.5+ & |
---|
| 165 | radw2*( & |
---|
| 166 | (0.5*(aluarw/(bluarw*(1.+fluarw2)))**(-0.5)) * & |
---|
| 167 | ( & |
---|
| 168 | (daluarw*(bluarw*(1.+fluarw2))- & |
---|
| 169 | aluarw*(dbluarw*(1.+fluarw2)+bluarw*2.*fluarw*dfluarw)) & |
---|
| 170 | /((bluarw*(1.+fluarw2))**2) & |
---|
| 171 | ) & |
---|
| 172 | ) |
---|
| 173 | |
---|
| 174 | wa=(fluarw*sigmawa) |
---|
| 175 | wb=(fluarw*sigmawb) |
---|
| 176 | dwa=dfluarw*sigmawa+fluarw*dsigmawa |
---|
| 177 | dwb=dfluarw*sigmawb+fluarw*dsigmawb |
---|
| 178 | |
---|
| 179 | deltawa=wold-wa |
---|
| 180 | deltawb=wold+wb |
---|
| 181 | wold2=wold*wold |
---|
| 182 | sigmawa2=sigmawa*sigmawa |
---|
| 183 | sigmawb2=sigmawb*sigmawb |
---|
| 184 | |
---|
| 185 | if (abs(deltawa).gt.6.*sigmawa.and.abs(deltawb).gt.6.*sigmawb) flagrein=1 |
---|
| 186 | |
---|
| 187 | pa=(usurad2p*(1./sigmawa))*(exp(-(0.5*((deltawa/sigmawa)**2.)))) |
---|
| 188 | pb=(usurad2p*(1./sigmawb))*(exp(-(0.5*((deltawb/sigmawb)**2.)))) |
---|
[23547f3] | 189 | |
---|
[8a65cb0] | 190 | ptot=dens*aluarw*pa+dens*bluarw*pb |
---|
| 191 | |
---|
| 192 | aperfa=deltawa*usurad2/sigmawa |
---|
| 193 | aperfb=deltawb*usurad2/sigmawb |
---|
| 194 | |
---|
| 195 | Phi=-0.5* & |
---|
| 196 | (aluarw*dens*dwa+dens*wa*daluarw+aluarw*wa*ddens)*erf(aperfa) & |
---|
| 197 | +sigmawa*(aluarw*dens*dsigmawa*(wold2/sigmawa2+1.)+ & |
---|
| 198 | sigmawa*dens*daluarw+sigmawa*ddens*aluarw+ & |
---|
| 199 | aluarw*wold*dens/sigmawa2*(sigmawa*dwa-wa*dsigmawa))*pa & |
---|
| 200 | +0.5* & |
---|
| 201 | (bluarw*dens*dwb+wb*dens*dbluarw+wb*bluarw*ddens)*erf(aperfb) & |
---|
| 202 | +sigmawb*(bluarw*dens*dsigmawb*(wold2/sigmawb2+1.)+ & |
---|
| 203 | sigmawb*dens*dbluarw+sigmawb*ddens*bluarw+ & |
---|
| 204 | bluarw*wold*dens/sigmawb2*(-sigmawb*dwb+wb*dsigmawb))*pb |
---|
| 205 | |
---|
| 206 | Q=timedir*((aluarw*dens*deltawa/sigmawa2)*pa+ & |
---|
| 207 | (bluarw*dens*deltawb/sigmawb2)*pb) |
---|
| 208 | |
---|
| 209 | ath=(1./ptot)*(-(C0/2.)*alfa*Q+phi) |
---|
| 210 | bth=sqrt(C0*alfa) |
---|
| 211 | !bth=sngl(sigmaw*sqrt(2.*tlw)) |
---|
| 212 | |
---|
| 213 | return |
---|
| 214 | end |
---|
| 215 | |
---|
| 216 | |
---|
| 217 | |
---|
| 218 | |
---|
| 219 | |
---|
| 220 | FUNCTION CUBEROOT (X) RESULT (Y) |
---|
| 221 | |
---|
| 222 | IMPLICIT NONE |
---|
| 223 | |
---|
| 224 | real, INTENT(IN) :: X |
---|
| 225 | real:: Y |
---|
| 226 | |
---|
| 227 | real, PARAMETER :: THIRD = 0.333333333 |
---|
| 228 | |
---|
| 229 | |
---|
| 230 | Y = SIGN((ABS(X))**THIRD, X) |
---|
| 231 | |
---|
| 232 | RETURN |
---|
| 233 | |
---|
| 234 | END FUNCTION CUBEROOT |
---|
| 235 | |
---|
| 236 | |
---|
| 237 | |
---|
| 238 | |
---|
| 239 | FUNCTION CUBEROOTD (X) RESULT (Y) |
---|
| 240 | |
---|
| 241 | IMPLICIT NONE |
---|
| 242 | |
---|
| 243 | DOUBLE PRECISION, INTENT(IN) :: X |
---|
| 244 | DOUBLE PRECISION :: Y |
---|
| 245 | |
---|
| 246 | DOUBLE PRECISION, PARAMETER :: THIRD = 0.33333333333333333333333333333333333333333333333333333333333333333333333333333333333D0 |
---|
| 247 | |
---|
| 248 | |
---|
| 249 | Y = SIGN((ABS(X))**THIRD, X) |
---|
| 250 | |
---|
| 251 | RETURN |
---|
| 252 | |
---|
| 253 | END FUNCTION CUBEROOTD |
---|