source: trunk/src/calcpar.f90 @ 28

Last change on this file since 28 was 4, checked in by mlanger, 11 years ago
File size: 9.1 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 calcpar(n,uuh,vvh,pvh)
23  !                   i  i   i   o
24  !*****************************************************************************
25  !                                                                            *
26  !     Computation of several boundary layer parameters needed for the        *
27  !     dispersion calculation and calculation of dry deposition velocities.   *
28  !     All parameters are calculated over the entire grid.                    *
29  !                                                                            *
30  !     Author: A. Stohl                                                       *
31  !                                                                            *
32  !     21 May 1995                                                            *
33  !                                                                            *
34  ! ------------------------------------------------------------------         *
35  !     Petra Seibert, Feb 2000:                                               *
36  !     convection scheme:                                                     *
37  !     new variables in call to richardson                                    *
38  !                                                                            *
39  !*****************************************************************************
40  !  Changes, Bernd C. Krueger, Feb. 2001:
41  !   Variables tth and qvh (on eta coordinates) in common block
42  !*****************************************************************************
43  !                                                                            *
44  ! Variables:                                                                 *
45  ! n                  temporal index for meteorological fields (1 to 3)       *
46  !                                                                            *
47  ! Constants:                                                                 *
48  !                                                                            *
49  !                                                                            *
50  ! Functions:                                                                 *
51  ! scalev             computation of ustar                                    *
52  ! obukhov            computatio of Obukhov length                            *
53  !                                                                            *
54  !*****************************************************************************
55
56  use par_mod
57  use com_mod
58
59  implicit none
60
61  integer :: n,ix,jy,i,kz,lz,kzmin
62  real :: ttlev(nuvzmax),qvlev(nuvzmax),obukhov,scalev,ol,hmixplus
63  real :: ulev(nuvzmax),vlev(nuvzmax),ew,rh,vd(maxspec),subsceff,ylat
64  real :: altmin,tvold,pold,zold,pint,tv,zlev(nuvzmax)
65  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
66  real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
67  real :: pvh(0:nxmax-1,0:nymax-1,nuvzmax)
68  real,parameter :: const=r_air/ga
69
70  !write(*,*) 'in calcpar writting snowheight'
71  !***********************************
72  ! for test: write out snow depths
73
74  ! open(4,file='slandusetest',form='formatted')
75  ! do 5 ix=0,nxmin1
76  !5       write (4,*) (sd(ix,jy,1,n),jy=0,nymin1)
77  !  close(4)
78
79
80  ! Loop over entire grid
81  !**********************
82
83  do jy=0,nymin1
84
85  ! Set minimum height for tropopause
86  !**********************************
87
88    ylat=ylat0+real(jy)*dy
89    if ((ylat.ge.-20.).and.(ylat.le.20.)) then
90      altmin = 5000.
91    else
92      if ((ylat.gt.20.).and.(ylat.lt.40.)) then
93        altmin=2500.+(40.-ylat)*125.
94      else if ((ylat.gt.-40.).and.(ylat.lt.-20.)) then
95        altmin=2500.+(40.+ylat)*125.
96      else
97        altmin=2500.
98      endif
99    endif
100
101    do ix=0,nxmin1
102
103  ! 1) Calculation of friction velocity
104  !************************************
105
106      ustar(ix,jy,1,n)=scalev(ps(ix,jy,1,n),tt2(ix,jy,1,n), &
107           td2(ix,jy,1,n),surfstr(ix,jy,1,n))
108      if (ustar(ix,jy,1,n).le.1.e-8) ustar(ix,jy,1,n)=1.e-8
109
110  ! 2) Calculation of inverse Obukhov length scale
111  !***********************************************
112
113      ol=obukhov(ps(ix,jy,1,n),tt2(ix,jy,1,n),td2(ix,jy,1,n), &
114           tth(ix,jy,2,n),ustar(ix,jy,1,n),sshf(ix,jy,1,n),akm,bkm)
115      if (ol.ne.0.) then
116        oli(ix,jy,1,n)=1./ol
117      else
118        oli(ix,jy,1,n)=99999.
119      endif
120
121
122  ! 3) Calculation of convective velocity scale and mixing height
123  !**************************************************************
124
125      do i=1,nuvz
126        ulev(i)=uuh(ix,jy,i)
127        vlev(i)=vvh(ix,jy,i)
128        ttlev(i)=tth(ix,jy,i,n)
129        qvlev(i)=qvh(ix,jy,i,n)
130      end do
131
132      call richardson(ps(ix,jy,1,n),ustar(ix,jy,1,n),ttlev,qvlev, &
133           ulev,vlev,nuvz,akz,bkz,sshf(ix,jy,1,n),tt2(ix,jy,1,n), &
134           td2(ix,jy,1,n),hmix(ix,jy,1,n),wstar(ix,jy,1,n),hmixplus)
135
136      if(lsubgrid.eq.1) then
137        subsceff=min(excessoro(ix,jy),hmixplus)
138      else
139        subsceff=0.0
140      endif
141  !
142  ! CALCULATE HMIX EXCESS ACCORDING TO SUBGRIDSCALE VARIABILITY AND STABILITY
143  !
144      hmix(ix,jy,1,n)=hmix(ix,jy,1,n)+subsceff
145      hmix(ix,jy,1,n)=max(hmixmin,hmix(ix,jy,1,n)) ! set minimum PBL height
146      hmix(ix,jy,1,n)=min(hmixmax,hmix(ix,jy,1,n)) ! set maximum PBL height
147
148  ! 4) Calculation of dry deposition velocities
149  !********************************************
150
151      if (DRYDEP) then
152  ! Sabine Eckhardt, Dec 06: use new index for z0 for water depending on
153  ! windspeed
154        z0(7)=0.016*ustar(ix,jy,1,n)*ustar(ix,jy,1,n)/ga
155
156  ! Calculate relative humidity at surface
157  !***************************************
158        rh=ew(td2(ix,jy,1,n))/ew(tt2(ix,jy,1,n))
159
160        call getvdep(n,ix,jy,ustar(ix,jy,1,n), &
161             tt2(ix,jy,1,n),ps(ix,jy,1,n),1./oli(ix,jy,1,n), &
162             ssr(ix,jy,1,n),rh,lsprec(ix,jy,1,n)+convprec(ix,jy,1,n), &
163             sd(ix,jy,1,n),vd)
164
165        do i=1,nspec
166          vdep(ix,jy,i,n)=vd(i)
167        end do
168
169      endif
170
171  !******************************************************
172  ! Calculate height of thermal tropopause (Hoinka, 1997)
173  !******************************************************
174
175  ! 1) Calculate altitudes of ECMWF model levels
176  !*********************************************
177
178      tvold=tt2(ix,jy,1,n)*(1.+0.378*ew(td2(ix,jy,1,n))/ &
179           ps(ix,jy,1,n))
180      pold=ps(ix,jy,1,n)
181      zold=0.
182      do kz=2,nuvz
183        pint=akz(kz)+bkz(kz)*ps(ix,jy,1,n)  ! pressure on model layers
184        tv=tth(ix,jy,kz,n)*(1.+0.608*qvh(ix,jy,kz,n))
185
186        if (abs(tv-tvold).gt.0.2) then
187         zlev(kz)=zold+const*log(pold/pint)*(tv-tvold)/log(tv/tvold)
188        else
189          zlev(kz)=zold+const*log(pold/pint)*tv
190        endif
191        tvold=tv
192        pold=pint
193        zold=zlev(kz)
194      end do
195
196  ! 2) Define a minimum level kzmin, from which upward the tropopause is
197  !    searched for. This is to avoid inversions in the lower troposphere
198  !    to be identified as the tropopause
199  !************************************************************************
200
201      do kz=1,nuvz
202        if (zlev(kz).ge.altmin) then
203          kzmin=kz
204          goto 45
205        endif
206      end do
20745    continue
208
209  ! 3) Search for first stable layer above minimum height that fulfills the
210  !    thermal tropopause criterion
211  !************************************************************************
212
213      do kz=kzmin,nuvz
214        do lz=kz+1,nuvz
215          if ((zlev(lz)-zlev(kz)).gt.2000.) then
216            if (((tth(ix,jy,kz,n)-tth(ix,jy,lz,n))/ &
217                 (zlev(lz)-zlev(kz))).lt.0.002) then
218              tropopause(ix,jy,1,n)=zlev(kz)
219              goto 51
220            endif
221            goto 50
222          endif
223        end do
22450      continue
225      end do
22651    continue
227
228
229    end do
230  end do
231
232  ! Calculation of potential vorticity on 3-d grid
233  !***********************************************
234
235  call calcpv(n,uuh,vvh,pvh)
236
237
238end subroutine calcpar
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG