22 use ctlblk_mod,
only: im, jsta_2l, jend_2u, jsta_m, jend_m, spval,&
54 subroutine dvdxdudy(uwnd,vwnd)
58 REAL,
dimension(ista_2l:iend_2u,jsta_2l:jend_2u),
intent(in) :: uwnd, vwnd
64 INTEGER,
allocatable :: ihe(:),ihw(:)
75 IF(gridtype ==
'A')
THEN
79 IF(vwnd(i+1,j)<spval.AND.vwnd(i-1,j)<spval.AND. &
80 uwnd(i,j+1)<spval.AND.uwnd(i,j-1)<spval.AND. &
81 uwnd(i,j)<spval.AND.vwnd(i,j)<spval.AND. &
82 abs(dx(i,j))>1.e-5.AND.abs(dy(i,j))>1.e-5)
THEN
83 r2dx = 1./(2.*dx(i,j))
84 r2dy = 1./(2.*dy(i,j))
85 ddvdx(i,j) = (vwnd(i+1,j)-vwnd(i-1,j))*r2dx
86 ddudy(i,j) = (uwnd(i,j+1)-uwnd(i,j-1))*r2dy
87 uuavg(i,j) = uwnd(i,j)
91 ELSE IF (gridtype ==
'E')
THEN
92 allocate(ihw(jsta_2l:jend_2u), ihe(jsta_2l:jend_2u))
101 IF(vwnd(i+ihe(j),j) < spval.AND.vwnd(i+ihw(j),j) < spval .AND. &
102 & uwnd(i,j+1) < spval .AND.uwnd(i,j-1) < spval)
THEN
103 r2dx = 1./(2.*dx(i,j))
104 r2dy = 1./(2.*dy(i,j))
105 ddvdx(i,j) = (vwnd(i+ihe(j),j)-vwnd(i+ihw(j),j))*r2dx
106 ddudy(i,j) = (uwnd(i,j+1)-uwnd(i,j-1))*r2dy
107 uuavg(i,j) = 0.25*(uwnd(i+ihe(j),j)+uwnd(i+ihw(j),j) &
108 & + uwnd(i,j+1)+uwnd(i,j-1))
113 ELSE IF (gridtype ==
'B')
THEN
119 if(vwnd(i, j)==spval .or. vwnd(i, j-1)==spval .or. &
120 vwnd(i-1,j)==spval .or. vwnd(i-1,j-1)==spval .or. &
121 uwnd(i, j)==spval .or. uwnd(i-1,j )==spval .or. &
122 uwnd(i,j-1)==spval .or. uwnd(i-1,j-1)==spval) cycle
123 ddvdx(i,j) = (0.5*(vwnd(i,j)+vwnd(i,j-1))-0.5*(vwnd(i-1,j) &
124 & + vwnd(i-1,j-1)))*r2dx
125 ddudy(i,j) = (0.5*(uwnd(i,j)+uwnd(i-1,j))-0.5*(uwnd(i,j-1) &
126 & + uwnd(i-1,j-1)))*r2dy
127 uuavg(i,j) = 0.25*(uwnd(i-1,j-1)+uwnd(i-1,j) &
128 & + uwnd(i, j-1)+uwnd(i, j))
137 subroutine h2u(ingrid,outgrid)
147 use ctlblk_mod,
only: spval, jsta, jend, jsta_m, jend_m, me, num_procs, jm,&
148 im, jsta_2l, jend_2u, ista, iend, ista_m, iend_m, ista_2l, iend_2u
149 use gridspec_mod,
only: gridtype
155 real,
dimension(ISTA_2L:IEND_2U,JSTA_2L:JEND_2U),
intent(in)::ingrid
156 real,
dimension(ISTA_2L:IEND_2U,JSTA_2L:JEND_2U),
intent(out)::outgrid
158 if(gridtype ==
'A')
THEN
161 outgrid(i,j)=ingrid(i,j)
164 else IF(gridtype ==
'E')
THEN
165 call exch(ingrid(ista_2l,jsta_2l))
170 outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
173 ELSE IF(gridtype ==
'B')
THEN
174 call exch(ingrid(ista_2l,jsta_2l))
177 outgrid(i,j)=(ingrid(i,j)+ingrid(i,j+1)+ingrid(i+1,j)+ingrid(i+1,j+1))/4.0
182 outgrid(iend,j)=outgrid(iend-1,j)
184 IF(me == (num_procs-1) .and. jend_2u >= jm)
then
186 outgrid(i,jend) = outgrid(i,jend-1)
189 ELSE IF(gridtype ==
'C')
THEN
190 call exch(ingrid(ista_2l,jsta_2l))
193 outgrid(i,j)=(ingrid(i,j)+ingrid(i+1,j))/2.0
202 subroutine h2v(ingrid,outgrid)
211 use ctlblk_mod,
only: spval, jsta, jend, jsta_m, jend_m, im, jsta_2l, jend_2u,&
212 ista, iend, ista_m, iend_m, ista_2l, iend_2u
213 use gridspec_mod,
only: gridtype
217 real,
dimension(ISTA_2L:IEND_2U,JSTA_2L:JEND_2U),
intent(in)::ingrid
218 real,
dimension(ISTA_2L:IEND_2U,JSTA_2L:JEND_2U),
intent(out)::outgrid
220 if(gridtype ==
'A')
THEN
223 outgrid(i,j)=ingrid(i,j)
226 else IF(gridtype ==
'E')
THEN
227 call exch(ingrid(ista_2l,jsta_2l))
232 outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
235 ELSE IF(gridtype ==
'B')
THEN
236 call exch(ingrid(ista_2l,jsta_2l))
239 outgrid(i,j)=(ingrid(i,j)+ingrid(i,j+1)+ingrid(i+1,j)+ingrid(i+1,j+1))/4.0
242 ELSE IF(gridtype ==
'C')
THEN
243 call exch(ingrid(ista_2l,jsta_2l))
246 outgrid(i,j)=(ingrid(i,j)+ingrid(i,j+1))/2.0
255 subroutine u2h(ingrid,outgrid)
264 use ctlblk_mod,
only: spval, jsta, jend, jsta_m, jend_m, im, jsta_2l, jend_2u,&
265 ista, iend, ista_m, iend_m, ista_2l, iend_2u
266 use gridspec_mod,
only: gridtype
270 real,
dimension(ISTA_2L:IEND_2U,JSTA_2L:JEND_2U),
intent(in)::ingrid
271 real,
dimension(ISTA_2L:IEND_2U,JSTA_2L:JEND_2U),
intent(out)::outgrid
273 if(gridtype ==
'A')
THEN
276 outgrid(i,j)=ingrid(i,j)
279 else IF(gridtype ==
'E')
THEN
280 call exch(ingrid(ista_2l,jsta_2l))
285 outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
288 ELSE IF(gridtype ==
'B')
THEN
289 call exch(ingrid(ista_2l,jsta_2l))
292 outgrid(i,j)=(ingrid(i-1,j-1)+ingrid(i,j-1)+ingrid(i-1,j)+ingrid(i,j))/4.0
295 ELSE IF(gridtype ==
'C')
THEN
296 call exch(ingrid(ista_2l,jsta_2l))
299 outgrid(i,j)=(ingrid(i-1,j)+ingrid(i,j))/2.0
308 subroutine v2h(ingrid,outgrid)
317 use ctlblk_mod,
only: spval, jsta, jend, jsta_m, jend_m, im, jsta_2l, jend_2u,&
318 ista, iend, ista_m, iend_m, ista_2l, iend_2u
319 use gridspec_mod,
only: gridtype
323 real,
dimension(ISTA_2L:IEND_2U,JSTA_2L:JEND_2U),
intent(in)::ingrid
324 real,
dimension(ISTA_2L:IEND_2U,JSTA_2L:JEND_2U),
intent(out)::outgrid
326 if(gridtype ==
'A')
THEN
329 outgrid(i,j)=ingrid(i,j)
332 else IF(gridtype ==
'E')
THEN
333 call exch(ingrid(ista_2l,jsta_2l))
338 outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
341 ELSE IF(gridtype ==
'B')
THEN
342 call exch(ingrid(ista_2l,jsta_2l))
345 outgrid(i,j)=(ingrid(i-1,j-1)+ingrid(i,j-1)+ingrid(i-1,j)+ingrid(i,j))/4.0
348 ELSE IF(gridtype ==
'C')
THEN
349 call exch(ingrid(ista_2l,jsta_2l))
352 outgrid(i,j)=(ingrid(i,j-1)+ingrid(i,j))/2.0