source: flexpart.git/src/pbl_profile.f90 @ 3481cc1

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

move license from headers to a different file

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