source: flexpart.git/src/richardson.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: 7.5 KB
Line 
1subroutine richardson(psurf,ust,ttlev,qvlev,ulev,vlev,nuvz, &
2       akz,bkz,hf,tt2,td2,h,wst,hmixplus,metdata_format)
3  !                        i    i    i     i    i    i    i
4  ! i   i  i   i   i  o  o     o
5  !****************************************************************************
6  !                                                                           *
7  !     Calculation of mixing height based on the critical Richardson number. *
8  !     Calculation of convective time scale.                                 *
9  !     For unstable conditions, one iteration is performed. An excess        *
10  !     temperature (dependent on hf and wst) is calculated, added to the     *
11  !     temperature at the lowest model level. Then the procedure is repeated.*
12  !                                                                           *
13  !     Author: A. Stohl                                                      *
14  !                                                                           *
15  !     22 August 1996                                                        *
16  !                                                                           *
17  !     Literature:                                                           *
18  !     Vogelezang DHP and Holtslag AAM (1996): Evaluation and model impacts  *
19  !     of alternative boundary-layer height formulations. Boundary-Layer     *
20  !     Meteor. 81, 245-269.                                                  *
21  !                                                                           *
22  !****************************************************************************
23  !                                                                           *
24  !     Update: 1999-02-01 by G. Wotawa                                       *
25  !                                                                           *
26  !     Two meter level (temperature, humidity) is taken as reference level   *
27  !     instead of first model level.                                         *
28  !     New input variables tt2, td2 introduced.                              *
29  !                                                                           *
30  !     CHANGE: 17/11/2005 Caroline Forster NCEP GFS version                  *
31  !                                                                           *
32  !     Unified ECMWF and GFS builds                                          *
33  !     Marian Harustak, 12.5.2017                                            *
34  !       - Merged richardson and richardson_gfs into one routine using       *
35  !         if-then for meteo-type dependent code                             *
36  !                                                                           *
37  !****************************************************************************
38  !                                                                           *
39  ! Variables:                                                                *
40  ! h                          mixing height [m]                              *
41  ! hf                         sensible heat flux                             *
42  ! psurf                      surface pressure at point (xt,yt) [Pa]         *
43  ! tv                         virtual temperature                            *
44  ! wst                        convective velocity scale                      *
45  ! metdata_format             format of metdata (ecmwf/gfs)                  *
46  !                                                                           *
47  ! Constants:                                                                *
48  ! ric                        critical Richardson number                     *
49  !                                                                           *
50  !****************************************************************************
51
52  use par_mod
53  use class_gribfile
54
55  implicit none
56
57  integer :: metdata_format
58  integer :: i,k,nuvz,iter,llev,loop_start
59  real :: tv,tvold,zref,z,zold,pint,pold,theta,thetaref,ri
60  real :: akz(nuvz),bkz(nuvz),ulev(nuvz),vlev(nuvz),hf,wst,tt2,td2,ew
61  real :: psurf,ust,ttlev(nuvz),qvlev(nuvz),h,excess
62  real :: thetaold,zl,ul,vl,thetal,ril,hmixplus,wspeed,bvfsq,bvf
63  real :: f_qvsat,rh,rhold,rhl,theta1,theta2,zl1,zl2,thetam
64  real,parameter    :: const=r_air/ga, ric=0.25, b=100., bs=8.5
65  integer,parameter :: itmax=3
66
67  excess=0.0
68  iter=0
69
70  if (metdata_format.eq.GRIBFILE_CENTRE_NCEP) then
71    ! NCEP version: find first model level above ground
72    !**************************************************
73
74     llev = 0
75     do i=1,nuvz
76       if (psurf.lt.akz(i)) llev=i
77     end do
78     llev = llev+1
79    ! sec llev should not be 1!
80     if (llev.eq.1) llev = 2
81     if (llev.gt.nuvz) llev = nuvz-1
82    ! NCEP version
83  end if
84
85
86  ! Compute virtual temperature and virtual potential temperature at
87  ! reference level (2 m)
88  !*****************************************************************
89
9030   iter=iter+1
91
92  pold=psurf
93  tvold=tt2*(1.+0.378*ew(td2)/psurf)
94  zold=2.0
95  zref=zold
96  rhold=ew(td2)/ew(tt2)
97
98  thetaref=tvold*(100000./pold)**(r_air/cpa)+excess
99  thetaold=thetaref
100
101
102  ! Integrate z up to one level above zt
103  !*************************************
104  if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then
105    loop_start=2
106  else
107    loop_start=llev
108  end if
109  do k=loop_start,nuvz
110    pint=akz(k)+bkz(k)*psurf  ! pressure on model layers
111    tv=ttlev(k)*(1.+0.608*qvlev(k))
112
113    if (abs(tv-tvold).gt.0.2) then
114      z=zold+const*log(pold/pint)*(tv-tvold)/log(tv/tvold)
115    else
116      z=zold+const*log(pold/pint)*tv
117    endif
118
119    theta=tv*(100000./pint)**(r_air/cpa)
120  ! Petra
121    rh = qvlev(k) / f_qvsat( pint, ttlev(k) )
122
123
124  ! Calculate Richardson number at each level
125  !****************************************
126
127    ri=ga/thetaref*(theta-thetaref)*(z-zref)/ &
128         max(((ulev(k)-ulev(2))**2+(vlev(k)-vlev(2))**2+b*ust**2),0.1)
129
130  !  addition of second condition: MH should not be placed in an
131  !  unstable layer (PS / Feb 2000)
132    if (ri.gt.ric .and. thetaold.lt.theta) goto 20
133
134    tvold=tv
135    pold=pint
136    rhold=rh
137    thetaold=theta
138    zold=z
139  end do
140  k=k-1 ! ESO: make sure k <= nuvz (ticket #139)
14120   continue
142
143  ! Determine Richardson number between the critical levels
144  !********************************************************
145
146  zl1=zold
147  theta1=thetaold
148  do i=1,20
149    zl=zold+real(i)/20.*(z-zold)
150    ul=ulev(k-1)+real(i)/20.*(ulev(k)-ulev(k-1))
151    vl=vlev(k-1)+real(i)/20.*(vlev(k)-vlev(k-1))
152    thetal=thetaold+real(i)/20.*(theta-thetaold)
153    rhl=rhold+real(i)/20.*(rh-rhold)
154    ril=ga/thetaref*(thetal-thetaref)*(zl-zref)/ &
155         max(((ul-ulev(2))**2+(vl-vlev(2))**2+b*ust**2),0.1)
156    zl2=zl
157    theta2=thetal
158    if (ril.gt.ric) goto 25
159    zl1=zl
160    theta1=thetal
161  end do
162
16325   continue
164  h=zl
165  thetam=0.5*(theta1+theta2)
166  wspeed=sqrt(ul**2+vl**2)                    ! Wind speed at z=hmix
167  bvfsq=(ga/thetam)*(theta2-theta1)/(zl2-zl1) ! Brunt-Vaisala frequency
168                                              ! at z=hmix
169
170  ! Under stable conditions, limit the maximum effect of the subgrid-scale topography
171  ! by the maximum lifting possible from the available kinetic energy
172  !*****************************************************************************
173
174  if(bvfsq.le.0.) then
175    hmixplus=9999.
176  else
177    bvf=sqrt(bvfsq)
178    hmixplus=wspeed/bvf*convke
179  endif
180
181
182  ! Calculate convective velocity scale
183  !************************************
184
185  if (hf.lt.0.) then
186    wst=(-h*ga/thetaref*hf/cpa)**0.333
187    excess=-bs*hf/cpa/wst
188    if (iter.lt.itmax) goto 30
189  else
190    wst=0.
191  endif
192
193end subroutine richardson
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG