[92fab65] | 1 | ! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt |
---|
| 2 | ! SPDX-License-Identifier: GPL-3.0-or-later |
---|
[332fbbd] | 3 | |
---|
[e200b7a] | 4 | subroutine hanna(z) |
---|
| 5 | ! i |
---|
| 6 | !***************************************************************************** |
---|
| 7 | ! * |
---|
| 8 | ! Computation of \sigma_i and \tau_L based on the scheme of Hanna (1982) * |
---|
| 9 | ! * |
---|
| 10 | ! Author: A. Stohl * |
---|
| 11 | ! * |
---|
| 12 | ! 4 December 1997 * |
---|
| 13 | ! * |
---|
| 14 | !***************************************************************************** |
---|
| 15 | ! * |
---|
| 16 | ! Variables: * |
---|
| 17 | ! dsigwdz [1/s] vertical gradient of sigw * |
---|
| 18 | ! ol [m] Obukhov length * |
---|
| 19 | ! sigu, sigv, sigw standard deviations of turbulent velocity fluctuations * |
---|
| 20 | ! tlu [s] Lagrangian time scale for the along wind component. * |
---|
| 21 | ! tlv [s] Lagrangian time scale for the cross wind component. * |
---|
| 22 | ! tlw [s] Lagrangian time scale for the vertical wind component. * |
---|
| 23 | ! ust, ustar [m/s] friction velocity * |
---|
| 24 | ! wst, wstar [m/s] convective velocity scale * |
---|
| 25 | ! * |
---|
| 26 | !***************************************************************************** |
---|
| 27 | |
---|
| 28 | use par_mod |
---|
| 29 | use com_mod |
---|
| 30 | use hanna_mod |
---|
| 31 | |
---|
| 32 | implicit none |
---|
| 33 | |
---|
| 34 | real :: corr,z |
---|
| 35 | |
---|
| 36 | |
---|
| 37 | |
---|
| 38 | !********************** |
---|
| 39 | ! 1. Neutral conditions |
---|
| 40 | !********************** |
---|
| 41 | |
---|
| 42 | if (h/abs(ol).lt.1.) then |
---|
| 43 | ust=max(1.e-4,ust) |
---|
| 44 | corr=z/ust |
---|
| 45 | sigu=1.e-2+2.0*ust*exp(-3.e-4*corr) |
---|
| 46 | sigw=1.3*ust*exp(-2.e-4*corr) |
---|
| 47 | dsigwdz=-2.e-4*sigw |
---|
| 48 | sigw=sigw+1.e-2 |
---|
| 49 | sigv=sigw |
---|
| 50 | tlu=0.5*z/sigw/(1.+1.5e-3*corr) |
---|
| 51 | tlv=tlu |
---|
| 52 | tlw=tlu |
---|
| 53 | |
---|
| 54 | |
---|
| 55 | !*********************** |
---|
| 56 | ! 2. Unstable conditions |
---|
| 57 | !*********************** |
---|
| 58 | |
---|
| 59 | else if (ol.lt.0.) then |
---|
| 60 | |
---|
| 61 | |
---|
| 62 | ! Determine sigmas |
---|
| 63 | !***************** |
---|
| 64 | |
---|
| 65 | sigu=1.e-2+ust*(12-0.5*h/ol)**0.33333 |
---|
| 66 | sigv=sigu |
---|
| 67 | sigw=sqrt(1.2*wst**2*(1.-.9*zeta)*zeta**0.66666+ & |
---|
| 68 | (1.8-1.4*zeta)*ust**2)+1.e-2 |
---|
| 69 | dsigwdz=0.5/sigw/h*(-1.4*ust**2+wst**2* & |
---|
| 70 | (0.8*max(zeta,1.e-3)**(-.33333)-1.8*zeta**0.66666)) |
---|
| 71 | |
---|
| 72 | |
---|
| 73 | ! Determine average Lagrangian time scale |
---|
| 74 | !**************************************** |
---|
| 75 | |
---|
| 76 | tlu=0.15*h/sigu |
---|
| 77 | tlv=tlu |
---|
| 78 | if (z.lt.abs(ol)) then |
---|
| 79 | tlw=0.1*z/(sigw*(0.55-0.38*abs(z/ol))) |
---|
| 80 | else if (zeta.lt.0.1) then |
---|
| 81 | tlw=0.59*z/sigw |
---|
| 82 | else |
---|
| 83 | tlw=0.15*h/sigw*(1.-exp(-5*zeta)) |
---|
| 84 | endif |
---|
| 85 | |
---|
| 86 | |
---|
| 87 | !********************* |
---|
| 88 | ! 3. Stable conditions |
---|
| 89 | !********************* |
---|
| 90 | |
---|
| 91 | else |
---|
| 92 | sigu=1.e-2+2.*ust*(1.-zeta) |
---|
| 93 | sigv=1.e-2+1.3*ust*(1.-zeta) |
---|
| 94 | sigw=sigv |
---|
| 95 | dsigwdz=-1.3*ust/h |
---|
| 96 | tlu=0.15*h/sigu*(sqrt(zeta)) |
---|
| 97 | tlv=0.467*tlu |
---|
| 98 | tlw=0.1*h/sigw*zeta**0.8 |
---|
| 99 | endif |
---|
| 100 | |
---|
| 101 | |
---|
| 102 | tlu=max(10.,tlu) |
---|
| 103 | tlv=max(10.,tlv) |
---|
| 104 | tlw=max(30.,tlw) |
---|
| 105 | |
---|
| 106 | if (dsigwdz.eq.0.) dsigwdz=1.e-10 |
---|
| 107 | |
---|
| 108 | end subroutine hanna |
---|