source: branches/jerome/src_flexwrf_v3.1/hanna.f90

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

sources for flexwrf v3.1

File size: 4.8 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 hanna(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 :: corr,z
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    ust=max(1.e-4,ust)
67    corr=z/ust
68    sigu=1.e-2+2.0*ust*exp(-3.e-4*corr)
69    sigw=1.3*ust*exp(-2.e-4*corr)
70    dsigwdz=-2.e-4*sigw
71    sigw=sigw+1.e-2
72    sigv=sigw
73    tlu=0.5*z/sigw/(1.+1.5e-3*corr)
74    tlv=tlu
75    tlw=tlu
76
77
78  !***********************
79  ! 2. Unstable conditions
80  !***********************
81
82  else if (ol.lt.0.) then
83
84
85  ! Determine sigmas
86  !*****************
87
88    sigu=1.e-2+ust*(12-0.5*h/ol)**0.33333
89    sigv=sigu
90    sigw=sqrt(1.2*wst**2*(1.-.9*zeta)*zeta**0.66666+ &
91         (1.8-1.4*zeta)*ust**2)+1.e-2
92    dsigwdz=0.5/sigw/h*(-1.4*ust**2+wst**2* &
93         (0.8*max(zeta,1.e-3)**(-.33333)-1.8*zeta**0.66666))
94
95
96  ! Determine average Lagrangian time scale
97  !****************************************
98
99    tlu=0.15*h/sigu
100    tlv=tlu
101    if (z.lt.abs(ol)) then
102      tlw=0.1*z/(sigw*(0.55-0.38*abs(z/ol)))
103    else if (zeta.lt.0.1) then
104      tlw=0.59*z/sigw
105    else
106      tlw=0.15*h/sigw*(1.-exp(-5*zeta))
107    endif
108
109
110  !*********************
111  ! 3. Stable conditions
112  !*********************
113
114  else
115    sigu=1.e-2+2.*ust*(1.-zeta)
116    sigv=1.e-2+1.3*ust*(1.-zeta)
117    sigw=sigv
118    dsigwdz=-1.3*ust/h
119    tlu=0.15*h/sigu*(sqrt(zeta))
120    tlv=0.467*tlu
121    tlw=0.1*h/sigw*zeta**0.8
122  endif
123
124
125  tlu=max(10.,tlu)
126  tlv=max(10.,tlv)
127  tlw=max(30.,tlw)
128
129  if (dsigwdz.eq.0.) dsigwdz=1.e-10
130
131end subroutine hanna
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG