source: flexpart.git/src/hanna.f90 @ 3481cc1

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

move license from headers to a different file

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