UPP  001
 All Data Structures Files Functions Pages
UPP_MATH.f
Go to the documentation of this file.
1 
5 
17  module upp_math
18 
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
22 !
23  implicit none
24 
25  private
26 
27  public :: ddvdx, ddudy, uuavg
28  public :: dvdxdudy
29  public :: h2u, h2v, u2h, v2h
30 
31  REAL, ALLOCATABLE :: ddvdx(:,:)
32  REAL, ALLOCATABLE :: ddudy(:,:)
33  REAL, ALLOCATABLE :: uuavg(:,:)
34 !
35 !-------------------------------------------------------------------------------------
36 !
37  contains
38 !
39 !-------------------------------------------------------------------------------------
40  subroutine dvdxdudy(uwnd,vwnd)
41 !
42  implicit none
43 
44  REAL, dimension(im,jsta_2l:jend_2u), intent(in) :: uwnd, vwnd
45 !
46 !-- local variables
47 !--
48  integer i, j
49  real r2dx, r2dy
50  INTEGER, allocatable :: ihe(:),ihw(:)
51 !
52  IF(gridtype == 'A')THEN
53 !$omp parallel do private(i,j,r2dx,r2dy)
54  DO j=jsta_m,jend_m
55  DO i=2,im-1
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)
63  END IF
64  END DO
65  END DO
66  ELSE IF (gridtype == 'E')THEN
67  allocate(ihw(jsta_2l:jend_2u), ihe(jsta_2l:jend_2u))
68 !$omp parallel do private(j)
69  DO j=jsta_2l,jend_2u
70  ihw(j) = -mod(j,2)
71  ihe(j) = ihw(j)+1
72  ENDDO
73 !$omp parallel do private(i,j,r2dx,r2dy)
74  DO j=jsta_m,jend_m
75  DO i=2,im-1
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))
84  END IF
85  END DO
86  END DO
87  deallocate(ihw, ihe)
88  ELSE IF (gridtype == 'B')THEN
89 !$omp parallel do private(i,j,r2dx,r2dy)
90  DO j=jsta_m,jend_m
91  DO i=2,im-1
92  r2dx = 1./dx(i,j)
93  r2dy = 1./dy(i,j)
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))
104  END DO
105  END DO
106  END IF
107 
108  end subroutine dvdxdudy
109 !
110 !-------------------------------------------------------------------------------------
111 !
112  subroutine h2u(ingrid,outgrid)
113 ! This subroutine interpolates from H points onto U points
114 ! Author: CHUANG, EMC, Dec. 2010
115 
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
119 
120  implicit none
121 
122  include "mpif.h"
123  integer:: i,j,ie,iw
124  real,dimension(IM,JSTA_2L:JEND_2U),intent(in)::ingrid
125  real,dimension(IM,JSTA_2L:JEND_2U),intent(out)::outgrid
126  outgrid=spval
127  if(gridtype == 'A')THEN
128  do j=jsta,jend
129  do i=1,im
130  outgrid(i,j)=ingrid(i,j)
131  end do
132  end do
133  else IF(gridtype == 'E')THEN
134  call exch(ingrid(1,jsta_2l))
135  DO j=jsta_m,jend_m
136  DO i=2,im-1
137  ie=i+mod(j,2)
138  iw=ie-1
139  outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
140  end do
141  end do
142  ELSE IF(gridtype == 'B')THEN
143  call exch(ingrid(1,jsta_2l))
144  DO j=jsta,jend_m
145  DO i=1,im-1
146  outgrid(i,j)=(ingrid(i,j)+ingrid(i,j+1)+ingrid(i+1,j)+ingrid(i+1,j+1))/4.0
147  end do
148  end do
149 ! Fill in boundary points because hysplit fails when 10 m wind has bitmaps
150  do j=jsta,jend_m
151  outgrid(im,j)=outgrid(im-1,j)
152  end do
153  IF(me == (num_procs-1) .and. jend_2u >= jm) then
154  DO i=1,im
155  outgrid(i,jm) = outgrid(i,jm-1)
156  END DO
157  END IF
158  ELSE IF(gridtype == 'C')THEN
159  DO j=jsta,jend
160  DO i=1,im-1
161  outgrid(i,j)=(ingrid(i,j)+ingrid(i+1,j))/2.0
162  end do
163  end do
164  end if
165 
166  end subroutine h2u
167 !
168 !-------------------------------------------------------------------------------------
169 !
170  subroutine h2v(ingrid,outgrid)
171 ! This subroutine interpolates from H points onto V points
172 ! Author: CHUANG, EMC, Dec. 2010
173  use ctlblk_mod, only: spval, jsta, jend, jsta_m, jend_m, im, jsta_2l, jend_2u
174  use gridspec_mod, only: gridtype
175  implicit none
176  include "mpif.h"
177  integer:: i,j,ie,iw
178  real,dimension(IM,JSTA_2L:JEND_2U),intent(in)::ingrid
179  real,dimension(IM,JSTA_2L:JEND_2U),intent(out)::outgrid
180  outgrid=spval
181  if(gridtype == 'A')THEN
182  do j=jsta,jend
183  do i=1,im
184  outgrid(i,j)=ingrid(i,j)
185  end do
186  end do
187  else IF(gridtype == 'E')THEN
188  call exch(ingrid(1,jsta_2l))
189  DO j=jsta_m,jend_m
190  DO i=2,im-1
191  ie=i+mod(j,2)
192  iw=ie-1
193  outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
194  end do
195  end do
196  ELSE IF(gridtype == 'B')THEN
197  call exch(ingrid(1,jsta_2l))
198  DO j=jsta,jend_m
199  DO i=1,im-1
200  outgrid(i,j)=(ingrid(i,j)+ingrid(i,j+1)+ingrid(i+1,j)+ingrid(i+1,j+1))/4.0
201  end do
202  end do
203  ELSE IF(gridtype == 'C')THEN
204  call exch(ingrid(1,jsta_2l))
205  DO j=jsta,jend_m
206  DO i=1,im
207  outgrid(i,j)=(ingrid(i,j)+ingrid(i,j+1))/2.0
208  end do
209  end do
210  end if
211 
212  end subroutine h2v
213 !
214 !-------------------------------------------------------------------------------------
215 !
216  subroutine u2h(ingrid,outgrid)
217 ! This subroutine interpolates from U points onto H points
218 ! Author: CHUANG, EMC, Dec. 2010
219  use ctlblk_mod, only: spval, jsta, jend, jsta_m, jend_m, im, jsta_2l, jend_2u
220  use gridspec_mod, only: gridtype
221  implicit none
222  include "mpif.h"
223  integer:: i,j,ie,iw
224  real,dimension(IM,JSTA_2L:JEND_2U),intent(in)::ingrid
225  real,dimension(IM,JSTA_2L:JEND_2U),intent(out)::outgrid
226  outgrid=spval
227  if(gridtype == 'A')THEN
228  do j=jsta,jend
229  do i=1,im
230  outgrid(i,j)=ingrid(i,j)
231  end do
232  end do
233  else IF(gridtype == 'E')THEN
234  call exch(ingrid(1,jsta_2l))
235  DO j=jsta_m,jend_m
236  DO i=2,im-1
237  ie=i+mod(j+1,2)
238  iw=ie-1
239  outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
240  end do
241  end do
242  ELSE IF(gridtype == 'B')THEN
243  call exch(ingrid(1,jsta_2l))
244  DO j=jsta_m,jend_m
245  DO i=2,im-1
246  outgrid(i,j)=(ingrid(i-1,j-1)+ingrid(i,j-1)+ingrid(i-1,j)+ingrid(i,j))/4.0
247  end do
248  end do
249  ELSE IF(gridtype == 'C')THEN
250  DO j=jsta,jend
251  DO i=2,im
252  outgrid(i,j)=(ingrid(i-1,j)+ingrid(i,j))/2.0
253  end do
254  end do
255  end if
256 
257  end subroutine u2h
258 !
259 !-------------------------------------------------------------------------------------
260 !
261  subroutine v2h(ingrid,outgrid)
262 ! This subroutine interpolates from V points onto H points
263 ! Author: CHUANG, EMC, Dec. 2010
264  use ctlblk_mod, only: spval, jsta, jend, jsta_m, jend_m, im, jsta_2l, jend_2u
265  use gridspec_mod, only: gridtype
266  implicit none
267  include "mpif.h"
268  integer:: i,j,ie,iw
269  real,dimension(IM,JSTA_2L:JEND_2U),intent(in)::ingrid
270  real,dimension(IM,JSTA_2L:JEND_2U),intent(out)::outgrid
271  outgrid=spval
272  if(gridtype == 'A')THEN
273  do j=jsta,jend
274  do i=1,im
275  outgrid(i,j)=ingrid(i,j)
276  end do
277  end do
278  else IF(gridtype == 'E')THEN
279  call exch(ingrid(1,jsta_2l))
280  DO j=jsta_m,jend_m
281  DO i=2,im-1
282  ie=i+mod(j,2)
283  iw=ie-1
284  outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
285  end do
286  end do
287  ELSE IF(gridtype == 'B')THEN
288  call exch(ingrid(1,jsta_2l))
289  DO j=jsta_m,jend_m
290  DO i=2,im-1
291  outgrid(i,j)=(ingrid(i-1,j-1)+ingrid(i,j-1)+ingrid(i-1,j)+ingrid(i,j))/4.0
292  end do
293  end do
294  ELSE IF(gridtype == 'C')THEN
295  call exch(ingrid(1,jsta_2l))
296  DO j=jsta_m,jend
297  DO i=1,im
298  outgrid(i,j)=(ingrid(i,j-1)+ingrid(i,j))/2.0
299  end do
300  end do
301  end if
302 
303  end subroutine v2h
304 !
305 !-------------------------------------------------------------------------------------
306 !
307  end module upp_math
Definition: MASKS_mod.f:1
dvdxdudy() computes dudy, dvdx, uwnd
Definition: UPP_MATH.f:17