CTBTO FLEXPART WO4 (2015-10-15)
 All Classes Files Functions Variables
hanna.f90
Go to the documentation of this file.
1 !**********************************************************************
2 ! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010 *
3 ! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa, *
4 ! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann *
5 ! *
6 ! This file is part of FLEXPART. *
7 ! *
8 ! FLEXPART is free software: you can redistribute it and/or modify *
9 ! it under the terms of the GNU General Public License as published by*
10 ! the Free Software Foundation, either version 3 of the License, or *
11 ! (at your option) any later version. *
12 ! *
13 ! FLEXPART is distributed in the hope that it will be useful, *
14 ! but WITHOUT ANY WARRANTY; without even the implied warranty of *
15 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
16 ! GNU General Public License for more details. *
17 ! *
18 ! You should have received a copy of the GNU General Public License *
19 ! along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. *
20 !**********************************************************************
21 
22 subroutine hanna(z)
23  ! i
24  !*****************************************************************************
25  ! *
26  ! Computation of \sigma_i and \tau_L based on the scheme of Hanna (1982) *
27  ! *
28  ! Author: A. Stohl *
29  ! *
30  ! 4 December 1997 *
31  ! *
32  !*****************************************************************************
33  ! *
34  ! Variables: *
35  ! dsigwdz [1/s] vertical gradient of sigw *
36  ! ol [m] Obukhov length *
37  ! sigu, sigv, sigw standard deviations of turbulent velocity fluctuations *
38  ! tlu [s] Lagrangian time scale for the along wind component. *
39  ! tlv [s] Lagrangian time scale for the cross wind component. *
40  ! tlw [s] Lagrangian time scale for the vertical wind component. *
41  ! ust, ustar [m/s] friction velocity *
42  ! wst, wstar [m/s] convective velocity scale *
43  ! *
44  !*****************************************************************************
45 
46  use par_mod
47  use com_mod
48  use hanna_mod
49 
50  implicit none
51 
52  real :: corr,z
53 
54 
55 
56  !**********************
57  ! 1. Neutral conditions
58  !**********************
59 
60  if (h/abs(ol).lt.1.) then
61  ust=max(1.e-4,ust)
62  corr=z/ust
63  sigu=1.e-2+2.0*ust*exp(-3.e-4*corr)
64  sigw=1.3*ust*exp(-2.e-4*corr)
65  dsigwdz=-2.e-4*sigw
66  sigw=sigw+1.e-2
67  sigv=sigw
68  tlu=0.5*z/sigw/(1.+1.5e-3*corr)
69  tlv=tlu
70  tlw=tlu
71 
72 
73  !***********************
74  ! 2. Unstable conditions
75  !***********************
76 
77  else if (ol.lt.0.) then
78 
79 
80  ! Determine sigmas
81  !*****************
82 
83  sigu=1.e-2+ust*(12-0.5*h/ol)**0.33333
84  sigv=sigu
85  sigw=sqrt(1.2*wst**2*(1.-.9*zeta)*zeta**0.66666+ &
86  (1.8-1.4*zeta)*ust**2)+1.e-2
87  dsigwdz=0.5/sigw/h*(-1.4*ust**2+wst**2* &
88  (0.8*max(zeta,1.e-3)**(-.33333)-1.8*zeta**0.66666))
89 
90 
91  ! Determine average Lagrangian time scale
92  !****************************************
93 
94  tlu=0.15*h/sigu
95  tlv=tlu
96  if (z.lt.abs(ol)) then
97  tlw=0.1*z/(sigw*(0.55-0.38*abs(z/ol)))
98  else if (zeta.lt.0.1) then
99  tlw=0.59*z/sigw
100  else
101  tlw=0.15*h/sigw*(1.-exp(-5*zeta))
102  endif
103 
104 
105  !*********************
106  ! 3. Stable conditions
107  !*********************
108 
109  else
110  sigu=1.e-2+2.*ust*(1.-zeta)
111  sigv=1.e-2+1.3*ust*(1.-zeta)
112  sigw=sigv
113  dsigwdz=-1.3*ust/h
114  tlu=0.15*h/sigu*(sqrt(zeta))
115  tlv=0.467*tlu
116  tlw=0.1*h/sigw*zeta**0.8
117  endif
118 
119 
120  tlu=max(10.,tlu)
121  tlv=max(10.,tlv)
122  tlw=max(30.,tlw)
123 
124  if (dsigwdz.eq.0.) dsigwdz=1.e-10
125 
126 end subroutine hanna
subroutine hanna(z)
Definition: hanna.f90:22