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