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