source: flexpart.git/src/hanna1.f90

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file was 92fab65, checked in by Ignacio Pisso <ip@…>, 4 years ago

add SPDX-License-Identifier to all .f90 files

  • Property mode set to 100644
File size: 3.8 KB
Line 
1! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
2! SPDX-License-Identifier: GPL-3.0-or-later
3
4subroutine hanna1(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 :: z,s1,s2
35
36
37
38  !**********************
39  ! 1. Neutral conditions
40  !**********************
41
42  if (h/abs(ol).lt.1.) then
43
44    ust=max(1.e-4,ust)
45    sigu=2.0*ust*exp(-3.e-4*z/ust)
46    sigu=max(sigu,1.e-5)
47    sigv=1.3*ust*exp(-2.e-4*z/ust)
48    sigv=max(sigv,1.e-5)
49    sigw=sigv
50    dsigw2dz=-6.76e-4*ust*exp(-4.e-4*z/ust)
51    tlu=0.5*z/sigw/(1.+1.5e-3*z/ust)
52    tlv=tlu
53    tlw=tlu
54
55
56  !***********************
57  ! 2. Unstable conditions
58  !***********************
59
60  else if (ol.lt.0.) then
61
62
63  ! Determine sigmas
64  !*****************
65
66    sigu=ust*(12-0.5*h/ol)**0.33333
67    sigu=max(sigu,1.e-6)
68    sigv=sigu
69
70    if (zeta.lt.0.03) then
71      sigw=0.96*wst*(3*zeta-ol/h)**0.33333
72      dsigw2dz=1.8432*wst*wst/h*(3*zeta-ol/h)**(-0.33333)
73    else if (zeta.lt.0.4) then
74      s1=0.96*(3*zeta-ol/h)**0.33333
75      s2=0.763*zeta**0.175
76      if (s1.lt.s2) then
77        sigw=wst*s1
78        dsigw2dz=1.8432*wst*wst/h*(3*zeta-ol/h)**(-0.33333)
79      else
80        sigw=wst*s2
81        dsigw2dz=0.203759*wst*wst/h*zeta**(-0.65)
82      endif
83    else if (zeta.lt.0.96) then
84      sigw=0.722*wst*(1-zeta)**0.207
85      dsigw2dz=-.215812*wst*wst/h*(1-zeta)**(-0.586)
86    else if (zeta.lt.1.00) then
87      sigw=0.37*wst
88      dsigw2dz=0.
89    endif
90    sigw=max(sigw,1.e-6)
91
92
93  ! Determine average Lagrangian time scale
94  !****************************************
95
96    tlu=0.15*h/sigu
97    tlv=tlu
98    if (z.lt.abs(ol)) then
99      tlw=0.1*z/(sigw*(0.55-0.38*abs(z/ol)))
100    else if (zeta.lt.0.1) then
101      tlw=0.59*z/sigw
102    else
103      tlw=0.15*h/sigw*(1.-exp(-5*zeta))
104    endif
105
106
107  !*********************
108  ! 3. Stable conditions
109  !*********************
110
111  else
112    sigu=2.*ust*(1.-zeta)
113    sigv=1.3*ust*(1.-zeta)
114    sigu=max(sigu,1.e-6)
115    sigv=max(sigv,1.e-6)
116    sigw=sigv
117    dsigw2dz=3.38*ust*ust*(zeta-1.)/h
118    tlu=0.15*h/sigu*(sqrt(zeta))
119    tlv=0.467*tlu
120    tlw=0.1*h/sigw*zeta**0.8
121  endif
122
123
124
125
126  tlu=max(10.,tlu)
127  tlv=max(10.,tlv)
128  tlw=max(30.,tlw)
129
130
131end subroutine hanna1
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG