source: branches/flexpart91_hasod/src_parallel/calcpv_nests.f90 @ 8

Last change on this file since 8 was 8, checked in by hasod, 11 years ago

Added parallel version of Flexpart91

File size: 11.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,ppmk
51  real :: ppml(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!$OMP PARALLEL DEFAULT(none) &
62!$OMP SHARED(l, nlck, nuvz, nxn, nyn, dxn, dyn, xlon0n, ylat0n,  &
63!$OMP   n, akz, bkz, psn, tthn, uuhn, vvhn, pvhn) &
64!$OMP PRIVATE(ix, jy, ii, i, phi, f, cosphi, tanphi, jyvp, jyvm, jumpy, juy, &
65!$OMP   ixvp, ixvm, jumpx, ivrp, ivrm, jux, ppml, ppmk, theta, klvrp, klvrm, klpt, &
66!$OMP   thetap, thetam, dthetadp, ivr, kup, kdn, kch, k, thdn, thup, dt1, dt2, dt, &
67!$OMP   vx, dvdx, jj, uy, dudy)
68
69#if (defined STATIC_SCHED)
70!$OMP DO SCHEDULE(static)
71#else
72!$OMP DO SCHEDULE(dynamic, max(1,(nyn(l)-1)/10))
73#endif
74  !
75  ! Loop over entire grid
76  !**********************
77  latloop : do jy=0,nyn(l)-1
78    phi = (ylat0n(l) + jy * dyn(l)) * pi / 180.
79    f = 0.00014585 * sin(phi)
80    tanphi = tan(phi)
81    cosphi = cos(phi)
82  ! Provide a virtual jy+1 and jy-1 in case we are on domain edge (Lat)
83      jyvp=jy+1
84      jyvm=jy-1
85      if (jy.eq.0) jyvm=0
86      if (jy.eq.nyn(l)-1) jyvp=nyn(l)-1
87  ! Define absolute gap length
88      jumpy=2
89      if (jy.eq.0.or.jy.eq.nyn(l)-1) jumpy=1
90      juy=jumpy
91  !
92    do ix=0,nxn(l)-1
93  ! Provide a virtual ix+1 and ix-1 in case we are on domain edge (Long)
94      ixvp=ix+1
95      ixvm=ix-1
96      jumpx=2
97      if (ix.eq.0) ixvm=0
98      if (ix.eq.nxn(l)-1) ixvp=nxn(l)-1
99      ivrp=ixvp
100      ivrm=ixvm
101  ! Define absolute gap length
102      if (ix.eq.0.or.ix.eq.nxn(l)-1) jumpx=1
103      jux=jumpx
104  ! Precalculate pressure values for efficiency
105      do kl=1,nuvz
106        ppml(kl)=akz(kl)+bkz(kl)*psn(ix,jy,1,n,l)
107      end do
108  !
109  ! Loop over the vertical
110  !***********************
111
112      do kl=1,nuvz
113        ppmk=akz(kl)+bkz(kl)*psn(ix,jy,1,n,l)
114        theta=tthn(ix,jy,kl,n,l)*(100000./ppmk)**kappa
115        klvrp=kl+1
116        klvrm=kl-1
117        klpt=kl
118  ! If top or bottom level, dthetadp is evaluated between the current
119  ! level and the level inside, otherwise between level+1 and level-1
120  !
121        if (klvrp.gt.nuvz) klvrp=nuvz
122        if (klvrm.lt.1) klvrm=1
123        ppmk=akz(klvrp)+bkz(klvrp)*psn(ix,jy,1,n,l)
124        thetap=tthn(ix,jy,klvrp,n,l)*(100000./ppmk)**kappa
125        ppmk=akz(klvrm)+bkz(klvrm)*psn(ix,jy,1,n,l)
126        thetam=tthn(ix,jy,klvrm,n,l)*(100000./ppmk)**kappa
127        dthetadp=(thetap-thetam)/(ppml(klvrp)-ppml(klvrm))
128
129  ! Compute vertical position at pot. temperature surface on subgrid
130  ! and the wind at that position
131  !*****************************************************************
132  ! a) in x direction
133        ii=0
134        do i=ixvm,ixvp,jumpx
135          ivr=i
136          ii=ii+1
137  ! Search adjacent levels for current theta value
138  ! Spiral out from current level for efficiency
139          kup=klpt-1
140          kdn=klpt
141          kch=0
14240        continue
143  ! Upward branch
144          kup=kup+1
145          if (kch.ge.nlck) goto 21     ! No more levels to check,
146  !                                       ! and no values found
147          if (kup.ge.nuvz) goto 41
148          kch=kch+1
149          k=kup
150          ppmk=akz(k)+bkz(k)*psn(ivr,jy,1,n,l)
151          thdn=tthn(ivr,jy,k,n,l)*(100000./ppmk)**kappa
152          ppmk=akz(k+1)+bkz(k+1)*psn(ivr,jy,1,n,l)
153          thup=tthn(ivr,jy,k+1,n,l)*(100000./ppmk)**kappa
154
155      if (((thdn.ge.theta).and.(thup.le.theta)).or. &
156           ((thdn.le.theta).and.(thup.ge.theta))) then
157              dt1=abs(theta-thdn)
158              dt2=abs(theta-thup)
159              dt=dt1+dt2
160              if (dt.lt.eps) then   ! Avoid division by zero error
161                dt1=0.5             ! G.W., 10.4.1996
162                dt2=0.5
163                dt=1.0
164              endif
165    vx(ii)=(vvhn(ivr,jy,k,l)*dt2+vvhn(ivr,jy,k+1,l)*dt1)/dt
166              goto 20
167            endif
16841        continue
169  ! Downward branch
170          kdn=kdn-1
171          if (kdn.lt.1) goto 40
172          kch=kch+1
173          k=kdn
174          ppmk=akz(k)+bkz(k)*psn(ivr,jy,1,n,l)
175          thdn=tthn(ivr,jy,k,n,l)*(100000./ppmk)**kappa
176          ppmk=akz(k+1)+bkz(k+1)*psn(ivr,jy,1,n,l)
177          thup=tthn(ivr,jy,k+1,n,l)*(100000./ppmk)**kappa
178      if (((thdn.ge.theta).and.(thup.le.theta)).or. &
179           ((thdn.le.theta).and.(thup.ge.theta))) then
180              dt1=abs(theta-thdn)
181              dt2=abs(theta-thup)
182              dt=dt1+dt2
183              if (dt.lt.eps) then   ! Avoid division by zero error
184                dt1=0.5             ! G.W., 10.4.1996
185                dt2=0.5
186                dt=1.0
187              endif
188    vx(ii)=(vvhn(ivr,jy,k,l)*dt2+vvhn(ivr,jy,k+1,l)*dt1)/dt
189              goto 20
190            endif
191            goto 40
192  ! This section used when no values were found
19321      continue
194  ! Must use vv at current level and long. jux becomes smaller by 1
195        vx(ii)=vvhn(ix,jy,kl,l)
196        jux=jux-1
197  ! Otherwise OK
19820        continue
199        end do
200      if (jux.gt.0) then
201      dvdx=(vx(2)-vx(1))/real(jux)/(dxn(l)*pi/180.)
202      else
203      dvdx=vvhn(ivrp,jy,kl,l)-vvhn(ivrm,jy,kl,l)
204      dvdx=dvdx/real(jumpx)/(dxn(l)*pi/180.)
205  ! Only happens if no equivalent theta value
206  ! can be found on either side, hence must use values
207  ! from either side, same pressure level.
208      end if
209
210  ! b) in y direction
211
212        jj=0
213        do j=jyvm,jyvp,jumpy
214          jj=jj+1
215  ! Search adjacent levels for current theta value
216  ! Spiral out from current level for efficiency
217          kup=klpt-1
218          kdn=klpt
219          kch=0
22070        continue
221  ! Upward branch
222          kup=kup+1
223          if (kch.ge.nlck) goto 51     ! No more levels to check,
224  !                                     ! and no values found
225          if (kup.ge.nuvz) goto 71
226          kch=kch+1
227          k=kup
228          ppmk=akz(k)+bkz(k)*psn(ix,j,1,n,l)
229          thdn=tthn(ix,j,k,n,l)*(100000./ppmk)**kappa
230          ppmk=akz(k+1)+bkz(k+1)*psn(ix,j,1,n,l)
231          thup=tthn(ix,j,k+1,n,l)*(100000./ppmk)**kappa
232      if (((thdn.ge.theta).and.(thup.le.theta)).or. &
233           ((thdn.le.theta).and.(thup.ge.theta))) then
234              dt1=abs(theta-thdn)
235              dt2=abs(theta-thup)
236              dt=dt1+dt2
237              if (dt.lt.eps) then   ! Avoid division by zero error
238                dt1=0.5             ! G.W., 10.4.1996
239                dt2=0.5
240                dt=1.0
241              endif
242        uy(jj)=(uuhn(ix,j,k,l)*dt2+uuhn(ix,j,k+1,l)*dt1)/dt
243              goto 50
244            endif
24571        continue
246  ! Downward branch
247          kdn=kdn-1
248          if (kdn.lt.1) goto 70
249          kch=kch+1
250          k=kdn
251          ppmk=akz(k)+bkz(k)*psn(ix,j,1,n,l)
252          thdn=tthn(ix,j,k,n,l)*(100000./ppmk)**kappa
253          ppmk=akz(k+1)+bkz(k+1)*psn(ix,j,1,n,l)
254          thup=tthn(ix,j,k+1,n,l)*(100000./ppmk)**kappa
255      if (((thdn.ge.theta).and.(thup.le.theta)).or. &
256           ((thdn.le.theta).and.(thup.ge.theta))) then
257              dt1=abs(theta-thdn)
258              dt2=abs(theta-thup)
259              dt=dt1+dt2
260              if (dt.lt.eps) then   ! Avoid division by zero error
261                dt1=0.5             ! G.W., 10.4.1996
262                dt2=0.5
263                dt=1.0
264              endif
265        uy(jj)=(uuhn(ix,j,k,l)*dt2+uuhn(ix,j,k+1,l)*dt1)/dt
266              goto 50
267            endif
268            goto 70
269  ! This section used when no values were found
27051      continue
271  ! Must use uu at current level and lat. juy becomes smaller by 1
272        uy(jj)=uuhn(ix,jy,kl,l)
273        juy=juy-1
274  ! Otherwise OK
27550        continue
276        end do
277      if (juy.gt.0) then
278      dudy=(uy(2)-uy(1))/real(juy)/(dyn(l)*pi/180.)
279      else
280      dudy=uuhn(ix,jyvp,kl,l)-uuhn(ix,jyvm,kl,l)
281      dudy=dudy/real(jumpy)/(dyn(l)*pi/180.)
282      end if
283  !
284      pvhn(ix,jy,kl,l)=dthetadp*(f+(dvdx/cosphi-dudy &
285           +uuhn(ix,jy,kl,l)*tanphi)/r_earth)*(-1.e6)*9.81
286
287  !
288  ! Resest jux and juy
289      jux=jumpx
290      juy=jumpy
291      end do
292    end do
293  end do latloop
294!$OMP END DO
295!$OMP END PARALLEL
296  !
297end subroutine calcpv_nests
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG