source: flexpart.git/src/calcpar_nests.f90 @ 7999df47

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

Completed handling of nested wind fields with cloud water (for wet deposition).

  • Property mode set to 100644
File size: 9.3 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_nests(n,uuhn,vvhn,pvhn)
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  !     This routine is similar to calcpar, but is used for the nested grids.  *
30  !                                                                            *
31  !     Author: A. Stohl                                                       *
32  !                                                                            *
33  !     8 February 1999                                                        *
34  !                                                                            *
35  ! ------------------------------------------------------------------         *
36  !     Petra Seibert, Feb 2000:                                               *
37  !     convection scheme:                                                     *
38  !     new variables in call to richardson                                    *
39  !                                                                            *
40  !*****************************************************************************
41  !  Changes, Bernd C. Krueger, Feb. 2001:
42  !   Variables tth and qvh (on eta coordinates) in common block
43  !*****************************************************************************
44  !                                                                            *
45  ! Variables:                                                                 *
46  ! n                  temporal index for meteorological fields (1 to 3)       *
47  !                                                                            *
48  ! Constants:                                                                 *
49  !                                                                            *
50  !                                                                            *
51  ! Functions:                                                                 *
52  ! scalev             computation of ustar                                    *
53  ! obukhov            computatio of Obukhov length                            *
54  !                                                                            *
55  !*****************************************************************************
56
57  use par_mod
58  use com_mod
59
60  implicit none
61
62  integer :: n,ix,jy,i,l,kz,lz,kzmin
63  real :: ttlev(nuvzmax),qvlev(nuvzmax),obukhov,scalev,ol,hmixplus
64  real :: ulev(nuvzmax),vlev(nuvzmax),ew,rh,vd(maxspec),subsceff,ylat
65  real :: altmin,tvold,pold,zold,pint,tv,zlev(nuvzmax)
66  real :: uuhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
67  real :: vvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
68  real :: pvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
69  real,parameter :: const=r_air/ga
70
71
72  ! Loop over all nests
73  !********************
74
75  do l=1,numbnests
76
77  ! Loop over entire grid
78  !**********************
79
80  do jy=0,nyn(l)-1
81
82  ! Set minimum height for tropopause
83  !**********************************
84
85    ylat=ylat0n(l)+real(jy)*dyn(l)
86    if ((ylat.ge.-20.).and.(ylat.le.20.)) then
87      altmin = 5000.
88    else
89      if ((ylat.gt.20.).and.(ylat.lt.40.)) then
90        altmin=2500.+(40.-ylat)*125.
91      else if ((ylat.gt.-40.).and.(ylat.lt.-20.)) then
92        altmin=2500.+(40.+ylat)*125.
93      else
94        altmin=2500.
95      endif
96    endif
97
98    do ix=0,nxn(l)-1
99
100  ! 1) Calculation of friction velocity
101  !************************************
102
103      ustarn(ix,jy,1,n,l)=scalev(psn(ix,jy,1,n,l),tt2n(ix,jy,1,n,l), &
104           td2n(ix,jy,1,n,l),surfstrn(ix,jy,1,n,l))
105      if (ustarn(ix,jy,1,n,l).le.1.e-8) ustarn(ix,jy,1,n,l)=1.e-8
106
107  ! 2) Calculation of inverse Obukhov length scale
108  !***********************************************
109
110      ol=obukhov(psn(ix,jy,1,n,l),tt2n(ix,jy,1,n,l), &
111           td2n(ix,jy,1,n,l),tthn(ix,jy,2,n,l),ustarn(ix,jy,1,n,l), &
112           sshfn(ix,jy,1,n,l),akm,bkm)
113      if (ol.ne.0.) then
114        olin(ix,jy,1,n,l)=1./ol
115      else
116        olin(ix,jy,1,n,l)=99999.
117      endif
118
119
120  ! 3) Calculation of convective velocity scale and mixing height
121  !**************************************************************
122
123      do i=1,nuvz
124        ulev(i)=uuhn(ix,jy,i,l)
125        vlev(i)=vvhn(ix,jy,i,l)
126        ttlev(i)=tthn(ix,jy,i,n,l)
127        qvlev(i)=qvhn(ix,jy,i,n,l)
128      end do
129
130      call richardson(psn(ix,jy,1,n,l),ustarn(ix,jy,1,n,l),ttlev, &
131           qvlev,ulev,vlev,nuvz,akz,bkz,sshfn(ix,jy,1,n,l), &
132           tt2n(ix,jy,1,n,l),td2n(ix,jy,1,n,l),hmixn(ix,jy,1,n,l), &
133           wstarn(ix,jy,1,n,l),hmixplus)
134
135      if(lsubgrid.eq.1) then
136        subsceff=min(excessoron(ix,jy,l),hmixplus)
137      else
138        subsceff=0.0
139      endif
140  !
141  ! CALCULATE HMIX EXCESS ACCORDING TO SUBGRIDSCALE VARIABILITY AND STABILITY
142  !
143      hmixn(ix,jy,1,n,l)=hmixn(ix,jy,1,n,l)+subsceff
144      hmixn(ix,jy,1,n,l)=max(hmixmin,hmixn(ix,jy,1,n,l)) ! minim PBL height
145      hmixn(ix,jy,1,n,l)=min(hmixmax,hmixn(ix,jy,1,n,l)) ! maxim PBL height
146
147
148  ! 4) Calculation of dry deposition velocities
149  !********************************************
150
151      if (DRYDEP) then
152        ! z0(4)=0.016*ustarn(ix,jy,1,n,l)*ustarn(ix,jy,1,n,l)/ga
153        ! z0(9)=0.016*ustarn(ix,jy,1,n,l)*ustarn(ix,jy,1,n,l)/ga
154        z0(7)=0.016*ustarn(ix,jy,1,n,l)*ustarn(ix,jy,1,n,l)/ga
155
156  ! Calculate relative humidity at surface
157  !***************************************
158        rh=ew(td2n(ix,jy,1,n,l))/ew(tt2n(ix,jy,1,n,l))
159
160        call getvdep_nests(n,ix,jy,ustarn(ix,jy,1,n,l), &
161             tt2n(ix,jy,1,n,l),psn(ix,jy,1,n,l),1./olin(ix,jy,1,n,l), &
162             ssrn(ix,jy,1,n,l),rh,lsprecn(ix,jy,1,n,l)+ &
163             convprecn(ix,jy,1,n,l),sdn(ix,jy,1,n,l),vd,l)
164
165        do i=1,nspec
166          vdepn(ix,jy,i,n,l)=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=tt2n(ix,jy,1,n,l)*(1.+0.378*ew(td2n(ix,jy,1,n,l))/ &
179           psn(ix,jy,1,n,l))
180      pold=psn(ix,jy,1,n,l)
181      zold=0.
182      do kz=2,nuvz
183        pint=akz(kz)+bkz(kz)*psn(ix,jy,1,n,l)  ! pressure on model layers
184        tv=tthn(ix,jy,kz,n,l)*(1.+0.608*qvhn(ix,jy,kz,n,l))
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 (((tthn(ix,jy,kz,n,l)-tthn(ix,jy,lz,n,l))/ &
217                 (zlev(lz)-zlev(kz))).lt.0.002) then
218              tropopausen(ix,jy,1,n,l)=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_nests(l,n,uuhn,vvhn,pvhn)
236
237  end do
238
239
240end subroutine calcpar_nests
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG