Line data Source code
1 : ! This code is part of Radiative Transfer for Energetics (RTE)
2 : !
3 : ! Contacts: Robert Pincus and Eli Mlawer
4 : ! email: rrtmgp@aer.com
5 : !
6 : ! Copyright 2015- Atmospheric and Environmental Research,
7 : ! Regents of the University of Colorado,
8 : ! Trustees of Columbia University in the City of New York
9 : ! All right reserved.
10 : !
11 : ! Use and duplication is permitted under the terms of the
12 : ! BSD 3-clause license, see http://opensource.org/licenses/BSD-3-Clause
13 : ! -------------------------------------------------------------------------------------------------
14 : module mo_rte_util_array_validation
15 : !
16 : !> Provide utilites for sanitizing input arrays:
17 : ! checking values and sizes
18 : ! These are in a module so code can be written for both CPUs and GPUs
19 : ! Used only by Fortran classes so routines don't need C bindings and can use assumed-shape
20 : !
21 : use mo_rte_kind, only: wp, wl
22 : implicit none
23 : !>
24 : !> Values less than a floor (including masked versions)
25 : !>
26 : interface any_vals_less_than
27 : module procedure any_vals_less_than_1D, any_vals_less_than_2D, any_vals_less_than_3D
28 : module procedure any_vals_less_than_1D_masked, any_vals_less_than_2D_masked, any_vals_less_than_3D_masked
29 : end interface
30 : !>
31 : !> Values outside a range (including masked versions)
32 : !>
33 : interface any_vals_outside
34 : module procedure any_vals_outside_1D, any_vals_outside_2D, any_vals_outside_3D
35 : module procedure any_vals_outside_1D_masked, any_vals_outside_2D_masked, any_vals_outside_3D_masked
36 : end interface
37 : !>
38 : !> Find the extents of an array
39 : !>
40 : interface extents_are
41 : module procedure extents_are_1D, extents_are_2D, extents_are_3D
42 : module procedure extents_are_4D, extents_are_5D, extents_are_6D
43 : module procedure extents_are_2d_int
44 : end interface extents_are
45 :
46 : private
47 : public :: any_vals_less_than, any_vals_outside, extents_are
48 : contains
49 : !-------------------------------------------------------------------------------------------------
50 : ! Values less than a floor
51 : !-------------------------------------------------------------------------------------------------
52 0 : logical function any_vals_less_than_1D(array, check_value)
53 : real(wp), dimension(:), intent(in) :: array
54 : real(wp), intent(in) :: check_value
55 :
56 : real(wp) :: minValue
57 :
58 : !$acc kernels copyin(array)
59 : !$omp target map(to:array) map(from:minValue)
60 0 : minValue = minval(array)
61 : !$acc end kernels
62 : !$omp end target
63 :
64 0 : any_vals_less_than_1D = (minValue < check_value)
65 :
66 0 : end function any_vals_less_than_1D
67 : !-------------------------------------------------------------------------------------------------
68 6069330 : logical function any_vals_less_than_2D(array, check_value)
69 : real(wp), dimension(:,:), intent(in) :: array
70 : real(wp), intent(in) :: check_value
71 :
72 : real(wp) :: minValue
73 :
74 : !$acc kernels copyin(array)
75 : !$omp target map(to:array) map(from:minValue)
76 3374427337 : minValue = minval(array)
77 : !$acc end kernels
78 : !$omp end target
79 :
80 6069330 : any_vals_less_than_2D = (minValue < check_value)
81 :
82 6069330 : end function any_vals_less_than_2D
83 : !-------------------------------------------------------------------------------------------------
84 4152333 : logical function any_vals_less_than_3D(array, check_value)
85 : real(wp), dimension(:,:,:), intent(in) :: array
86 : real(wp), intent(in) :: check_value
87 :
88 : real(wp) :: minValue
89 :
90 : #ifdef _OPENMP
91 : integer :: dim1, dim2, dim3, i, j, k
92 : dim1 = size(array,1)
93 : dim2 = size(array,2)
94 : dim3 = size(array,3)
95 : minValue = check_value + epsilon(check_value) ! initialize to some value
96 : !$omp target teams map(to:array) &
97 : !$omp defaultmap(tofrom:scalar) reduction(min:minValue)
98 : !$omp distribute parallel do simd reduction(min:minValue)
99 : do i = 1, dim1
100 : do j = 1, dim2
101 : do k = 1, dim3
102 : minValue = min(minValue,array(i,j,k))
103 : enddo
104 : enddo
105 : enddo
106 : !$omp end target teams
107 : #else
108 : !$acc kernels copyin(array)
109 >66614*10^7 : minValue = minval(array)
110 : !$acc end kernels
111 : #endif
112 :
113 4152333 : any_vals_less_than_3D = (minValue < check_value)
114 :
115 4152333 : end function any_vals_less_than_3D
116 : !-------------------------------------------------------------------------------------------------
117 : ! Masked versions
118 : !-------------------------------------------------------------------------------------------------
119 0 : logical function any_vals_less_than_1D_masked(array, mask, check_value)
120 : real (wp), dimension(:), intent(in) :: array
121 : logical(wl), dimension(:), intent(in) :: mask
122 : real (wp), intent(in) :: check_value
123 :
124 : real(wp) :: minValue
125 :
126 : !$acc kernels copyin(array, mask)
127 : !$omp target map(to: array, mask) map(from:minValue)
128 0 : minValue = minval(array, mask=mask)
129 : !$acc end kernels
130 : !$omp end target
131 :
132 0 : any_vals_less_than_1D_masked = (minValue < check_value)
133 :
134 0 : end function any_vals_less_than_1D_masked
135 : !-------------------------------------------------------------------------------------------------
136 0 : logical function any_vals_less_than_2D_masked(array, mask, check_value)
137 : real (wp), dimension(:,:), intent(in) :: array
138 : logical(wl), dimension(:,:), intent(in) :: mask
139 : real (wp), intent(in) :: check_value
140 :
141 : real(wp) :: minValue
142 :
143 : !$acc kernels copyin(array, mask)
144 : !$omp target map(to: array, mask) map(from:minValue)
145 0 : minValue = minval(array, mask=mask)
146 : !$acc end kernels
147 : !$omp end target
148 :
149 0 : any_vals_less_than_2D_masked = (minValue < check_value)
150 :
151 0 : end function any_vals_less_than_2D_masked
152 : !-------------------------------------------------------------------------------------------------
153 0 : logical function any_vals_less_than_3D_masked(array, mask, check_value)
154 : real (wp), dimension(:,:,:), intent(in) :: array
155 : logical(wl), dimension(:,:,:), intent(in) :: mask
156 : real (wp), intent(in) :: check_value
157 :
158 : real(wp) :: minValue
159 :
160 : !$acc kernels copyin(array, mask)
161 : !$omp target map(to: array, mask) map(from:minValue)
162 0 : minValue = minval(array, mask=mask)
163 : !$acc end kernels
164 : !$omp end target
165 :
166 0 : any_vals_less_than_3D_masked = (minValue < check_value)
167 :
168 0 : end function any_vals_less_than_3D_masked
169 : !-------------------------------------------------------------------------------------------------
170 : ! Values outside a range
171 : !-------------------------------------------------------------------------------------------------
172 746136 : logical function any_vals_outside_1D(array, checkMin, checkMax)
173 : real(wp), dimension(:), intent(in) :: array
174 : real(wp), intent(in) :: checkMin, checkMax
175 :
176 : real(wp) :: minValue, maxValue
177 :
178 : !$acc kernels copyin(array)
179 : !$omp target map(to:array) map(from:minValue, maxValue)
180 13204872 : minValue = minval(array)
181 13204872 : maxValue = maxval(array)
182 : !$acc end kernels
183 : !$omp end target
184 746136 : any_vals_outside_1D = minValue < checkMin .or. maxValue > checkMax
185 :
186 746136 : end function any_vals_outside_1D
187 : ! ----------------------------------------------------------
188 15181840 : logical function any_vals_outside_2D(array, checkMin, checkMax)
189 : real(wp), dimension(:,:), intent(in) :: array
190 : real(wp), intent(in) :: checkMin, checkMax
191 :
192 : real(wp) :: minValue, maxValue
193 :
194 : !$acc kernels copyin(array)
195 : !$omp target map(to:array) map(from:minValue, maxValue)
196 19520995144 : minValue = minval(array)
197 19520995144 : maxValue = maxval(array)
198 : !$acc end kernels
199 : !$omp end target
200 15181840 : any_vals_outside_2D = minValue < checkMin .or. maxValue > checkMax
201 :
202 15181840 : end function any_vals_outside_2D
203 : ! ----------------------------------------------------------
204 2335578 : logical function any_vals_outside_3D(array, checkMin, checkMax)
205 : real(wp), dimension(:,:,:), intent(in) :: array
206 : real(wp), intent(in) :: checkMin, checkMax
207 :
208 : ! Compact version using intrinsics below
209 : ! but an explicit loop is the only current solution on GPUs
210 : real(wp) :: minValue, maxValue
211 :
212 :
213 : #ifdef _OPENMP
214 : integer :: dim1, dim2, dim3, i, j, k
215 : dim1 = size(array,1)
216 : dim2 = size(array,2)
217 : dim3 = size(array,3)
218 : minValue = checkMin + epsilon(checkMin) ! initialize to some value
219 : maxValue = checkMax - epsilon(checkMax) ! initialize to some value
220 : !$omp target teams map(to:array) &
221 : !$omp defaultmap(tofrom:scalar) reduction(min:minValue) reduction(max:maxValue)
222 : !$omp distribute parallel do simd reduction(min:minValue) reduction(max:maxValue)
223 : do i= 1, dim1
224 : do j = 1, dim2
225 : do k = 1, dim3
226 : minValue = min(minValue,array(i,j,k))
227 : maxValue = max(maxValue,array(i,j,k))
228 : enddo
229 : enddo
230 : enddo
231 : !$omp end target teams
232 : #else
233 : !$acc kernels copyin(array)
234 >39478*10^7 : minValue = minval(array)
235 >39478*10^7 : maxValue = maxval(array)
236 : !$acc end kernels
237 : #endif
238 :
239 2335578 : any_vals_outside_3D = minValue < checkMin .or. maxValue > checkMax
240 :
241 2335578 : end function any_vals_outside_3D
242 : ! ----------------------------------------------------------
243 : ! Masked versions
244 : ! ----------------------------------------------------------
245 0 : logical function any_vals_outside_1D_masked(array, mask, checkMin, checkMax)
246 : real (wp), dimension(:), intent(in) :: array
247 : logical(wl), dimension(:), intent(in) :: mask
248 : real(wp), intent(in) :: checkMin, checkMax
249 :
250 : real(wp) :: minValue, maxValue
251 :
252 : !$acc kernels copyin(array, mask)
253 : !$omp target map(to: array, mask) map(from:minValue, maxValue)
254 0 : minValue = minval(array, mask=mask)
255 0 : maxValue = maxval(array, mask=mask)
256 : !$acc end kernels
257 : !$omp end target
258 0 : any_vals_outside_1D_masked = minValue < checkMin .or. maxValue > checkMax
259 :
260 0 : end function any_vals_outside_1D_masked
261 : ! ----------------------------------------------------------
262 0 : logical function any_vals_outside_2D_masked(array, mask, checkMin, checkMax)
263 : real (wp), dimension(:,:), intent(in) :: array
264 : logical(wl), dimension(:,:), intent(in) :: mask
265 : real(wp), intent(in) :: checkMin, checkMax
266 :
267 : real(wp) :: minValue, maxValue
268 :
269 : !$acc kernels copyin(array, mask)
270 : !$omp target map(to: array, mask) map(from:minValue, maxValue)
271 0 : minValue = minval(array, mask=mask)
272 0 : maxValue = maxval(array, mask=mask)
273 : !$acc end kernels
274 : !$omp end target
275 0 : any_vals_outside_2D_masked = minValue < checkMin .or. maxValue > checkMax
276 :
277 0 : end function any_vals_outside_2D_masked
278 : ! ----------------------------------------------------------
279 0 : logical function any_vals_outside_3D_masked(array, mask, checkMin, checkMax)
280 : real (wp), dimension(:,:,:), intent(in) :: array
281 : logical(wl), dimension(:,:,:), intent(in) :: mask
282 : real(wp), intent(in) :: checkMin, checkMax
283 :
284 : real(wp) :: minValue, maxValue
285 :
286 : !$acc kernels copyin(array, mask)
287 : !$omp target map(to: array, mask) map(from:minValue, maxValue)
288 0 : minValue = minval(array, mask=mask)
289 0 : maxValue = maxval(array, mask=mask)
290 : !$acc end kernels
291 : !$omp end target
292 0 : any_vals_outside_3D_masked = minValue < checkMin .or. maxValue > checkMax
293 :
294 0 : end function any_vals_outside_3D_masked
295 : !--------------------------------------------------------------------------------------------------------------------
296 : !
297 : ! Extents
298 : !
299 : ! --------------------------------------------------------------------------------------
300 746136 : function extents_are_1d(array, n1)
301 : real(wp), dimension(: ), intent(in) :: array
302 : integer, intent(in) :: n1
303 : logical(wl) :: extents_are_1d
304 :
305 746136 : extents_are_1d = (size(array,1) == n1)
306 746136 : end function extents_are_1d
307 : ! --------------------------------------------------------------------------------------
308 12197296 : function extents_are_2d(array, n1, n2)
309 : real(wp), dimension(:,: ), intent(in) :: array
310 : integer, intent(in) :: n1, n2
311 : logical(wl) :: extents_are_2d
312 :
313 : extents_are_2d = (size(array,1) == n1 .and. &
314 12197296 : size(array,2) == n2 )
315 12197296 : end function extents_are_2d
316 : ! --------------------------------------------------------------------------------------
317 5417292 : function extents_are_3d(array, n1, n2, n3)
318 : real(wp), dimension(:,:,: ), intent(in) :: array
319 : integer, intent(in) :: n1, n2, n3
320 : logical(wl) :: extents_are_3d
321 :
322 : extents_are_3d = (size(array,1) == n1 .and. &
323 : size(array,2) == n2 .and. &
324 5417292 : size(array,3) == n3)
325 5417292 : end function extents_are_3d
326 : ! --------------------------------------------------------------------------------------
327 0 : function extents_are_4d(array, n1, n2, n3, n4)
328 : real(wp), dimension(:,:,:,: ), intent(in) :: array
329 : integer, intent(in) :: n1, n2, n3, n4
330 : logical(wl) :: extents_are_4d
331 :
332 : extents_are_4d = (size(array,1) == n1 .and. &
333 : size(array,2) == n2 .and. &
334 : size(array,3) == n3 .and. &
335 0 : size(array,4) == n4)
336 0 : end function extents_are_4d
337 : ! --------------------------------------------------------------------------------------
338 0 : function extents_are_5d(array, n1, n2, n3, n4, n5)
339 : real(wp), dimension(:,:,:,:,: ), intent(in) :: array
340 : integer, intent(in) :: n1, n2, n3, n4, n5
341 : logical(wl) :: extents_are_5d
342 :
343 : extents_are_5d = (size(array,1) == n1 .and. &
344 : size(array,2) == n2 .and. &
345 : size(array,3) == n3 .and. &
346 : size(array,4) == n4 .and. &
347 0 : size(array,5) == n5 )
348 0 : end function extents_are_5d
349 : ! --------------------------------------------------------------------------------------
350 0 : function extents_are_6d(array, n1, n2, n3, n4, n5, n6)
351 : real(wp), dimension(:,:,:,:,:,:), intent(in) :: array
352 : integer, intent(in) :: n1, n2, n3, n4, n5, n6
353 : logical(wl) :: extents_are_6d
354 :
355 : extents_are_6d = (size(array,1) == n1 .and. &
356 : size(array,2) == n2 .and. &
357 : size(array,3) == n3 .and. &
358 : size(array,4) == n4 .and. &
359 : size(array,5) == n5 .and. &
360 0 : size(array,6) == n6 )
361 0 : end function extents_are_6d
362 : ! --------------------------------------------------------------------------------------
363 3020006 : function extents_are_2d_int(array, n1, n2)
364 : integer, dimension(:,: ), intent(in) :: array
365 : integer, intent(in) :: n1, n2
366 : logical(wl) :: extents_are_2d_int
367 :
368 : extents_are_2d_int = (size(array,1) == n1 .and. &
369 3020006 : size(array,2) == n2 )
370 3020006 : end function extents_are_2d_int
371 :
372 : end module mo_rte_util_array_validation
|