UPP (develop)
Loading...
Searching...
No Matches
UPP_MATH.f
Go to the documentation of this file.
1
19 module upp_math
20
21 use masks, only: dx, dy
22 use ctlblk_mod, only: im, jsta_2l, jend_2u, jsta_m, jend_m, spval,&
23 ista_2l, iend_2u, ista_m, iend_m
24 use gridspec_mod, only: gridtype
25!
26 implicit none
27
28 private
29
30 public :: ddvdx, ddudy, uuavg
31 public :: dvdxdudy
32 public :: h2u, h2v, u2h, v2h
33
35 REAL, ALLOCATABLE :: ddvdx(:,:)
36
38 REAL, ALLOCATABLE :: ddudy(:,:)
39
41 REAL, ALLOCATABLE :: uuavg(:,:)
42!
43!-------------------------------------------------------------------------------------
44!
45 contains
46!
47!-------------------------------------------------------------------------------------
52!
53
54 subroutine dvdxdudy(uwnd,vwnd)
55!
56 implicit none
57
58 REAL, dimension(ista_2l:iend_2u,jsta_2l:jend_2u), intent(in) :: uwnd, vwnd
59!
60!-- local variables
61!--
62 integer i, j
63 real r2dx, r2dy
64 INTEGER, allocatable :: ihe(:),ihw(:)
65
66!Initializing
67 DO j=jsta_m,jend_m
68 DO i=ista_m,iend_m
69 ddvdx(i,j)=spval
70 ddudy(i,j)=spval
71 uuavg(i,j)=spval
72 ENDDO
73 ENDDO
74
75 IF(gridtype == 'A')THEN
76!$omp parallel do private(i,j,r2dx,r2dy)
77 DO j=jsta_m,jend_m
78 DO i=ista_m,iend_m
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)
88 END IF
89 END DO
90 END DO
91 ELSE IF (gridtype == 'E')THEN
92 allocate(ihw(jsta_2l:jend_2u), ihe(jsta_2l:jend_2u))
93!$omp parallel do private(j)
94 DO j=jsta_2l,jend_2u
95 ihw(j) = -mod(j,2)
96 ihe(j) = ihw(j)+1
97 ENDDO
98!$omp parallel do private(i,j,r2dx,r2dy)
99 DO j=jsta_m,jend_m
100 DO i=ista_m,iend_m
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))
109 END IF
110 END DO
111 END DO
112 deallocate(ihw, ihe)
113 ELSE IF (gridtype == 'B')THEN
114!$omp parallel do private(i,j,r2dx,r2dy)
115 DO j=jsta_m,jend_m
116 DO i=ista_m,iend_m
117 r2dx = 1./dx(i,j)
118 r2dy = 1./dy(i,j)
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))
129 END DO
130 END DO
131 END IF
132
133 end subroutine dvdxdudy
134!
135!-------------------------------------------------------------------------------------
136!
137 subroutine h2u(ingrid,outgrid)
138!
145!
146
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
150
151 implicit none
152
153 include "mpif.h"
154 integer:: i,j,ie,iw
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
157 outgrid=spval
158 if(gridtype == 'A')THEN
159 do j=jsta,jend
160 do i=ista,iend
161 outgrid(i,j)=ingrid(i,j)
162 end do
163 end do
164 else IF(gridtype == 'E')THEN
165 call exch(ingrid(ista_2l,jsta_2l))
166 DO j=jsta_m,jend_m
167 DO i=ista_m,iend_m
168 ie=i+mod(j,2)
169 iw=ie-1
170 outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
171 end do
172 end do
173 ELSE IF(gridtype == 'B')THEN
174 call exch(ingrid(ista_2l,jsta_2l))
175 DO j=jsta,jend_m
176 DO i=ista,iend_m
177 outgrid(i,j)=(ingrid(i,j)+ingrid(i,j+1)+ingrid(i+1,j)+ingrid(i+1,j+1))/4.0
178 end do
179 end do
180! Fill in boundary points because hysplit fails when 10 m wind has bitmaps
181 do j=jsta,jend_m
182 outgrid(iend,j)=outgrid(iend-1,j)
183 end do
184 IF(me == (num_procs-1) .and. jend_2u >= jm) then
185 DO i=ista,iend
186 outgrid(i,jend) = outgrid(i,jend-1)
187 END DO
188 END IF
189 ELSE IF(gridtype == 'C')THEN
190 call exch(ingrid(ista_2l,jsta_2l))
191 DO j=jsta,jend
192 DO i=ista,iend_m
193 outgrid(i,j)=(ingrid(i,j)+ingrid(i+1,j))/2.0
194 end do
195 end do
196 end if
197
198 end subroutine h2u
199!
200!-------------------------------------------------------------------------------------
201!
202 subroutine h2v(ingrid,outgrid)
203!
210!
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
214 implicit none
215 include "mpif.h"
216 integer:: i,j,ie,iw
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
219 outgrid=spval
220 if(gridtype == 'A')THEN
221 do j=jsta,jend
222 do i=ista,iend
223 outgrid(i,j)=ingrid(i,j)
224 end do
225 end do
226 else IF(gridtype == 'E')THEN
227 call exch(ingrid(ista_2l,jsta_2l))
228 DO j=jsta_m,jend_m
229 DO i=ista_m,iend_m
230 ie=i+mod(j,2)
231 iw=ie-1
232 outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
233 end do
234 end do
235 ELSE IF(gridtype == 'B')THEN
236 call exch(ingrid(ista_2l,jsta_2l))
237 DO j=jsta,jend_m
238 DO i=ista,iend_m
239 outgrid(i,j)=(ingrid(i,j)+ingrid(i,j+1)+ingrid(i+1,j)+ingrid(i+1,j+1))/4.0
240 end do
241 end do
242 ELSE IF(gridtype == 'C')THEN
243 call exch(ingrid(ista_2l,jsta_2l))
244 DO j=jsta,jend_m
245 DO i=ista,iend
246 outgrid(i,j)=(ingrid(i,j)+ingrid(i,j+1))/2.0
247 end do
248 end do
249 end if
250
251 end subroutine h2v
252!
253!-------------------------------------------------------------------------------------
254!
255 subroutine u2h(ingrid,outgrid)
256!
263!
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
267 implicit none
268 include "mpif.h"
269 integer:: i,j,ie,iw
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
272 outgrid=spval
273 if(gridtype == 'A')THEN
274 do j=jsta,jend
275 do i=ista,iend
276 outgrid(i,j)=ingrid(i,j)
277 end do
278 end do
279 else IF(gridtype == 'E')THEN
280 call exch(ingrid(ista_2l,jsta_2l))
281 DO j=jsta_m,jend_m
282 DO i=ista_m,iend_m
283 ie=i+mod(j+1,2)
284 iw=ie-1
285 outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
286 end do
287 end do
288 ELSE IF(gridtype == 'B')THEN
289 call exch(ingrid(ista_2l,jsta_2l))
290 DO j=jsta_m,jend_m
291 DO i=ista_m,iend_m
292 outgrid(i,j)=(ingrid(i-1,j-1)+ingrid(i,j-1)+ingrid(i-1,j)+ingrid(i,j))/4.0
293 end do
294 end do
295 ELSE IF(gridtype == 'C')THEN
296 call exch(ingrid(ista_2l,jsta_2l))
297 DO j=jsta,jend
298 DO i=ista_m,iend
299 outgrid(i,j)=(ingrid(i-1,j)+ingrid(i,j))/2.0
300 end do
301 end do
302 end if
303
304 end subroutine u2h
305!
306!-------------------------------------------------------------------------------------
307!
308 subroutine v2h(ingrid,outgrid)
309!
316!
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
320 implicit none
321 include "mpif.h"
322 integer:: i,j,ie,iw
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
325 outgrid=spval
326 if(gridtype == 'A')THEN
327 do j=jsta,jend
328 do i=ista,iend
329 outgrid(i,j)=ingrid(i,j)
330 end do
331 end do
332 else IF(gridtype == 'E')THEN
333 call exch(ingrid(ista_2l,jsta_2l))
334 DO j=jsta_m,jend_m
335 DO i=ista_m,iend_m
336 ie=i+mod(j,2)
337 iw=ie-1
338 outgrid(i,j)=(ingrid(iw,j)+ingrid(ie,j)+ingrid(i,j+1)+ingrid(i,j-1))/4.0
339 end do
340 end do
341 ELSE IF(gridtype == 'B')THEN
342 call exch(ingrid(ista_2l,jsta_2l))
343 DO j=jsta_m,jend_m
344 DO i=ista_m,iend_m
345 outgrid(i,j)=(ingrid(i-1,j-1)+ingrid(i,j-1)+ingrid(i-1,j)+ingrid(i,j))/4.0
346 end do
347 end do
348 ELSE IF(gridtype == 'C')THEN
349 call exch(ingrid(ista_2l,jsta_2l))
350 DO j=jsta_m,jend
351 DO i=ista,iend
352 outgrid(i,j)=(ingrid(i,j-1)+ingrid(i,j))/2.0
353 end do
354 end do
355 end if
356
357 end subroutine v2h
358!
359!-------------------------------------------------------------------------------------
360!
361 end module upp_math
subroutine exch(a)
exch() Subroutine that exchanges one halo row.
Definition EXCH.f:27