source: branches/jerome/src_flexwrf_v3.1/hanna1.f90 @ 16

Last change on this file since 16 was 16, checked in by jebri, 11 years ago

sources for flexwrf v3.1

File size: 5.3 KB
Line 
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
22subroutine hanna1(z, &
23  ust,wst,ol,h,zeta,sigu,sigv,tlu,tlv,tlw, &
24  sigw,dsigwdz,dsigw2dz)
25
26  !                  i
27  !*****************************************************************************
28  !                                                                            *
29  !   Computation of \sigma_i and \tau_L based on the scheme of Hanna (1982)   *
30  !                                                                            *
31  !   Author: A. Stohl                                                         *
32  !                                                                            *
33  !   4 December 1997                                                          *
34  !                                                                            *
35  !*****************************************************************************
36  !                                                                            *
37  ! Variables:                                                                 *
38  ! dsigwdz [1/s]     vertical gradient of sigw                                *
39  ! ol [m]            Obukhov length                                           *
40  ! sigu, sigv, sigw  standard deviations of turbulent velocity fluctuations   *
41  ! tlu [s]           Lagrangian time scale for the along wind component.      *
42  ! tlv [s]           Lagrangian time scale for the cross wind component.      *
43  ! tlw [s]           Lagrangian time scale for the vertical wind component.   *
44  ! ust, ustar [m/s]  friction velocity                                        *
45  ! wst, wstar [m/s]  convective velocity scale                                *
46  !                                                                            *
47  !*****************************************************************************
48
49  use par_mod
50  use com_mod
51!  use hanna_mod
52
53  implicit none
54
55  real :: z,s1,s2
56  real :: ust,wst,ol,h,zeta,sigu,sigv,tlu,tlv,tlw
57  real :: sigw,dsigwdz,dsigw2dz
58
59
60
61  !**********************
62  ! 1. Neutral conditions
63  !**********************
64
65  if (h/abs(ol).lt.1.) then
66
67    ust=max(1.e-4,ust)
68    sigu=2.0*ust*exp(-3.e-4*z/ust)
69    sigu=max(sigu,1.e-5)
70    sigv=1.3*ust*exp(-2.e-4*z/ust)
71    sigv=max(sigv,1.e-5)
72    sigw=sigv
73    dsigw2dz=-6.76e-4*ust*exp(-4.e-4*z/ust)
74    tlu=0.5*z/sigw/(1.+1.5e-3*z/ust)
75    tlv=tlu
76    tlw=tlu
77
78
79  !***********************
80  ! 2. Unstable conditions
81  !***********************
82
83  else if (ol.lt.0.) then
84
85
86  ! Determine sigmas
87  !*****************
88
89    sigu=ust*(12-0.5*h/ol)**0.33333
90    sigu=max(sigu,1.e-6)
91    sigv=sigu
92
93    if (zeta.lt.0.03) then
94      sigw=0.96*wst*(3*zeta-ol/h)**0.33333
95      dsigw2dz=1.8432*wst*wst/h*(3*zeta-ol/h)**(-0.33333)
96    else if (zeta.lt.0.4) then
97      s1=0.96*(3*zeta-ol/h)**0.33333
98      s2=0.763*zeta**0.175
99      if (s1.lt.s2) then
100        sigw=wst*s1
101        dsigw2dz=1.8432*wst*wst/h*(3*zeta-ol/h)**(-0.33333)
102      else
103        sigw=wst*s2
104        dsigw2dz=0.203759*wst*wst/h*zeta**(-0.65)
105      endif
106    else if (zeta.lt.0.96) then
107      sigw=0.722*wst*(1-zeta)**0.207
108      dsigw2dz=-.215812*wst*wst/h*(1-zeta)**(-0.586)
109    else if (zeta.lt.1.00) then
110      sigw=0.37*wst
111      dsigw2dz=0.
112    endif
113    sigw=max(sigw,1.e-6)
114
115
116  ! Determine average Lagrangian time scale
117  !****************************************
118
119    tlu=0.15*h/sigu
120    tlv=tlu
121    if (z.lt.abs(ol)) then
122      tlw=0.1*z/(sigw*(0.55-0.38*abs(z/ol)))
123    else if (zeta.lt.0.1) then
124      tlw=0.59*z/sigw
125    else
126      tlw=0.15*h/sigw*(1.-exp(-5*zeta))
127    endif
128
129
130  !*********************
131  ! 3. Stable conditions
132  !*********************
133
134  else
135    sigu=2.*ust*(1.-zeta)
136    sigv=1.3*ust*(1.-zeta)
137    sigu=max(sigu,1.e-6)
138    sigv=max(sigv,1.e-6)
139    sigw=sigv
140    dsigw2dz=3.38*ust*ust*(zeta-1.)/h
141    tlu=0.15*h/sigu*(sqrt(zeta))
142    tlv=0.467*tlu
143    tlw=0.1*h/sigw*zeta**0.8
144  endif
145
146
147
148
149  tlu=max(10.,tlu)
150  tlv=max(10.,tlv)
151  tlw=max(30.,tlw)
152
153
154end subroutine hanna1
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG