source: trunk/src/pbl_profile.f90 @ 28

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