19 use masks, only: dx, dy
20 use ctlblk_mod
, only: im, jsta_2l, jend_2u, jsta_m, jend_m, spval
21 use gridspec_mod
, only: gridtype
27 public :: ddvdx, ddudy, uuavg
29 public :: h2u, h2v, u2h, v2h
31 REAL,
ALLOCATABLE :: ddvdx(:,:)
32 REAL,
ALLOCATABLE :: ddudy(:,:)
33 REAL,
ALLOCATABLE :: uuavg(:,:)
40 subroutine dvdxdudy(uwnd,vwnd)
44 REAL,
dimension(im,jsta_2l:jend_2u),
intent(in) :: uwnd, vwnd
50 INTEGER,
allocatable :: ihe(:),ihw(:)
52 IF(gridtype ==
'A')
THEN
56 IF(vwnd(i+1,j)<spval.AND.vwnd(i-1,j)<spval.AND. &
57 & uwnd(i,j+1)<spval.AND.uwnd(i,j-1)<spval)
THEN
58 r2dx = 1./(2.*dx(i,j))
59 r2dy = 1./(2.*dy(i,j))
60 ddvdx(i,j) = (vwnd(i+1,j)-vwnd(i-1,j))*r2dx
61 ddudy(i,j) = (uwnd(i,j+1)-uwnd(i,j-1))*r2dy
62 uuavg(i,j) = uwnd(i,j)
66 ELSE IF (gridtype ==
'E')
THEN
67 allocate(ihw(jsta_2l:jend_2u), ihe(jsta_2l:jend_2u))
76 IF(vwnd(i+ihe(j),j) < spval.AND.vwnd(i+ihw(j),j) < spval .AND. &
77 & uwnd(i,j+1) < spval .AND.uwnd(i,j-1) < spval)
THEN
78 r2dx = 1./(2.*dx(i,j))
79 r2dy = 1./(2.*dy(i,j))
80 ddvdx(i,j) = (vwnd(i+ihe(j),j)-vwnd(i+ihw(j),j))*r2dx
81 ddudy(i,j) = (uwnd(i,j+1)-uwnd(i,j-1))*r2dy
82 uuavg(i,j) = 0.25*(uwnd(i+ihe(j),j)+uwnd(i+ihw(j),j) &
83 & + uwnd(i,j+1)+uwnd(i,j-1))
88 ELSE IF (gridtype ==
'B')
THEN
94 if(vwnd(i, j)==spval .or. vwnd(i, j-1)==spval .or. &
95 vwnd(i-1,j)==spval .or. vwnd(i-1,j-1)==spval .or. &
96 uwnd(i, j)==spval .or. uwnd(i-1,j )==spval .or. &
97 uwnd(i,j-1)==spval .or. uwnd(i-1,j-1)==spval) cycle
98 ddvdx(i,j) = (0.5*(vwnd(i,j)+vwnd(i,j-1))-0.5*(vwnd(i-1,j) &
99 & + vwnd(i-1,j-1)))*r2dx
100 ddudy(i,j) = (0.5*(uwnd(i,j)+uwnd(i-1,j))-0.5*(uwnd(i,j-1) &
101 & + uwnd(i-1,j-1)))*r2dy
102 uuavg(i,j) = 0.25*(uwnd(i-1,j-1)+uwnd(i-1,j) &
103 & + uwnd(i, j-1)+uwnd(i, j))
108 end subroutine dvdxdudy
112 subroutine h2u(ingrid,outgrid)
116 use ctlblk_mod
, only: spval, jsta, jend, jsta_m, jend, me, num_procs, jm,&
117 im, jsta_2l, jend_2u , jend_m
118 use gridspec_mod
, only: gridtype
124 real,
dimension(IM,JSTA_2L:JEND_2U),
intent(in)::ingrid
125 real,
dimension(IM,JSTA_2L:JEND_2U),
intent(out)::outgrid
127 if(gridtype ==
'A')
THEN
130 outgrid(i,j)=ingrid(i,j)
133 else IF(gridtype ==
'E')
THEN
134 call exch(ingrid(1,jsta_2l))
139 outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
142 ELSE IF(gridtype ==
'B')
THEN
143 call exch(ingrid(1,jsta_2l))
146 outgrid(i,j)=(ingrid(i,j)+ingrid(i,j+1)+ingrid(i+1,j)+ingrid(i+1,j+1))/4.0
151 outgrid(im,j)=outgrid(im-1,j)
153 IF(me == (num_procs-1) .and. jend_2u >= jm)
then
155 outgrid(i,jm) = outgrid(i,jm-1)
158 ELSE IF(gridtype ==
'C')
THEN
161 outgrid(i,j)=(ingrid(i,j)+ingrid(i+1,j))/2.0
170 subroutine h2v(ingrid,outgrid)
173 use ctlblk_mod
, only: spval, jsta, jend, jsta_m, jend_m, im, jsta_2l, jend_2u
174 use gridspec_mod
, only: gridtype
178 real,
dimension(IM,JSTA_2L:JEND_2U),
intent(in)::ingrid
179 real,
dimension(IM,JSTA_2L:JEND_2U),
intent(out)::outgrid
181 if(gridtype ==
'A')
THEN
184 outgrid(i,j)=ingrid(i,j)
187 else IF(gridtype ==
'E')
THEN
188 call exch(ingrid(1,jsta_2l))
193 outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
196 ELSE IF(gridtype ==
'B')
THEN
197 call exch(ingrid(1,jsta_2l))
200 outgrid(i,j)=(ingrid(i,j)+ingrid(i,j+1)+ingrid(i+1,j)+ingrid(i+1,j+1))/4.0
203 ELSE IF(gridtype ==
'C')
THEN
204 call exch(ingrid(1,jsta_2l))
207 outgrid(i,j)=(ingrid(i,j)+ingrid(i,j+1))/2.0
216 subroutine u2h(ingrid,outgrid)
219 use ctlblk_mod
, only: spval, jsta, jend, jsta_m, jend_m, im, jsta_2l, jend_2u
220 use gridspec_mod
, only: gridtype
224 real,
dimension(IM,JSTA_2L:JEND_2U),
intent(in)::ingrid
225 real,
dimension(IM,JSTA_2L:JEND_2U),
intent(out)::outgrid
227 if(gridtype ==
'A')
THEN
230 outgrid(i,j)=ingrid(i,j)
233 else IF(gridtype ==
'E')
THEN
234 call exch(ingrid(1,jsta_2l))
239 outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
242 ELSE IF(gridtype ==
'B')
THEN
243 call exch(ingrid(1,jsta_2l))
246 outgrid(i,j)=(ingrid(i-1,j-1)+ingrid(i,j-1)+ingrid(i-1,j)+ingrid(i,j))/4.0
249 ELSE IF(gridtype ==
'C')
THEN
252 outgrid(i,j)=(ingrid(i-1,j)+ingrid(i,j))/2.0
261 subroutine v2h(ingrid,outgrid)
264 use ctlblk_mod
, only: spval, jsta, jend, jsta_m, jend_m, im, jsta_2l, jend_2u
265 use gridspec_mod
, only: gridtype
269 real,
dimension(IM,JSTA_2L:JEND_2U),
intent(in)::ingrid
270 real,
dimension(IM,JSTA_2L:JEND_2U),
intent(out)::outgrid
272 if(gridtype ==
'A')
THEN
275 outgrid(i,j)=ingrid(i,j)
278 else IF(gridtype ==
'E')
THEN
279 call exch(ingrid(1,jsta_2l))
284 outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
287 ELSE IF(gridtype ==
'B')
THEN
288 call exch(ingrid(1,jsta_2l))
291 outgrid(i,j)=(ingrid(i-1,j-1)+ingrid(i,j-1)+ingrid(i-1,j)+ingrid(i,j))/4.0
294 ELSE IF(gridtype ==
'C')
THEN
295 call exch(ingrid(1,jsta_2l))
298 outgrid(i,j)=(ingrid(i,j-1)+ingrid(i,j))/2.0
dvdxdudy() computes dudy, dvdx, uwnd