source: flexpart.git/src/calcpv.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 9.7 KB
Line 
1subroutine calcpv(n,uuh,vvh,pvh)
2  !               i  i   i   o
3  !*****************************************************************************
4  !                                                                            *
5  !  Calculation of potential vorticity on 3-d grid.                           *
6  !                                                                            *
7  !     Author: P. James                                                       *
8  !     3 February 2000                                                        *
9  !                                                                            *
10  !     Adaptation to FLEXPART, A. Stohl, 1 May 2000                           *
11  !                                                                            *
12  !*****************************************************************************
13  !                                                                            *
14  ! Variables:                                                                 *
15  ! n                  temporal index for meteorological fields (1 to 2)       *
16  !                                                                            *
17  ! Constants:                                                                 *
18  !                                                                            *
19  !*****************************************************************************
20
21  use par_mod
22  use com_mod
23
24  implicit none
25
26  integer :: n,ix,jy,i,j,k,kl,ii,jj,klvrp,klvrm,klpt,kup,kdn,kch
27  integer :: jyvp,jyvm,ixvp,ixvm,jumpx,jumpy,jux,juy,ivrm,ivrp,ivr
28  integer :: nlck
29  real :: vx(2),uy(2),phi,tanphi,cosphi,dvdx,dudy,f
30  real :: theta,thetap,thetam,dthetadp,dt1,dt2,dt
31  real :: pvavr,ppml(0:nxmax-1,0:nymax-1,nuvzmax),ppmk(0:nxmax-1,0:nymax-1,nuvzmax)
32  real :: thup,thdn
33  real,parameter :: eps=1.e-5, p0=101325
34  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
35  real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
36  real :: pvh(0:nxmax-1,0:nymax-1,nuvzmax)
37
38  ! Set number of levels to check for adjacent theta
39  nlck=nuvz/3
40  !
41  ! Loop over entire grid
42  !**********************
43  do kl=1,nuvz
44    do jy=0,nymin1
45      do ix=0,nxmin1
46         ppml(ix,jy,kl)=akz(kl)+bkz(kl)*ps(ix,jy,1,n)
47      enddo
48    enddo
49  enddo
50
51!  ppmk(:,:,1:nuvz)=(100000./ppml(:,:,1:nuvz))**kappa
52  ppmk(0:nxmin1,0:nymin1,1:nuvz)=(100000./ppml(0:nxmin1,0:nymin1,1:nuvz))**kappa
53
54  do jy=0,nymin1
55    if (sglobal.and.jy.eq.0) goto 10
56    if (nglobal.and.jy.eq.nymin1) goto 10
57    phi = (ylat0 + jy * dy) * pi / 180.
58    f = 0.00014585 * sin(phi)
59    tanphi = tan(phi)
60    cosphi = cos(phi)
61  ! Provide a virtual jy+1 and jy-1 in case we are on domain edge (Lat)
62      jyvp=jy+1
63      jyvm=jy-1
64      if (jy.eq.0) jyvm=0
65      if (jy.eq.nymin1) jyvp=nymin1
66  ! Define absolute gap length
67      jumpy=2
68      if (jy.eq.0.or.jy.eq.nymin1) jumpy=1
69      if (sglobal.and.jy.eq.1) then
70         jyvm=1
71         jumpy=1
72      end if
73      if (nglobal.and.jy.eq.ny-2) then
74         jyvp=ny-2
75         jumpy=1
76      end if
77      juy=jumpy
78  !
79    do ix=0,nxmin1
80  ! Provide a virtual ix+1 and ix-1 in case we are on domain edge (Long)
81      ixvp=ix+1
82      ixvm=ix-1
83      jumpx=2
84      if (xglobal) then
85         ivrp=ixvp
86         ivrm=ixvm
87         if (ixvm.lt.0) ivrm=ixvm+nxmin1
88         if (ixvp.ge.nx) ivrp=ixvp-nx+1
89      else
90        if (ix.eq.0) ixvm=0
91        if (ix.eq.nxmin1) ixvp=nxmin1
92        ivrp=ixvp
93        ivrm=ixvm
94  ! Define absolute gap length
95        if (ix.eq.0.or.ix.eq.nxmin1) jumpx=1
96      end if
97      jux=jumpx
98  !
99  ! Loop over the vertical
100  !***********************
101
102      do kl=1,nuvz
103        theta=tth(ix,jy,kl,n)*ppmk(ix,jy,kl)
104        klvrp=kl+1
105        klvrm=kl-1
106        klpt=kl
107  ! If top or bottom level, dthetadp is evaluated between the current
108  ! level and the level inside, otherwise between level+1 and level-1
109  !
110        if (klvrp.gt.nuvz) klvrp=nuvz
111        if (klvrm.lt.1) klvrm=1
112        thetap=tth(ix,jy,klvrp,n)*ppmk(ix,jy,klvrp)
113        thetam=tth(ix,jy,klvrm,n)*ppmk(ix,jy,klvrm)
114        dthetadp=(thetap-thetam)/(ppml(ix,jy,klvrp)-ppml(ix,jy,klvrm))
115
116  ! Compute vertical position at pot. temperature surface on subgrid
117  ! and the wind at that position
118  !*****************************************************************
119  ! a) in x direction
120        ii=0
121        do i=ixvm,ixvp,jumpx
122          ivr=i
123          if (xglobal) then
124             if (i.lt.0) ivr=ivr+nxmin1
125             if (i.ge.nx) ivr=ivr-nx+1
126          end if
127          ii=ii+1
128  ! Search adjacent levels for current theta value
129  ! Spiral out from current level for efficiency
130          kup=klpt-1
131          kdn=klpt
132          kch=0
13340        continue
134  ! Upward branch
135          kup=kup+1
136          if (kch.ge.nlck) goto 21     ! No more levels to check,
137  !                                       ! and no values found
138          if (kup.ge.nuvz) goto 41
139          kch=kch+1
140          k=kup
141          thdn=tth(ivr,jy,k,n)*ppmk(ivr,jy,k)
142          thup=tth(ivr,jy,k+1,n)*ppmk(ivr,jy,k+1)
143
144
145          if (((thdn.ge.theta).and.(thup.le.theta)).or. &
146          ((thdn.le.theta).and.(thup.ge.theta))) then
147            dt1=abs(theta-thdn)
148            dt2=abs(theta-thup)
149            dt=dt1+dt2
150            if (dt.lt.eps) then   ! Avoid division by zero error
151              dt1=0.5             ! G.W., 10.4.1996
152              dt2=0.5
153              dt=1.0
154            endif
155            vx(ii)=(vvh(ivr,jy,k)*dt2+vvh(ivr,jy,k+1)*dt1)/dt
156            goto 20
157          endif
15841        continue
159  ! Downward branch
160          kdn=kdn-1
161          if (kdn.lt.1) goto 40
162          kch=kch+1
163          k=kdn
164          thdn=tth(ivr,jy,k,n)*ppmk(ivr,jy,k)
165          thup=tth(ivr,jy,k+1,n)*ppmk(ivr,jy,k+1)
166
167          if (((thdn.ge.theta).and.(thup.le.theta)).or. &
168          ((thdn.le.theta).and.(thup.ge.theta))) then
169            dt1=abs(theta-thdn)
170            dt2=abs(theta-thup)
171            dt=dt1+dt2
172            if (dt.lt.eps) then   ! Avoid division by zero error
173              dt1=0.5             ! G.W., 10.4.1996
174              dt2=0.5
175              dt=1.0
176            endif
177            vx(ii)=(vvh(ivr,jy,k)*dt2+vvh(ivr,jy,k+1)*dt1)/dt
178            goto 20
179          endif
180          goto 40
181  ! This section used when no values were found
18221      continue
183  ! Must use vv at current level and long. jux becomes smaller by 1
184        vx(ii)=vvh(ix,jy,kl)
185        jux=jux-1
186  ! Otherwise OK
18720        continue
188        end do
189      if (jux.gt.0) then
190      dvdx=(vx(2)-vx(1))/real(jux)/(dx*pi/180.)
191      else
192      dvdx=vvh(ivrp,jy,kl)-vvh(ivrm,jy,kl)
193      dvdx=dvdx/real(jumpx)/(dx*pi/180.)
194  ! Only happens if no equivalent theta value
195  ! can be found on either side, hence must use values
196  ! from either side, same pressure level.
197      end if
198
199  ! b) in y direction
200
201        jj=0
202        do j=jyvm,jyvp,jumpy
203          jj=jj+1
204  ! Search adjacent levels for current theta value
205  ! Spiral out from current level for efficiency
206          kup=klpt-1
207          kdn=klpt
208          kch=0
20970        continue
210  ! Upward branch
211          kup=kup+1
212          if (kch.ge.nlck) goto 51     ! No more levels to check,
213  !                                     ! and no values found
214          if (kup.ge.nuvz) goto 71
215          kch=kch+1
216          k=kup
217          thdn=tth(ix,j,k,n)*ppmk(ix,j,k)
218          thup=tth(ix,j,k+1,n)*ppmk(ix,j,k+1)
219          if (((thdn.ge.theta).and.(thup.le.theta)).or. &
220          ((thdn.le.theta).and.(thup.ge.theta))) then
221            dt1=abs(theta-thdn)
222            dt2=abs(theta-thup)
223            dt=dt1+dt2
224            if (dt.lt.eps) then   ! Avoid division by zero error
225              dt1=0.5             ! G.W., 10.4.1996
226              dt2=0.5
227              dt=1.0
228            endif
229            uy(jj)=(uuh(ix,j,k)*dt2+uuh(ix,j,k+1)*dt1)/dt
230            goto 50
231          endif
23271        continue
233  ! Downward branch
234          kdn=kdn-1
235          if (kdn.lt.1) goto 70
236          kch=kch+1
237          k=kdn
238          thdn=tth(ix,j,k,n)*ppmk(ix,j,k)
239          thup=tth(ix,j,k+1,n)*ppmk(ix,j,k+1)
240          if (((thdn.ge.theta).and.(thup.le.theta)).or. &
241          ((thdn.le.theta).and.(thup.ge.theta))) then
242            dt1=abs(theta-thdn)
243            dt2=abs(theta-thup)
244            dt=dt1+dt2
245            if (dt.lt.eps) then   ! Avoid division by zero error
246              dt1=0.5             ! G.W., 10.4.1996
247              dt2=0.5
248              dt=1.0
249            endif
250            uy(jj)=(uuh(ix,j,k)*dt2+uuh(ix,j,k+1)*dt1)/dt
251            goto 50
252          endif
253          goto 70
254  ! This section used when no values were found
25551      continue
256  ! Must use uu at current level and lat. juy becomes smaller by 1
257          uy(jj)=uuh(ix,jy,kl)
258          juy=juy-1
259  ! Otherwise OK
26050        continue
261        end do
262      if (juy.gt.0) then
263      dudy=(uy(2)-uy(1))/real(juy)/(dy*pi/180.)
264      else
265      dudy=uuh(ix,jyvp,kl)-uuh(ix,jyvm,kl)
266      dudy=dudy/real(jumpy)/(dy*pi/180.)
267      end if
268  !
269      pvh(ix,jy,kl)=dthetadp*(f+(dvdx/cosphi-dudy &
270           +uuh(ix,jy,kl)*tanphi)/r_earth)*(-1.e6)*9.81
271
272
273  !
274  ! Resest jux and juy
275      jux=jumpx
276      juy=jumpy
277      end do
278    end do
27910  continue
280  end do
281  !
282  ! Fill in missing PV values on poles, if present
283  ! Use mean PV of surrounding latitude ring
284  !
285      if (sglobal) then
286         do kl=1,nuvz
287            pvavr=0.
288            do ix=0,nxmin1
289               pvavr=pvavr+pvh(ix,1,kl)
290            end do
291            pvavr=pvavr/real(nx)
292            jy=0
293            do ix=0,nxmin1
294               pvh(ix,jy,kl)=pvavr
295            end do
296         end do
297      end if
298      if (nglobal) then
299         do kl=1,nuvz
300            pvavr=0.
301            do ix=0,nxmin1
302               pvavr=pvavr+pvh(ix,ny-2,kl)
303            end do
304            pvavr=pvavr/real(nx)
305            jy=nymin1
306            do ix=0,nxmin1
307               pvh(ix,jy,kl)=pvavr
308            end do
309         end do
310      end if
311
312end subroutine calcpv
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG