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