source: flexpart.git/src/calcpv_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: 10.0 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 calcpv_nests(l,n,uuhn,vvhn,pvhn)
23  !                     i i  i    i    o
24  !*****************************************************************************
25  !                                                                            *
26  !  Calculation of potential vorticity on 3-d nested grid                     *
27  !                                                                            *
28  !     Author: P. James                                                       *
29  !     22 February 2000                                                       *
30  !                                                                            *
31  !*****************************************************************************
32  !                                                                            *
33  ! Variables:                                                                 *
34  ! n                  temporal index for meteorological fields (1 to 2)       *
35  ! l                  index of current nest                                   *
36  !                                                                            *
37  ! Constants:                                                                 *
38  !                                                                            *
39  !*****************************************************************************
40
41  use par_mod
42  use com_mod
43
44  implicit none
45
46  integer :: n,ix,jy,i,j,k,kl,ii,jj,klvrp,klvrm,klpt,kup,kdn,kch
47  integer :: jyvp,jyvm,ixvp,ixvm,jumpx,jumpy,jux,juy,ivrm,ivrp,ivr
48  integer :: nlck,l
49  real :: vx(2),uy(2),phi,tanphi,cosphi,dvdx,dudy,f
50  real :: theta,thetap,thetam,dthetadp,dt1,dt2,dt
51  real :: ppml(0:nxmaxn-1,0:nymaxn-1,nuvzmax),ppmk(0:nxmaxn-1,0:nymaxn-1,nuvzmax)
52  real :: thup,thdn
53  real,parameter :: eps=1.e-5,p0=101325
54  real :: uuhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
55  real :: vvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
56  real :: pvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
57
58  ! Set number of levels to check for adjacent theta
59  nlck=nuvz/3
60  !
61  ! Loop over entire grid
62  !**********************
63  do kl=1,nuvz
64    do jy=0,nyn(l)-1
65      do ix=0,nxn(l)-1
66         ppml(ix,jy,kl)=akz(kl)+bkz(kl)*psn(ix,jy,1,n,l)
67      enddo
68    enddo
69  enddo
70!  ppmk=(100000./ppml)**kappa
71  ppmk(0:nxn(l)-1,0:nyn(l)-1,1:nuvz)=(100000./ppml(0:nxn(l)-1,0:nyn(l)-1,1:nuvz))**kappa
72
73  do jy=0,nyn(l)-1
74    phi = (ylat0n(l) + jy * dyn(l)) * pi / 180.
75    f = 0.00014585 * sin(phi)
76    tanphi = tan(phi)
77    cosphi = cos(phi)
78  ! Provide a virtual jy+1 and jy-1 in case we are on domain edge (Lat)
79      jyvp=jy+1
80      jyvm=jy-1
81      if (jy.eq.0) jyvm=0
82      if (jy.eq.nyn(l)-1) jyvp=nyn(l)-1
83  ! Define absolute gap length
84      jumpy=2
85      if (jy.eq.0.or.jy.eq.nyn(l)-1) jumpy=1
86      juy=jumpy
87  !
88    do ix=0,nxn(l)-1
89  ! Provide a virtual ix+1 and ix-1 in case we are on domain edge (Long)
90      ixvp=ix+1
91      ixvm=ix-1
92      jumpx=2
93      if (ix.eq.0) ixvm=0
94      if (ix.eq.nxn(l)-1) ixvp=nxn(l)-1
95      ivrp=ixvp
96      ivrm=ixvm
97  ! Define absolute gap length
98      if (ix.eq.0.or.ix.eq.nxn(l)-1) jumpx=1
99      jux=jumpx
100  !
101  ! Loop over the vertical
102  !***********************
103
104      do kl=1,nuvz
105        theta=tthn(ix,jy,kl,n,l)*ppmk(ix,jy,kl)
106        klvrp=kl+1
107        klvrm=kl-1
108        klpt=kl
109  ! If top or bottom level, dthetadp is evaluated between the current
110  ! level and the level inside, otherwise between level+1 and level-1
111  !
112        if (klvrp.gt.nuvz) klvrp=nuvz
113        if (klvrm.lt.1) klvrm=1
114        thetap=tthn(ix,jy,klvrp,n,l)*ppmk(ix,jy,klvrp)
115        thetam=tthn(ix,jy,klvrm,n,l)*ppmk(ix,jy,klvrm)
116        dthetadp=(thetap-thetam)/(ppml(ix,jy,klvrp)-ppml(ix,jy,klvrm))
117
118  ! Compute vertical position at pot. temperature surface on subgrid
119  ! and the wind at that position
120  !*****************************************************************
121  ! a) in x direction
122        ii=0
123        do i=ixvm,ixvp,jumpx
124          ivr=i
125          ii=ii+1
126  ! Search adjacent levels for current theta value
127  ! Spiral out from current level for efficiency
128          kup=klpt-1
129          kdn=klpt
130          kch=0
13140        continue
132  ! Upward branch
133          kup=kup+1
134          if (kch.ge.nlck) goto 21     ! No more levels to check,
135  !                                       ! and no values found
136          if (kup.ge.nuvz) goto 41
137          kch=kch+1
138          k=kup
139          thdn=tthn(ivr,jy,k,n,l)*ppmk(ivr,jy,k)
140          thup=tthn(ivr,jy,k+1,n,l)*ppmk(ivr,jy,k+1)
141
142      if (((thdn.ge.theta).and.(thup.le.theta)).or. &
143           ((thdn.le.theta).and.(thup.ge.theta))) then
144              dt1=abs(theta-thdn)
145              dt2=abs(theta-thup)
146              dt=dt1+dt2
147              if (dt.lt.eps) then   ! Avoid division by zero error
148                dt1=0.5             ! G.W., 10.4.1996
149                dt2=0.5
150                dt=1.0
151              endif
152    vx(ii)=(vvhn(ivr,jy,k,l)*dt2+vvhn(ivr,jy,k+1,l)*dt1)/dt
153              goto 20
154            endif
15541        continue
156  ! Downward branch
157          kdn=kdn-1
158          if (kdn.lt.1) goto 40
159          kch=kch+1
160          k=kdn
161          thdn=tthn(ivr,jy,k,n,l)*ppmk(ivr,jy,k)
162          thup=tthn(ivr,jy,k+1,n,l)*ppmk(ivr,jy,k+1)
163      if (((thdn.ge.theta).and.(thup.le.theta)).or. &
164           ((thdn.le.theta).and.(thup.ge.theta))) then
165              dt1=abs(theta-thdn)
166              dt2=abs(theta-thup)
167              dt=dt1+dt2
168              if (dt.lt.eps) then   ! Avoid division by zero error
169                dt1=0.5             ! G.W., 10.4.1996
170                dt2=0.5
171                dt=1.0
172              endif
173    vx(ii)=(vvhn(ivr,jy,k,l)*dt2+vvhn(ivr,jy,k+1,l)*dt1)/dt
174              goto 20
175            endif
176            goto 40
177  ! This section used when no values were found
17821      continue
179  ! Must use vv at current level and long. jux becomes smaller by 1
180        vx(ii)=vvhn(ix,jy,kl,l)
181        jux=jux-1
182  ! Otherwise OK
18320        continue
184        end do
185      if (jux.gt.0) then
186      dvdx=(vx(2)-vx(1))/real(jux)/(dxn(l)*pi/180.)
187      else
188      dvdx=vvhn(ivrp,jy,kl,l)-vvhn(ivrm,jy,kl,l)
189      dvdx=dvdx/real(jumpx)/(dxn(l)*pi/180.)
190  ! Only happens if no equivalent theta value
191  ! can be found on either side, hence must use values
192  ! from either side, same pressure level.
193      end if
194
195  ! b) in y direction
196
197        jj=0
198        do j=jyvm,jyvp,jumpy
199          jj=jj+1
200  ! Search adjacent levels for current theta value
201  ! Spiral out from current level for efficiency
202          kup=klpt-1
203          kdn=klpt
204          kch=0
20570        continue
206  ! Upward branch
207          kup=kup+1
208          if (kch.ge.nlck) goto 51     ! No more levels to check,
209  !                                     ! and no values found
210          if (kup.ge.nuvz) goto 71
211          kch=kch+1
212          k=kup
213          thdn=tthn(ix,j,k,n,l)*ppmk(ix,j,k)
214          thup=tthn(ix,j,k+1,n,l)*ppmk(ix,j,k+1)
215      if (((thdn.ge.theta).and.(thup.le.theta)).or. &
216           ((thdn.le.theta).and.(thup.ge.theta))) then
217              dt1=abs(theta-thdn)
218              dt2=abs(theta-thup)
219              dt=dt1+dt2
220              if (dt.lt.eps) then   ! Avoid division by zero error
221                dt1=0.5             ! G.W., 10.4.1996
222                dt2=0.5
223                dt=1.0
224              endif
225        uy(jj)=(uuhn(ix,j,k,l)*dt2+uuhn(ix,j,k+1,l)*dt1)/dt
226              goto 50
227            endif
22871        continue
229  ! Downward branch
230          kdn=kdn-1
231          if (kdn.lt.1) goto 70
232          kch=kch+1
233          k=kdn
234          thdn=tthn(ix,j,k,n,l)*ppmk(ix,j,k)
235          thup=tthn(ix,j,k+1,n,l)*ppmk(ix,j,k+1)
236      if (((thdn.ge.theta).and.(thup.le.theta)).or. &
237           ((thdn.le.theta).and.(thup.ge.theta))) then
238              dt1=abs(theta-thdn)
239              dt2=abs(theta-thup)
240              dt=dt1+dt2
241              if (dt.lt.eps) then   ! Avoid division by zero error
242                dt1=0.5             ! G.W., 10.4.1996
243                dt2=0.5
244                dt=1.0
245              endif
246        uy(jj)=(uuhn(ix,j,k,l)*dt2+uuhn(ix,j,k+1,l)*dt1)/dt
247              goto 50
248            endif
249            goto 70
250  ! This section used when no values were found
25151      continue
252  ! Must use uu at current level and lat. juy becomes smaller by 1
253        uy(jj)=uuhn(ix,jy,kl,l)
254        juy=juy-1
255  ! Otherwise OK
25650        continue
257        end do
258      if (juy.gt.0) then
259      dudy=(uy(2)-uy(1))/real(juy)/(dyn(l)*pi/180.)
260      else
261      dudy=uuhn(ix,jyvp,kl,l)-uuhn(ix,jyvm,kl,l)
262      dudy=dudy/real(jumpy)/(dyn(l)*pi/180.)
263      end if
264  !
265      pvhn(ix,jy,kl,l)=dthetadp*(f+(dvdx/cosphi-dudy &
266           +uuhn(ix,jy,kl,l)*tanphi)/r_earth)*(-1.e6)*9.81
267
268  !
269  ! Resest jux and juy
270      jux=jumpx
271      juy=jumpy
272      end do
273    end do
274  end do
275  !
276end subroutine calcpv_nests
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG