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