source: flexpart.git/src/calcpar.f90 @ c0884a8

univie
Last change on this file since c0884a8 was c0884a8, checked in by pesei <petra seibert at univie ac at>, 6 years ago

replace CTBTO code for checking type of GRIB

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