Line data Source code
1 : module cmparray_mod
2 :
3 : use shr_kind_mod, only : r8 => shr_kind_r8
4 :
5 : implicit none
6 : private
7 : save
8 :
9 : public expdaynite, cmpdaynite
10 :
11 : interface CmpDayNite
12 : module procedure CmpDayNite_1d_R
13 : module procedure CmpDayNite_2d_R
14 : module procedure CmpDayNite_3d_R
15 : module procedure CmpDayNite_1d_R_Copy
16 : module procedure CmpDayNite_2d_R_Copy
17 : module procedure CmpDayNite_3d_R_Copy
18 : module procedure CmpDayNite_1d_I
19 : module procedure CmpDayNite_2d_I
20 : module procedure CmpDayNite_3d_I
21 : end interface ! CmpDayNite
22 :
23 : interface ExpDayNite
24 : module procedure ExpDayNite_1d_R
25 : module procedure ExpDayNite_2d_R
26 : module procedure ExpDayNite_3d_R
27 : module procedure ExpDayNite_1d_I
28 : module procedure ExpDayNite_2d_I
29 : module procedure ExpDayNite_3d_I
30 : end interface ! ExpDayNite
31 :
32 : interface cmparray
33 : module procedure cmparray_1d_R
34 : module procedure cmparray_2d_R
35 : module procedure cmparray_3d_R
36 : end interface ! cmparray
37 :
38 : interface chksum
39 : module procedure chksum_1d_R
40 : module procedure chksum_2d_R
41 : module procedure chksum_3d_R
42 : module procedure chksum_1d_I
43 : module procedure chksum_2d_I
44 : module procedure chksum_3d_I
45 : end interface ! chksum
46 :
47 : contains
48 :
49 0 : subroutine CmpDayNite_1d_R(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1)
50 : integer, intent(in) :: Nday, Nnite
51 : integer, intent(in) :: il1, iu1
52 : integer, intent(in), dimension(Nday) :: IdxDay
53 : integer, intent(in), dimension(Nnite) :: IdxNite
54 : real(r8), intent(inout), dimension(il1:iu1) :: Array
55 :
56 0 : call CmpDayNite_3d_R(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, 1, 1, 1, 1)
57 :
58 0 : return
59 : end subroutine CmpDayNite_1d_R
60 :
61 391783 : subroutine CmpDayNite_2d_R(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, il2, iu2)
62 : integer, intent(in) :: Nday, Nnite
63 : integer, intent(in) :: il1, iu1
64 : integer, intent(in) :: il2, iu2
65 : integer, intent(in), dimension(Nday) :: IdxDay
66 : integer, intent(in), dimension(Nnite) :: IdxNite
67 : real(r8), intent(inout), dimension(il1:iu1,il2:iu2) :: Array
68 :
69 391783 : call CmpDayNite_3d_R(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, il2, iu2, 1, 1)
70 :
71 391783 : return
72 : end subroutine CmpDayNite_2d_R
73 :
74 391783 : subroutine CmpDayNite_3d_R(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, il2,iu2, il3, iu3)
75 : integer, intent(in) :: Nday, Nnite
76 : integer, intent(in) :: il1, iu1
77 : integer, intent(in) :: il2, iu2
78 : integer, intent(in) :: il3, iu3
79 : integer, intent(in), dimension(Nday) :: IdxDay
80 : integer, intent(in), dimension(Nnite) :: IdxNite
81 : real(r8), intent(inout), dimension(il1:iu1,il2:iu2,il3:iu3) :: Array
82 :
83 783566 : real(r8), dimension(il1:iu1) :: tmp
84 : integer :: i, j, k
85 :
86 :
87 783566 : do k = il3, iu3
88 11361707 : do j = il2, iu2
89 :
90 17887230 : tmp(1:Nnite) = Array(IdxNite(1:Nnite),j,k)
91 338708682 : Array(il1:il1+Nday-1,j,k) = Array(IdxDay(1:Nday),j,k)
92 18279013 : Array(il1+Nday:il1+Nday+Nnite-1,j,k) = tmp(1:Nnite)
93 :
94 : end do
95 : end do
96 :
97 391783 : return
98 : end subroutine CmpDayNite_3d_R
99 :
100 2742481 : subroutine CmpDayNite_1d_R_Copy(InArray, OutArray, Nday, IdxDay, Nnite, IdxNite, il1, iu1)
101 : integer, intent(in) :: Nday, Nnite
102 : integer, intent(in) :: il1, iu1
103 : integer, intent(in), dimension(Nday) :: IdxDay
104 : integer, intent(in), dimension(Nnite) :: IdxNite
105 : real(r8), intent(in), dimension(il1:iu1) :: InArray
106 : real(r8), intent(out), dimension(il1:iu1) :: OutArray
107 :
108 2742481 : call CmpDayNite_3d_R_Copy(InArray, OutArray, Nday, IdxDay, Nnite, IdxNite, il1, iu1, 1, 1, 1, 1)
109 :
110 2742481 : return
111 : end subroutine CmpDayNite_1d_R_Copy
112 :
113 3526047 : subroutine CmpDayNite_2d_R_Copy(InArray, OutArray, Nday, IdxDay, Nnite, IdxNite, il1, iu1, il2, iu2)
114 : integer, intent(in) :: Nday, Nnite
115 : integer, intent(in) :: il1, iu1
116 : integer, intent(in) :: il2, iu2
117 : integer, intent(in), dimension(Nday) :: IdxDay
118 : integer, intent(in), dimension(Nnite) :: IdxNite
119 : real(r8), intent(in), dimension(il1:iu1,il2:iu2) :: InArray
120 : real(r8), intent(out), dimension(il1:iu1,il2:iu2) :: OutArray
121 :
122 3526047 : call CmpDayNite_3d_R_Copy(InArray, OutArray, Nday, IdxDay, Nnite, IdxNite, il1, iu1, il2, iu2, 1, 1)
123 :
124 3526047 : return
125 : end subroutine CmpDayNite_2d_R_Copy
126 :
127 7835660 : subroutine CmpDayNite_3d_R_Copy(InArray, OutArray, Nday, IdxDay, Nnite, IdxNite, il1, iu1, il2,iu2, il3, iu3)
128 : integer, intent(in) :: Nday, Nnite
129 : integer, intent(in) :: il1, iu1
130 : integer, intent(in) :: il2, iu2
131 : integer, intent(in) :: il3, iu3
132 : integer, intent(in), dimension(Nday) :: IdxDay
133 : integer, intent(in), dimension(Nnite) :: IdxNite
134 : real(r8), intent(in), dimension(il1:iu1,il2:iu2,il3:iu3) :: InArray
135 : real(r8), intent(out), dimension(il1:iu1,il2:iu2,il3:iu3) :: OutArray
136 :
137 : integer :: i, j, k
138 :
139 :
140 43879696 : do k = il3, iu3
141 942629898 : do j = il2, iu2
142 :
143 14388846602 : do i=il1,il1+Nday-1
144 14388846602 : OutArray(i,j,k) = InArray(IdxDay(i-il1+1),j,k)
145 : enddo
146 1555796096 : do i=il1+Nday,il1+Nday+Nnite-1
147 1519752060 : OutArray(i,j,k) = InArray(IdxNite(i-(il1+Nday)+1),j,k)
148 : enddo
149 :
150 :
151 : end do
152 : end do
153 :
154 7835660 : return
155 : end subroutine CmpDayNite_3d_R_Copy
156 :
157 391783 : subroutine CmpDayNite_1d_I(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1)
158 : integer, intent(in) :: Nday, Nnite
159 : integer, intent(in) :: il1, iu1
160 : integer, intent(in), dimension(Nday) :: IdxDay
161 : integer, intent(in), dimension(Nnite) :: IdxNite
162 : integer, intent(inout), dimension(il1:iu1) :: Array
163 :
164 391783 : call CmpDayNite_3d_I(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, 1, 1, 1, 1)
165 :
166 391783 : return
167 : end subroutine CmpDayNite_1d_I
168 :
169 0 : subroutine CmpDayNite_2d_I(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, il2, iu2)
170 : integer, intent(in) :: Nday, Nnite
171 : integer, intent(in) :: il1, iu1
172 : integer, intent(in) :: il2, iu2
173 : integer, intent(in), dimension(Nday) :: IdxDay
174 : integer, intent(in), dimension(Nnite) :: IdxNite
175 : integer, intent(inout), dimension(il1:iu1,il2:iu2) :: Array
176 :
177 0 : call CmpDayNite_3d_I(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, il2, iu2, 1, 1)
178 :
179 0 : return
180 : end subroutine CmpDayNite_2d_I
181 :
182 391783 : subroutine CmpDayNite_3d_I(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, il2,iu2, il3, iu3)
183 : integer, intent(in) :: Nday, Nnite
184 : integer, intent(in) :: il1, iu1
185 : integer, intent(in) :: il2, iu2
186 : integer, intent(in) :: il3, iu3
187 : integer, intent(in), dimension(Nday) :: IdxDay
188 : integer, intent(in), dimension(Nnite) :: IdxNite
189 : integer, intent(inout), dimension(il1:iu1,il2:iu2,il3:iu3) :: Array
190 :
191 783566 : integer, dimension(il1:iu1) :: tmp
192 : integer :: i, j, k
193 :
194 :
195 783566 : do k = il3, iu3
196 1175349 : do j = il2, iu2
197 :
198 662490 : tmp(1:Nnite) = Array(IdxNite(1:Nnite),j,k)
199 12544766 : Array(il1:il1+Nday-1,j,k) = Array(IdxDay(1:Nday),j,k)
200 1054273 : Array(il1+Nday:il1+Nday+Nnite-1,j,k) = tmp(1:Nnite)
201 :
202 : end do
203 : end do
204 :
205 391783 : return
206 : end subroutine CmpDayNite_3d_I
207 :
208 7052094 : subroutine ExpDayNite_1d_R(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1)
209 : integer, intent(in) :: Nday, Nnite
210 : integer, intent(in) :: il1, iu1
211 : integer, intent(in), dimension(Nday) :: IdxDay
212 : integer, intent(in), dimension(Nnite) :: IdxNite
213 : real(r8), intent(inout), dimension(il1:iu1) :: Array
214 :
215 7052094 : call ExpDayNite_3d_R(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, 1, 1, 1, 1)
216 :
217 7052094 : return
218 : end subroutine ExpDayNite_1d_R
219 :
220 2742481 : subroutine ExpDayNite_2d_R(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, il2, iu2)
221 : integer, intent(in) :: Nday, Nnite
222 : integer, intent(in) :: il1, iu1
223 : integer, intent(in) :: il2, iu2
224 : integer, intent(in), dimension(Nday) :: IdxDay
225 : integer, intent(in), dimension(Nnite) :: IdxNite
226 : real(r8), intent(inout), dimension(il1:iu1,il2:iu2) :: Array
227 :
228 2742481 : call ExpDayNite_3d_R(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, il2, iu2, 1, 1)
229 :
230 2742481 : return
231 : end subroutine ExpDayNite_2d_R
232 :
233 9794575 : subroutine ExpDayNite_3d_R(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, il2,iu2, il3, iu3)
234 : integer, intent(in) :: Nday, Nnite
235 : integer, intent(in) :: il1, iu1
236 : integer, intent(in) :: il2, iu2
237 : integer, intent(in) :: il3, iu3
238 : integer, intent(in), dimension(Nday) :: IdxDay
239 : integer, intent(in), dimension(Nnite) :: IdxNite
240 : real(r8), intent(inout), dimension(il1:iu1,il2:iu2,il3:iu3) :: Array
241 :
242 19589150 : real(r8), dimension(il1:iu1) :: tmp
243 : integer :: i, j, k
244 :
245 :
246 19589150 : do k = il3, iu3
247 99121099 : do j = il2, iu2
248 :
249 1273293749 : tmp(1:Nday) = Array(1:Nday,j,k)
250 268970940 : Array(IdxNite(1:Nnite),j,k) = Array(il1+Nday:il1+Nday+Nnite-1,j,k)
251 1283088324 : Array(IdxDay(1:Nday),j,k) = tmp(1:Nday)
252 :
253 : end do
254 : end do
255 :
256 9794575 : return
257 : end subroutine ExpDayNite_3d_R
258 :
259 391783 : subroutine ExpDayNite_1d_I(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1)
260 : integer, intent(in) :: Nday, Nnite
261 : integer, intent(in) :: il1, iu1
262 : integer, intent(in), dimension(Nday) :: IdxDay
263 : integer, intent(in), dimension(Nnite) :: IdxNite
264 : integer, intent(inout), dimension(il1:iu1) :: Array
265 :
266 391783 : call ExpDayNite_3d_I(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, 1, 1, 1, 1)
267 :
268 391783 : return
269 : end subroutine ExpDayNite_1d_I
270 :
271 0 : subroutine ExpDayNite_2d_I(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, il2, iu2)
272 : integer, intent(in) :: Nday, Nnite
273 : integer, intent(in) :: il1, iu1
274 : integer, intent(in) :: il2, iu2
275 : integer, intent(in), dimension(Nday) :: IdxDay
276 : integer, intent(in), dimension(Nnite) :: IdxNite
277 : integer, intent(inout), dimension(il1:iu1,il2:iu2) :: Array
278 :
279 0 : call ExpDayNite_3d_I(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, il2, iu2, 1, 1)
280 :
281 0 : return
282 : end subroutine ExpDayNite_2d_I
283 :
284 391783 : subroutine ExpDayNite_3d_I(Array, Nday, IdxDay, Nnite, IdxNite, il1, iu1, il2,iu2, il3, iu3)
285 : integer, intent(in) :: Nday, Nnite
286 : integer, intent(in) :: il1, iu1
287 : integer, intent(in) :: il2, iu2
288 : integer, intent(in) :: il3, iu3
289 : integer, intent(in), dimension(Nday) :: IdxDay
290 : integer, intent(in), dimension(Nnite) :: IdxNite
291 : integer, intent(inout), dimension(il1:iu1,il2:iu2,il3:iu3) :: Array
292 :
293 783566 : integer, dimension(il1:iu1) :: tmp
294 : integer :: i, j, k
295 :
296 :
297 783566 : do k = il3, iu3
298 1175349 : do j = il2, iu2
299 :
300 6272383 : tmp(1:Nday) = Array(1:Nday,j,k)
301 1324980 : Array(IdxNite(1:Nnite),j,k) = Array(il1+Nday:il1+Nday+Nnite-1,j,k)
302 6664166 : Array(IdxDay(1:Nday),j,k) = tmp(1:Nday)
303 :
304 : end do
305 : end do
306 :
307 391783 : return
308 : end subroutine ExpDayNite_3d_I
309 :
310 : !******************************************************************************!
311 : ! !
312 : ! DEBUG !
313 : ! !
314 : !******************************************************************************!
315 :
316 : subroutine cmparray_1d_R(name, Ref, New, id1, is1, ie1)
317 : character(*), intent(in) :: name
318 : integer, intent(in) :: id1, is1, ie1
319 : real(r8), intent(in), dimension(id1) :: Ref
320 : real(r8), intent(in), dimension(id1) :: New
321 :
322 : call cmparray_3d_R(name, Ref, New, id1, is1, ie1, 1, 1, 1, 1, 1, 1)
323 : end subroutine cmparray_1d_R
324 :
325 : subroutine cmparray_2d_R(name, Ref, New, id1, is1, ie1, id2, is2, ie2)
326 : character(*), intent(in) :: name
327 : integer, intent(in) :: id1, is1, ie1
328 : integer, intent(in) :: id2, is2, ie2
329 : real(r8), intent(in), dimension(id1, id2) :: Ref
330 : real(r8), intent(in), dimension(id1, id2) :: New
331 :
332 : call cmparray_3d_R(name, Ref, New, id1, is1, ie1, id2, is2, ie2, 1, 1, 1)
333 : end subroutine cmparray_2d_R
334 :
335 : subroutine cmparray_3d_R(name, Ref, New, id1, is1, ie1, id2, is2, ie2, id3, is3, ie3)
336 : character(*), intent(in) :: name
337 : integer, intent(in) :: id1, is1, ie1
338 : integer, intent(in) :: id2, is2, ie2
339 : integer, intent(in) :: id3, is3, ie3
340 : real(r8), intent(in), dimension(id1, id2, id3) :: Ref
341 : real(r8), intent(in), dimension(id1, id2, id3) :: New
342 :
343 : integer :: i, j, k
344 : integer :: nerr
345 : logical :: found
346 : real(r8):: rdiff
347 : real(r8), parameter :: rtol = 1.0e-13_r8
348 :
349 : nerr = 0
350 :
351 : do k = is3, ie3
352 : do j = is2, ie2
353 :
354 : found = .false.
355 : do i = is1, ie1
356 : rdiff = abs(New(i,j,k)-Ref(i,j,k))
357 : rdiff = rdiff / merge(abs(Ref(i,j,k)), 1.0_r8, Ref(i,j,k) /= 0.0_r8)
358 : if ( rdiff > rtol ) then
359 : found = .true.
360 : exit
361 : end if
362 : end do
363 :
364 : if ( found ) then
365 : do i = is1, ie1
366 : rdiff = abs(New(i,j,k)-Ref(i,j,k))
367 : rdiff = rdiff / merge(abs(Ref(i,j,k)), 1.0_r8, Ref(i,j,k) /= 0.0_r8)
368 : if ( rdiff > rtol ) then
369 : print 666, name, i, j, k, Ref(i, j, k), New(i, j, k), rdiff
370 : nerr = nerr + 1
371 : if ( nerr > 10 ) stop
372 : end if
373 : end do
374 : end if
375 :
376 : end do
377 : end do
378 :
379 : return
380 : 666 format('cmp3d: ', a10, 3(1x, i4), 3(1x, e20.14))
381 :
382 : end subroutine cmparray_3d_R
383 :
384 : subroutine chksum_1d_R(name, Ref, id1, is1, ie1)
385 : character(*), intent(in) :: name
386 : integer, intent(in) :: id1, is1, ie1
387 : real(r8), intent(in), dimension(id1) :: Ref
388 :
389 : call chksum_3d_R(name, Ref, id1, is1, ie1, 1, 1, 1, 1, 1, 1)
390 : end subroutine chksum_1d_R
391 :
392 : subroutine chksum_1d_I(name, Ref, id1, is1, ie1)
393 : character(*), intent(in) :: name
394 : integer, intent(in) :: id1, is1, ie1
395 : integer, intent(in), dimension(id1) :: Ref
396 :
397 : call chksum_3d_I(name, Ref, id1, is1, ie1, 1, 1, 1, 1, 1, 1)
398 : end subroutine chksum_1d_I
399 :
400 : subroutine chksum_2d_R(name, Ref, id1, is1, ie1, id2, is2, ie2)
401 : character(*), intent(in) :: name
402 : integer, intent(in) :: id1, is1, ie1
403 : integer, intent(in) :: id2, is2, ie2
404 : real(r8), intent(in), dimension(id1, id2) :: Ref
405 :
406 : call chksum_3d_R(name, Ref, id1, is1, ie1, id2, is2, ie2, 1, 1, 1)
407 : end subroutine chksum_2d_R
408 :
409 : subroutine chksum_2d_I(name, Ref, id1, is1, ie1, id2, is2, ie2)
410 : character(*), intent(in) :: name
411 : integer, intent(in) :: id1, is1, ie1
412 : integer, intent(in) :: id2, is2, ie2
413 : integer, intent(in), dimension(id1, id2) :: Ref
414 :
415 : call chksum_3d_I(name, Ref, id1, is1, ie1, id2, is2, ie2, 1, 1, 1)
416 : end subroutine chksum_2d_I
417 :
418 : subroutine chksum_3d_R(name, Ref, id1, is1, ie1, id2, is2, ie2, id3, is3, ie3)
419 : character(*), intent(in) :: name
420 : integer, intent(in) :: id1, is1, ie1
421 : integer, intent(in) :: id2, is2, ie2
422 : integer, intent(in) :: id3, is3, ie3
423 : !orig real(r8), intent(in), dimension(id1, id2, id3) :: Ref
424 : real(r8), intent(in), dimension(is1:ie1, is2:ie2, is3:ie3) :: Ref
425 :
426 : real(r8) :: chksum
427 : real(r8) :: rmin, rmax
428 : integer :: i, j, k
429 : integer :: imin, jmin, kmin
430 : integer :: imax, jmax, kmax
431 :
432 : imin = is1 ; jmin = is2 ; kmin = is3
433 : imax = is1 ; jmax = is2 ; kmax = is3
434 : rmin = Ref(is1, is2, is3) ; rmax = rmin
435 :
436 : chksum = 0.0_r8
437 :
438 : do k = is3, ie3
439 : do j = is2, ie2
440 : do i = is1, ie1
441 : chksum = chksum + abs(Ref(i,j,k))
442 : if ( Ref(i,j,k) < rmin ) then
443 : rmin = Ref(i,j,k)
444 : imin = i ; jmin = j ; kmin = k
445 : end if
446 : if ( Ref(i,j,k) > rmax ) then
447 : rmax = Ref(i,j,k)
448 : imax = i ; jmax = j ; kmax = k
449 : end if
450 : end do
451 : end do
452 : end do
453 :
454 : print 666, name, chksum, imin, jmin, kmin, imax, jmax, kmax
455 : 666 format('chksum: ', a8, 1x, e20.14, 6(1x, i4))
456 :
457 : end subroutine chksum_3d_R
458 :
459 : subroutine chksum_3d_I(name, Ref, id1, is1, ie1, id2, is2, ie2, id3, is3, ie3)
460 : character(*), intent(in) :: name
461 : integer, intent(in) :: id1, is1, ie1
462 : integer, intent(in) :: id2, is2, ie2
463 : integer, intent(in) :: id3, is3, ie3
464 : integer, intent(in), dimension(id1, id2, id3) :: Ref
465 :
466 : integer :: i, j, k
467 : integer :: chksum
468 : chksum = 0
469 :
470 : do k = is3, ie3
471 : do j = is2, ie2
472 : do i = is1, ie1
473 : chksum = chksum + abs(Ref(i,j,k))
474 : end do
475 : end do
476 : end do
477 :
478 : print 666, name, chksum
479 : 666 format('chksum: ', a8, 1x, i8)
480 :
481 : end subroutine chksum_3d_I
482 :
483 : end module cmparray_mod
|