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