source: flexpart.git/src/calcpar.f90 @ 027e844

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 027e844 was 027e844, checked in by Espen Sollum ATMOS <eso@…>, 6 years ago

Bugfix for calcpar.

  • Property mode set to 100644
File size: 11.3 KB
RevLine 
[e200b7a]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
[6ecb30a]22subroutine calcpar(n,uuh,vvh,pvh,metdata_format)
[e200b7a]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  !                                                                            *
[6ecb30a]39  !   CHANGE 17/11/2005 Caroline Forster NCEP GFS version                      *
40  !                                                                            *
41  !   Changes, Bernd C. Krueger, Feb. 2001:                                    *
42  !    Variables tth and qvh (on eta coordinates) in common block              *
43  !                                                                            *
44  !   Unified ECMWF and GFS builds                                             *
45  !   Marian Harustak, 12.5.2017                                               *
46  !     - Merged calcpar and calcpar_gfs into one routine using if-then        *
47  !       for meteo-type dependent code                                        *
[e200b7a]48  !*****************************************************************************
[6ecb30a]49
[e200b7a]50  !*****************************************************************************
51  !                                                                            *
52  ! Variables:                                                                 *
53  ! n                  temporal index for meteorological fields (1 to 3)       *
[6ecb30a]54  ! uuh                                                                        *
55  ! vvh                                                                        *
56  ! pvh                                                                        *
57  ! metdata_format     format of metdata (ecmwf/gfs)                           *
[e200b7a]58  !                                                                            *
59  ! Constants:                                                                 *
60  !                                                                            *
61  !                                                                            *
62  ! Functions:                                                                 *
63  ! scalev             computation of ustar                                    *
64  ! obukhov            computatio of Obukhov length                            *
65  !                                                                            *
66  !*****************************************************************************
67
68  use par_mod
69  use com_mod
[6ecb30a]70  use class_gribfile
[e200b7a]71
72  implicit none
73
[6ecb30a]74  integer :: metdata_format
75  integer :: n,ix,jy,i,kz,lz,kzmin,llev,loop_start
[e200b7a]76  real :: ttlev(nuvzmax),qvlev(nuvzmax),obukhov,scalev,ol,hmixplus
77  real :: ulev(nuvzmax),vlev(nuvzmax),ew,rh,vd(maxspec),subsceff,ylat
[027e844]78  real :: altmin,tvold,pold,zold,pint,tv,zlev(nuvzmax),hmixdummy,akzdummy
[e200b7a]79  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
80  real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
81  real :: pvh(0:nxmax-1,0:nymax-1,nuvzmax)
82  real,parameter :: const=r_air/ga
83
84  !write(*,*) 'in calcpar writting snowheight'
85  !***********************************
86  ! for test: write out snow depths
87
88  ! open(4,file='slandusetest',form='formatted')
89  ! do 5 ix=0,nxmin1
90  !5       write (4,*) (sd(ix,jy,1,n),jy=0,nymin1)
91  !  close(4)
92
93
94  ! Loop over entire grid
95  !**********************
96
97  do jy=0,nymin1
98
99  ! Set minimum height for tropopause
100  !**********************************
101
102    ylat=ylat0+real(jy)*dy
103    if ((ylat.ge.-20.).and.(ylat.le.20.)) then
104      altmin = 5000.
105    else
106      if ((ylat.gt.20.).and.(ylat.lt.40.)) then
107        altmin=2500.+(40.-ylat)*125.
108      else if ((ylat.gt.-40.).and.(ylat.lt.-20.)) then
109        altmin=2500.+(40.+ylat)*125.
110      else
111        altmin=2500.
112      endif
113    endif
114
115    do ix=0,nxmin1
116
117  ! 1) Calculation of friction velocity
118  !************************************
119
120      ustar(ix,jy,1,n)=scalev(ps(ix,jy,1,n),tt2(ix,jy,1,n), &
121           td2(ix,jy,1,n),surfstr(ix,jy,1,n))
122      if (ustar(ix,jy,1,n).le.1.e-8) ustar(ix,jy,1,n)=1.e-8
123
124  ! 2) Calculation of inverse Obukhov length scale
125  !***********************************************
126
[6ecb30a]127      if (metdata_format.eq.GRIBFILE_CENTRE_NCEP) then
128        ! NCEP version: find first level above ground
129        llev = 0
130        do i=1,nuvz
131          if (ps(ix,jy,1,n).lt.akz(i)) llev=i
132        end do
133        llev = llev+1
134        if (llev.gt.nuvz) llev = nuvz-1
135        ! NCEP version
136
137        ! calculate inverse Obukhov length scale with tth(llev)
[027e844]138        ol=obukhov(ps(ix,jy,1,n),tt2(ix,jy,1,n),td2(ix,jy,1,n), &
139             tth(ix,jy,llev,n),ustar(ix,jy,1,n),sshf(ix,jy,1,n),akm,bkm,akz(llev),metdata_format)
[6ecb30a]140      else
141        llev=0
142        ol=obukhov(ps(ix,jy,1,n),tt2(ix,jy,1,n),td2(ix,jy,1,n), &
[027e844]143            tth(ix,jy,2,n),ustar(ix,jy,1,n),sshf(ix,jy,1,n),akm,bkm,akzdummy,metdata_format)
[6ecb30a]144      end if
145
[e200b7a]146      if (ol.ne.0.) then
147        oli(ix,jy,1,n)=1./ol
148      else
149        oli(ix,jy,1,n)=99999.
150      endif
151
152
153  ! 3) Calculation of convective velocity scale and mixing height
154  !**************************************************************
155
156      do i=1,nuvz
157        ulev(i)=uuh(ix,jy,i)
158        vlev(i)=vvh(ix,jy,i)
159        ttlev(i)=tth(ix,jy,i,n)
160        qvlev(i)=qvh(ix,jy,i,n)
161      end do
162
[6ecb30a]163      if (metdata_format.eq.GRIBFILE_CENTRE_NCEP) then
164        ! NCEP version hmix has been read in in readwind.f, is therefore not calculated here
[e200b7a]165      call richardson(ps(ix,jy,1,n),ustar(ix,jy,1,n),ttlev,qvlev, &
166           ulev,vlev,nuvz,akz,bkz,sshf(ix,jy,1,n),tt2(ix,jy,1,n), &
[6ecb30a]167             td2(ix,jy,1,n),hmixdummy,wstar(ix,jy,1,n),hmixplus,metdata_format)
168      else
169        call richardson(ps(ix,jy,1,n),ustar(ix,jy,1,n),ttlev,qvlev, &
170             ulev,vlev,nuvz,akz,bkz,sshf(ix,jy,1,n),tt2(ix,jy,1,n), &
171             td2(ix,jy,1,n),hmix(ix,jy,1,n),wstar(ix,jy,1,n),hmixplus,metdata_format)
172      end if
[e200b7a]173
174      if(lsubgrid.eq.1) then
175        subsceff=min(excessoro(ix,jy),hmixplus)
176      else
177        subsceff=0.0
178      endif
179  !
180  ! CALCULATE HMIX EXCESS ACCORDING TO SUBGRIDSCALE VARIABILITY AND STABILITY
181  !
182      hmix(ix,jy,1,n)=hmix(ix,jy,1,n)+subsceff
183      hmix(ix,jy,1,n)=max(hmixmin,hmix(ix,jy,1,n)) ! set minimum PBL height
184      hmix(ix,jy,1,n)=min(hmixmax,hmix(ix,jy,1,n)) ! set maximum PBL height
185
186  ! 4) Calculation of dry deposition velocities
187  !********************************************
188
189      if (DRYDEP) then
190  ! Sabine Eckhardt, Dec 06: use new index for z0 for water depending on
191  ! windspeed
192        z0(7)=0.016*ustar(ix,jy,1,n)*ustar(ix,jy,1,n)/ga
193
194  ! Calculate relative humidity at surface
195  !***************************************
196        rh=ew(td2(ix,jy,1,n))/ew(tt2(ix,jy,1,n))
197
198        call getvdep(n,ix,jy,ustar(ix,jy,1,n), &
199             tt2(ix,jy,1,n),ps(ix,jy,1,n),1./oli(ix,jy,1,n), &
200             ssr(ix,jy,1,n),rh,lsprec(ix,jy,1,n)+convprec(ix,jy,1,n), &
201             sd(ix,jy,1,n),vd)
202
203        do i=1,nspec
204          vdep(ix,jy,i,n)=vd(i)
205        end do
206
207      endif
208
209  !******************************************************
210  ! Calculate height of thermal tropopause (Hoinka, 1997)
211  !******************************************************
212
[6ecb30a]213  ! 1) Calculate altitudes of model levels
214  !***************************************
[e200b7a]215
216      tvold=tt2(ix,jy,1,n)*(1.+0.378*ew(td2(ix,jy,1,n))/ &
217           ps(ix,jy,1,n))
218      pold=ps(ix,jy,1,n)
219      zold=0.
[6ecb30a]220      if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then
221        loop_start=2
222      else
223        loop_start=llev
224      end if
225      do kz=loop_start,nuvz
[e200b7a]226        pint=akz(kz)+bkz(kz)*ps(ix,jy,1,n)  ! pressure on model layers
227        tv=tth(ix,jy,kz,n)*(1.+0.608*qvh(ix,jy,kz,n))
228
229        if (abs(tv-tvold).gt.0.2) then
230         zlev(kz)=zold+const*log(pold/pint)*(tv-tvold)/log(tv/tvold)
231        else
232          zlev(kz)=zold+const*log(pold/pint)*tv
233        endif
234        tvold=tv
235        pold=pint
236        zold=zlev(kz)
237      end do
238
239  ! 2) Define a minimum level kzmin, from which upward the tropopause is
240  !    searched for. This is to avoid inversions in the lower troposphere
241  !    to be identified as the tropopause
242  !************************************************************************
243
[6ecb30a]244      if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then
245        loop_start=1
246      else
247        loop_start=llev
248      end if
249
250      do kz=loop_start,nuvz
[e200b7a]251        if (zlev(kz).ge.altmin) then
252          kzmin=kz
253          goto 45
254        endif
255      end do
25645    continue
257
258  ! 3) Search for first stable layer above minimum height that fulfills the
259  !    thermal tropopause criterion
260  !************************************************************************
261
262      do kz=kzmin,nuvz
263        do lz=kz+1,nuvz
264          if ((zlev(lz)-zlev(kz)).gt.2000.) then
265            if (((tth(ix,jy,kz,n)-tth(ix,jy,lz,n))/ &
266                 (zlev(lz)-zlev(kz))).lt.0.002) then
267              tropopause(ix,jy,1,n)=zlev(kz)
268              goto 51
269            endif
270            goto 50
271          endif
272        end do
27350      continue
274      end do
27551    continue
276
277
278    end do
279  end do
280
281  ! Calculation of potential vorticity on 3-d grid
282  !***********************************************
283
284  call calcpv(n,uuh,vvh,pvh)
285
286
287end subroutine calcpar
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG