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 | |
---|
24 | subroutine pbl_profile(ps,td2m,zml1,t2m,tml1,u10m,uml1,stress,hf) |
---|
25 | |
---|
26 | !******************************************************************** |
---|
27 | ! * |
---|
28 | ! G. WOTAWA, 1995-07-07 * |
---|
29 | ! * |
---|
30 | !******************************************************************** |
---|
31 | ! * |
---|
32 | ! DESCRIPTION: CALCULATION OF FRICTION VELOCITY AND SURFACE SENS- * |
---|
33 | ! IBLE HEAT FLUX USING THE PROFILE METHOD (BERKOVICZ * |
---|
34 | ! AND PRAHM, 1982) * |
---|
35 | ! * |
---|
36 | ! Output now is surface stress instead of ustar * |
---|
37 | ! * |
---|
38 | ! * |
---|
39 | !******************************************************************** |
---|
40 | ! * |
---|
41 | ! INPUT: * |
---|
42 | ! * |
---|
43 | ! * |
---|
44 | ! ps surface pressure(Pa) * |
---|
45 | ! td2m two metre dew point(K) * |
---|
46 | ! zml1 heigth of first model level (m) * |
---|
47 | ! t2m two metre temperature (K) * |
---|
48 | ! tml1 temperature first model level (K) * |
---|
49 | ! u10m ten metre wind speed (ms-1) * |
---|
50 | ! uml1 wind speed first model level (ms-1) * |
---|
51 | ! * |
---|
52 | !******************************************************************** |
---|
53 | ! * |
---|
54 | ! OUTPUT: * |
---|
55 | ! * |
---|
56 | ! stress surface stress (i.e., friction velocity (ms-1) squared * |
---|
57 | ! multiplied with air density) * |
---|
58 | ! hf surface sensible heat flux (Wm-2) * |
---|
59 | ! * |
---|
60 | !******************************************************************** |
---|
61 | ! ustar friction velocity (ms-1) * |
---|
62 | ! maxiter maximum number of iterations * |
---|
63 | !******************************************************************** |
---|
64 | |
---|
65 | use par_mod |
---|
66 | |
---|
67 | implicit none |
---|
68 | |
---|
69 | integer :: iter |
---|
70 | real :: ps,td2m,rhoa,zml1,t2m,tml1,u10m,uml1,ustar,hf |
---|
71 | real :: al,alold,aldiff,tmean,crit |
---|
72 | real :: deltau,deltat,thetastar,psim,psih,e,ew,tv,stress |
---|
73 | integer,parameter :: maxiter=10 |
---|
74 | real,parameter :: r1=0.74 |
---|
75 | |
---|
76 | e=ew(td2m) ! vapor pressure |
---|
77 | tv=t2m*(1.+0.378*e/ps) ! virtual temperature |
---|
78 | rhoa=ps/(r_air*tv) ! air density |
---|
79 | |
---|
80 | deltau=uml1-u10m !! Wind Speed difference between |
---|
81 | !! Model level 1 and 10 m |
---|
82 | |
---|
83 | if(deltau.le.0.001) then !! Monin-Obukhov Theory not |
---|
84 | al=9999. !! applicable --> Set dummy values |
---|
85 | ustar=0.01 |
---|
86 | stress=ustar*ustar*rhoa |
---|
87 | hf=0.0 |
---|
88 | return |
---|
89 | endif |
---|
90 | deltat=tml1-t2m+0.0098*(zml1-2.) !! Potential temperature difference |
---|
91 | !! between model level 1 and 10 m |
---|
92 | |
---|
93 | if(abs(deltat).le.0.03) then !! Neutral conditions |
---|
94 | hf=0.0 |
---|
95 | al=9999. |
---|
96 | ustar=(vonkarman*deltau)/ & |
---|
97 | (log(zml1/10.)-psim(zml1,al)+psim(10.,al)) |
---|
98 | stress=ustar*ustar*rhoa |
---|
99 | return |
---|
100 | endif |
---|
101 | |
---|
102 | tmean=0.5*(t2m+tml1) |
---|
103 | crit=(0.0219*tmean*(zml1-2.0)*deltau**2)/ & |
---|
104 | (deltat*(zml1-10.0)**2) |
---|
105 | if((deltat.gt.0).and.(crit.le.1.)) then |
---|
106 | !! Successive approximation will |
---|
107 | al=50. !! not converge |
---|
108 | ustar=(vonkarman*deltau)/ & |
---|
109 | (log(zml1/10.)-psim(zml1,al)+psim(10.,al)) |
---|
110 | thetastar=(vonkarman*deltat/r1)/ & |
---|
111 | (log(zml1/2.)-psih(zml1,al)+psih(2.,al)) |
---|
112 | hf=rhoa*cpa*ustar*thetastar |
---|
113 | stress=ustar*ustar*rhoa |
---|
114 | return |
---|
115 | endif |
---|
116 | |
---|
117 | al=9999. ! Start iteration assuming neutral conditions |
---|
118 | do iter=1,maxiter |
---|
119 | alold=al |
---|
120 | ustar=(vonkarman*deltau)/ & |
---|
121 | (log(zml1/10.)-psim(zml1,al)+psim(10.,al)) |
---|
122 | thetastar=(vonkarman*deltat/r1)/ & |
---|
123 | (log(zml1/2.)-psih(zml1,al)+psih(2.,al)) |
---|
124 | al=(tmean*ustar**2)/(ga*vonkarman*thetastar) |
---|
125 | aldiff=abs((al-alold)/alold) |
---|
126 | if(aldiff.lt.0.01) goto 30 !! Successive approximation successful |
---|
127 | end do |
---|
128 | 30 hf=rhoa*cpa*ustar*thetastar |
---|
129 | if(al.gt.9999.) al=9999. |
---|
130 | if(al.lt.-9999.) al=-9999. |
---|
131 | |
---|
132 | stress=ustar*ustar*rhoa |
---|
133 | |
---|
134 | end subroutine pbl_profile |
---|