[16] | 1 | !*********************************************************************** |
---|
| 2 | !* Copyright 2012,2013 * |
---|
| 3 | !* Jerome Brioude, Delia Arnold, Andreas Stohl, Wayne Angevine, * |
---|
| 4 | !* John Burkhart, Massimo Cassiani, Adam Dingwell, Richard C Easter, Sabine Eckhardt,* |
---|
| 5 | !* Stephanie Evan, Jerome D Fast, Don Morton, Ignacio Pisso, * |
---|
| 6 | !* Petra Seibert, Gerard Wotawa, Caroline Forster, Harald Sodemann, * |
---|
| 7 | !* * |
---|
| 8 | !* This file is part of FLEXPART WRF * |
---|
| 9 | !* * |
---|
| 10 | !* FLEXPART is free software: you can redistribute it and/or modify * |
---|
| 11 | !* it under the terms of the GNU General Public License as published by* |
---|
| 12 | !* the Free Software Foundation, either version 3 of the License, or * |
---|
| 13 | !* (at your option) any later version. * |
---|
| 14 | !* * |
---|
| 15 | !* FLEXPART is distributed in the hope that it will be useful, * |
---|
| 16 | !* but WITHOUT ANY WARRANTY; without even the implied warranty of * |
---|
| 17 | !* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * |
---|
| 18 | !* GNU General Public License for more details. * |
---|
| 19 | !* * |
---|
| 20 | !* You should have received a copy of the GNU General Public License * |
---|
| 21 | !* along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. * |
---|
| 22 | !*********************************************************************** |
---|
| 23 | subroutine richardson(psurf,ust,ttlev,qvlev,ulev,vlev,nuvz, & |
---|
| 24 | pplev,hf,tt2,td2,h,wst,hmixplus,ierr,sfc_option) |
---|
| 25 | ! subroutine richardson(psurf,ust,ttlev,qvlev,ulev,vlev,nuvz, |
---|
| 26 | ! +akz,bkz,hf,tt2,td2,h,wst,hmixplus,ierr) |
---|
| 27 | |
---|
| 28 | ! i i i i i i i |
---|
| 29 | ! i i i i i o o o |
---|
| 30 | !**************************************************************************** |
---|
| 31 | ! * |
---|
| 32 | ! Note: This is the FLEXPART_WRF version of subroutine richardson. * |
---|
| 33 | ! * |
---|
| 34 | ! Calculation of mixing height based on the critical Richardson number. * |
---|
| 35 | ! Calculation of convective time scale. * |
---|
| 36 | ! For unstable conditions, one iteration is performed. An excess * |
---|
| 37 | ! temperature (dependent on hf and wst) is calculated, added to the * |
---|
| 38 | ! temperature at the lowest model level. Then the procedure is repeated.* |
---|
| 39 | ! * |
---|
| 40 | ! Author: A. Stohl * |
---|
| 41 | ! * |
---|
| 42 | ! 22 August 1996 * |
---|
| 43 | ! * |
---|
| 44 | ! Literature: * |
---|
| 45 | ! Vogelezang DHP and Holtslag AAM (1996): Evaluation and model impacts * |
---|
| 46 | ! of alternative boundary-layer height formulations. Boundary-Layer * |
---|
| 47 | ! Meteor. 81, 245-269. * |
---|
| 48 | ! * |
---|
| 49 | ! Update: 1999-02-01 by G. Wotawa * |
---|
| 50 | ! * |
---|
| 51 | ! Two meter level (temperature, humidity) is taken as reference level * |
---|
| 52 | ! instead of first model level. * |
---|
| 53 | ! New input variables tt2, td2 introduced. * |
---|
| 54 | ! * |
---|
| 55 | ! 17 Oct 2005 - R. Easter - added ierr status flag (0=ok, -1=fail) * |
---|
| 56 | ! 15 Nov 2005 - R. Easter - use pplev instead of akz,bkz * |
---|
| 57 | ! * |
---|
| 58 | !**************************************************************************** |
---|
| 59 | ! * |
---|
| 60 | ! Variables: * |
---|
| 61 | ! h mixing height [m] * |
---|
| 62 | ! hf sensible heat flux * |
---|
| 63 | ! psurf surface pressure at point (xt,yt) [Pa] * |
---|
| 64 | ! tv virtual temperature * |
---|
| 65 | ! wst convective velocity scale * |
---|
| 66 | ! * |
---|
| 67 | ! Constants: * |
---|
| 68 | ! ric critical Richardson number * |
---|
| 69 | ! * |
---|
| 70 | !**************************************************************************** |
---|
| 71 | |
---|
| 72 | ! include 'includepar' |
---|
| 73 | |
---|
| 74 | ! integer ierr |
---|
| 75 | ! integer i,k,nuvz,itmax,iter |
---|
| 76 | ! real tv,tvold,zref,z,zold,pint,pold,theta,thetaref,ri |
---|
| 77 | !! real akz(nuvz),bkz(nuvz),ulev(nuvz),vlev(nuvz),hf,wst,tt2,td2,ew |
---|
| 78 | ! real pplev(nuvz),ulev(nuvz),vlev(nuvz),hf,wst,tt2,td2,ew |
---|
| 79 | ! real psurf,ust,ttlev(nuvz),qvlev(nuvz),h,const,ric,b,excess,bs |
---|
| 80 | ! real thetaold,zl,ul,vl,thetal,ril,hmixplus,wspeed,bvfsq,bvf |
---|
| 81 | ! real f_qvsat,rh,rhold,rhl,theta1,theta2,zl1,zl2,thetam |
---|
| 82 | ! parameter(const=r_air/ga,ric=0.25,b=100.,bs=8.5,itmax=3) |
---|
| 83 | |
---|
| 84 | use par_mod |
---|
| 85 | |
---|
| 86 | implicit none |
---|
| 87 | |
---|
| 88 | integer :: i,k,nuvz,iter,ierr |
---|
| 89 | real :: tv,tvold,zref,z,zold,pint,pold,theta,thetaref,ri |
---|
| 90 | real :: pplev(nuvz),ulev(nuvz),vlev(nuvz),hf,wst,tt2,td2,ew |
---|
| 91 | real :: psurf,ust,ttlev(nuvz),qvlev(nuvz),h,excess |
---|
| 92 | real :: thetaold,zl,ul,vl,thetal,ril,hmixplus,wspeed,bvfsq,bvf |
---|
| 93 | real :: f_qvsat,rh,rhold,rhl,theta1,theta2,zl1,zl2,thetam |
---|
| 94 | |
---|
| 95 | real,parameter :: const=r_air/ga, ric=0.25, b=100., bs=8.5 |
---|
| 96 | integer,parameter :: itmax=3 |
---|
| 97 | |
---|
| 98 | real :: duma |
---|
| 99 | integer :: sfc_option |
---|
| 100 | |
---|
| 101 | excess=0.0 |
---|
| 102 | iter=0 |
---|
| 103 | |
---|
| 104 | ! Compute virtual temperature and virtual potential temperature at |
---|
| 105 | ! reference level (2 m) |
---|
| 106 | !***************************************************************** |
---|
| 107 | |
---|
| 108 | 30 iter=iter+1 |
---|
| 109 | |
---|
| 110 | pold=psurf |
---|
| 111 | tvold=tt2*(1.+0.378*ew(td2)/psurf) |
---|
| 112 | zold=2.0 |
---|
| 113 | zref=zold |
---|
| 114 | rhold=ew(td2)/ew(tt2) |
---|
| 115 | |
---|
| 116 | thetaref=tvold*(100000./pold)**(r_air/cpa)+excess |
---|
| 117 | thetaold=thetaref |
---|
| 118 | |
---|
| 119 | |
---|
| 120 | ! Integrate z up to one level above zt |
---|
| 121 | !************************************* |
---|
| 122 | |
---|
| 123 | do k=2,nuvz |
---|
| 124 | ! pint=akz(k)+bkz(k)*psurf ! pressure on model layers |
---|
| 125 | pint=pplev(k) ! pressure on model layers |
---|
| 126 | tv=ttlev(k)*(1.+0.608*qvlev(k)) |
---|
| 127 | |
---|
| 128 | if (abs(tv-tvold).gt.0.2) then |
---|
| 129 | z=zold+const*log(pold/pint)*(tv-tvold)/log(tv/tvold) |
---|
| 130 | else |
---|
| 131 | z=zold+const*log(pold/pint)*tv |
---|
| 132 | endif |
---|
| 133 | |
---|
| 134 | theta=tv*(100000./pint)**(r_air/cpa) |
---|
| 135 | ! Petra |
---|
| 136 | rh = qvlev(k) / f_qvsat( pint, ttlev(k) ) |
---|
| 137 | |
---|
| 138 | |
---|
| 139 | !alculate Richardson number at each level |
---|
| 140 | !**************************************** |
---|
| 141 | |
---|
| 142 | ri=ga/thetaref*(theta-thetaref)*(z-zref)/ & |
---|
| 143 | max(((ulev(k)-ulev(2))**2+(vlev(k)-vlev(2))**2+b*ust**2),0.1) |
---|
| 144 | |
---|
| 145 | ! addition of second condition: MH should not be placed in an |
---|
| 146 | ! unstable layer (PS / Feb 2000) |
---|
| 147 | if (ri.gt.ric .and. thetaold.lt.theta) goto 20 |
---|
| 148 | |
---|
| 149 | tvold=tv |
---|
| 150 | pold=pint |
---|
| 151 | rhold=rh |
---|
| 152 | thetaold=theta |
---|
| 153 | zold=z |
---|
| 154 | end do |
---|
| 155 | |
---|
| 156 | if (k .ge. nuvz) then |
---|
| 157 | write(*,*) 'richardson not working -- k = nuvz' |
---|
| 158 | ierr = -10 |
---|
| 159 | goto 7000 |
---|
| 160 | end if |
---|
| 161 | |
---|
| 162 | 20 continue |
---|
| 163 | |
---|
| 164 | ! Determine Richardson number between the critical levels |
---|
| 165 | !******************************************************** |
---|
| 166 | |
---|
| 167 | zl1=zold |
---|
| 168 | theta1=thetaold |
---|
| 169 | do i=1,20 |
---|
| 170 | zl=zold+float(i)/20.*(z-zold) |
---|
| 171 | ul=ulev(k-1)+float(i)/20.*(ulev(k)-ulev(k-1)) |
---|
| 172 | vl=vlev(k-1)+float(i)/20.*(vlev(k)-vlev(k-1)) |
---|
| 173 | thetal=thetaold+float(i)/20.*(theta-thetaold) |
---|
| 174 | rhl=rhold+float(i)/20.*(rh-rhold) |
---|
| 175 | ril=ga/thetaref*(thetal-thetaref)*(zl-zref)/ & |
---|
| 176 | max(((ul-ulev(2))**2+(vl-vlev(2))**2+b*ust**2),0.1) |
---|
| 177 | zl2=zl |
---|
| 178 | theta2=thetal |
---|
| 179 | if (ril.gt.ric) goto 25 |
---|
| 180 | zl1=zl |
---|
| 181 | theta1=thetal |
---|
| 182 | enddo |
---|
| 183 | |
---|
| 184 | 25 continue |
---|
| 185 | ! if sfc_option = sfc_option_wrf, |
---|
| 186 | ! pbl heights are read from WRF met. files and put into hmix (=h) |
---|
| 187 | !JB |
---|
| 188 | ! h=zl |
---|
| 189 | if(sfc_option .eq. sfc_option_diagnosed) h=zl |
---|
| 190 | if (h .le. 0.0) then |
---|
| 191 | write(*,*) 'richardson not working -- bad h =', h |
---|
| 192 | ierr = -20 |
---|
| 193 | goto 7000 |
---|
| 194 | ! else if (h .lt. 10.0) then |
---|
| 195 | ! write(*,*) 'richardson not working -- too small h =', h |
---|
| 196 | ! ierr = +20 |
---|
| 197 | ! return |
---|
| 198 | end if |
---|
| 199 | |
---|
| 200 | thetam=0.5*(theta1+theta2) |
---|
| 201 | wspeed=sqrt(ul**2+vl**2) ! Wind speed at z=hmix |
---|
| 202 | bvfsq=(ga/thetam)*(theta2-theta1)/(zl2-zl1) ! Brunt-Vaisala frequency |
---|
| 203 | ! at z=hmix |
---|
| 204 | |
---|
| 205 | ! Under stable conditions, limit the maximum effect of the subgrid-scale topography |
---|
| 206 | ! by the maximum lifting possible from the available kinetic energy |
---|
| 207 | !********************************************************************************** |
---|
| 208 | |
---|
| 209 | if(bvfsq.le.0.) then |
---|
| 210 | hmixplus=9999. |
---|
| 211 | else |
---|
| 212 | bvf=sqrt(bvfsq) |
---|
| 213 | hmixplus=wspeed/bvf*convke ! keconv = kinetic energy |
---|
| 214 | endif ! used for lifting |
---|
| 215 | |
---|
| 216 | |
---|
| 217 | ! Calculate convective velocity scale |
---|
| 218 | !************************************ |
---|
| 219 | |
---|
| 220 | if (hf.lt.0.) then |
---|
| 221 | wst=(-h*ga/thetaref*hf/cpa)**0.333 |
---|
| 222 | excess=-bs*hf/cpa/wst |
---|
| 223 | if (iter.lt.itmax) goto 30 |
---|
| 224 | else |
---|
| 225 | wst=0. |
---|
| 226 | endif |
---|
| 227 | |
---|
| 228 | ierr = 0 |
---|
| 229 | return |
---|
| 230 | |
---|
| 231 | ! Fatal error -- print the inputs |
---|
| 232 | 7000 continue |
---|
| 233 | write(*,'(a )') 'nuvz' |
---|
| 234 | write(*,'(i5 )') nuvz |
---|
| 235 | write(*,'(a )') 'psurf,ust,hf,tt2,td2,h,wst,hmixplus' |
---|
| 236 | write(*,'(1p,4e18.10)') psurf,ust,hf,tt2,td2,h,wst,hmixplus |
---|
| 237 | write(*,'(a )') 'ttlev' |
---|
| 238 | write(*,'(1p,4e18.10)') ttlev |
---|
| 239 | write(*,'(a )') 'qvlev' |
---|
| 240 | write(*,'(1p,4e18.10)') qvlev |
---|
| 241 | write(*,'(a )') 'ulev' |
---|
| 242 | write(*,'(1p,4e18.10)') ulev |
---|
| 243 | write(*,'(a )') 'vlev' |
---|
| 244 | write(*,'(1p,4e18.10)') vlev |
---|
| 245 | write(*,'(a )') 'pplev' |
---|
| 246 | write(*,'(1p,4e18.10)') pplev |
---|
| 247 | return |
---|
| 248 | |
---|
| 249 | end subroutine richardson |
---|
| 250 | |
---|