source: flexpart.git/src/pbl_profile.f90

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file was 92fab65, checked in by Ignacio Pisso <ip@…>, 4 years ago

add SPDX-License-Identifier to all .f90 files

  • Property mode set to 100644
File size: 5.0 KB
Line 
1! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
2! SPDX-License-Identifier: GPL-3.0-or-later
3
4subroutine pbl_profile(ps,td2m,zml1,t2m,tml1,u10m,uml1,stress,hf)
5
6  !********************************************************************
7  !                                                                   *
8  !                    G. WOTAWA, 1995-07-07                          *
9  !                                                                   *
10  !********************************************************************
11  !                                                                   *
12  ! DESCRIPTION: CALCULATION OF FRICTION VELOCITY AND SURFACE SENS-   *
13  !              IBLE HEAT FLUX USING THE PROFILE METHOD (BERKOVICZ   *
14  !              AND PRAHM, 1982)                                     *
15  !                                                                   *
16  ! Output now is surface stress instead of ustar                     *
17  !                                                                   *
18  !                                                                   *
19  !********************************************************************
20  !                                                                   *
21  ! INPUT:                                                            *
22  !                                                                   *
23  !                                                                   *
24  ! ps      surface pressure(Pa)                                      *
25  ! td2m    two metre dew point(K)                                    *
26  ! zml1    heigth of first model level (m)                           *
27  ! t2m     two metre temperature (K)                                 *
28  ! tml1    temperature first model level (K)                         *
29  ! u10m    ten metre wind speed (ms-1)                               *
30  ! uml1    wind speed first model level (ms-1)                       *
31  !                                                                   *
32  !********************************************************************
33  !                                                                   *
34  ! OUTPUT:                                                           *
35  !                                                                   *
36  ! stress  surface stress (i.e., friction velocity (ms-1) squared    *
37  !                         multiplied with air density)              *
38  ! hf      surface sensible heat flux (Wm-2)                         *
39  !                                                                   *
40  !********************************************************************
41  ! ustar   friction velocity (ms-1)                                  *
42  ! maxiter maximum number of iterations                              *
43  !********************************************************************
44
45  use par_mod
46
47  implicit none
48
49  integer :: iter
50  real :: ps,td2m,rhoa,zml1,t2m,tml1,u10m,uml1,ustar,hf
51  real :: al,alold,aldiff,tmean,crit
52  real :: deltau,deltat,thetastar,psim,psih,e,ew,tv,stress
53  integer,parameter :: maxiter=10
54  real,parameter    :: r1=0.74
55
56  e=ew(td2m)               ! vapor pressure
57  tv=t2m*(1.+0.378*e/ps)   ! virtual temperature
58  rhoa=ps/(r_air*tv)       ! air density
59
60  deltau=uml1-u10m         !! Wind Speed difference between
61                           !! Model level 1 and 10 m
62
63  if(deltau.le.0.001) then    !! Monin-Obukhov Theory not
64    al=9999.               !! applicable --> Set dummy values
65    ustar=0.01
66    stress=ustar*ustar*rhoa
67    hf=0.0
68    return
69  endif
70  deltat=tml1-t2m+0.0098*(zml1-2.)  !! Potential temperature difference
71                                    !! between model level 1 and 10 m
72
73  if(abs(deltat).le.0.03) then    !! Neutral conditions
74    hf=0.0
75    al=9999.
76    ustar=(vonkarman*deltau)/ &
77         (log(zml1/10.)-psim(zml1,al)+psim(10.,al))
78    stress=ustar*ustar*rhoa
79    return
80  endif
81
82  tmean=0.5*(t2m+tml1)
83  crit=(0.0219*tmean*(zml1-2.0)*deltau**2)/ &
84       (deltat*(zml1-10.0)**2)
85  if((deltat.gt.0).and.(crit.le.1.)) then
86                                    !! Successive approximation will
87    al=50.                          !! not converge
88    ustar=(vonkarman*deltau)/ &
89         (log(zml1/10.)-psim(zml1,al)+psim(10.,al))
90    thetastar=(vonkarman*deltat/r1)/ &
91         (log(zml1/2.)-psih(zml1,al)+psih(2.,al))
92    hf=rhoa*cpa*ustar*thetastar
93    stress=ustar*ustar*rhoa
94    return
95  endif
96
97  al=9999.                 ! Start iteration assuming neutral conditions
98  do iter=1,maxiter
99    alold=al
100    ustar=(vonkarman*deltau)/ &
101         (log(zml1/10.)-psim(zml1,al)+psim(10.,al))
102    thetastar=(vonkarman*deltat/r1)/ &
103         (log(zml1/2.)-psih(zml1,al)+psih(2.,al))
104    al=(tmean*ustar**2)/(ga*vonkarman*thetastar)
105    aldiff=abs((al-alold)/alold)
106    if(aldiff.lt.0.01) goto 30  !! Successive approximation successful
107  end do
10830   hf=rhoa*cpa*ustar*thetastar
109  if(al.gt.9999.) al=9999.
110  if(al.lt.-9999.) al=-9999.
111
112  stress=ustar*ustar*rhoa
113
114end subroutine pbl_profile
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG