SCALE-RM
scale_debug.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_io
18  use scale_prof
19 #ifdef _OPENACC
20  use openacc
21 #endif
22  !-----------------------------------------------------------------------------
23  implicit none
24  private
25  !-----------------------------------------------------------------------------
26  !
27  !++ Public procedure
28  !
29  public :: check
30  public :: valcheck
31 
32  interface valcheck
33  module procedure valcheck_1d
34  module procedure valcheck_2d
35  module procedure valcheck_3d
36  end interface valcheck
37 
38  !-----------------------------------------------------------------------------
39  !
40  !++ Public parameters & variables
41  integer, public :: debug_domain_num = 0
42 
43  !
44  !-----------------------------------------------------------------------------
45  !
46  !++ Private procedure
47  !
48  !-----------------------------------------------------------------------------
49  !
50  !++ Private parameters & variables
51  !
52  !-----------------------------------------------------------------------------
53 contains
54  !-----------------------------------------------------------------------------
56  subroutine check( &
57  current_line, &
58  v )
59  use scale_prc, only: &
60  prc_abort
61  use scale_const, only: &
63  implicit none
64 
65  integer, intent(in) :: current_line
66  real(rp), intent(in) :: v
67  !---------------------------------------------------------------------------
68 
69  call prof_rapstart('Debug', 1)
70 
71  if ( .NOT. ( abs(v) < abs(const_undef) ) ) then
72  log_error("CHECK",*) 'uninitialized value at line:', current_line
73  call abort
74  end if
75 
76  call prof_rapend ('Debug', 1)
77 
78  return
79  end subroutine check
80 
81  !-----------------------------------------------------------------------------
83  subroutine valcheck_1d( &
84  KA, KS, KE, &
85  var, &
86  valmin, &
87  valmax, &
88  varname, &
89  current_file, &
90  current_line, &
91  mask )
92  use scale_prc, only: &
93  prc_abort, &
95  implicit none
96  integer , intent(in) :: KA, KS, KE
97  real(RP), intent(in) :: var(KA)
98  real(RP), intent(in) :: valmin
99  real(RP), intent(in) :: valmax
100  character(len=*), intent(in) :: varname
101  character(len=*), intent(in) :: current_file
102  integer, intent(in) :: current_line
103 
104  logical, intent(in), optional :: mask(KA)
105 
106  logical :: invalid_value
107  integer :: k
108  !---------------------------------------------------------------------------
109 
110  call prof_rapstart('Debug', 1)
111 
112  !$acc update host(var) if(acc_is_present(var))
113 
114  invalid_value = .false.
115  if ( present(mask) ) then
116  !$acc update host(mask) if(acc_is_present(mask))
117  do k = ks, ke
118  if ( .not. mask(k) ) cycle
119  if ( var(k)*0.0_rp /= 0.0_rp &
120  .OR. var(k) < valmin &
121  .OR. var(k) > valmax ) then
122  invalid_value = .true.
123  exit
124  endif
125  enddo
126  else
127  do k = ks, ke
128  if ( var(k)*0.0_rp /= 0.0_rp &
129  .OR. var(k) < valmin &
130  .OR. var(k) > valmax ) then
131  invalid_value = .true.
132  exit
133  endif
134  enddo
135  end if
136 
137  if ( invalid_value ) then
138  log_error("VALCHECK_1D",*) 'invalid value): ', trim(varname), &
139  '(', k, ')=', var(k)
140  log_error_cont(*) 'in file : ', trim(current_file), ', at line : ', current_line
141  call prc_abort
142  endif
143 
144  call prof_rapend ('Debug', 1)
145 
146  return
147  end subroutine valcheck_1d
148 
149  !-----------------------------------------------------------------------------
151  subroutine valcheck_2d( &
152  IA, IS, IE, JA, JS, JE, &
153  var, &
154  valmin, &
155  valmax, &
156  varname, &
157  current_file, &
158  current_line, &
159  mask )
160  use scale_prc, only: &
161  prc_abort, &
162  prc_myrank
163  implicit none
164  integer , intent(in) :: IA, IS, IE
165  integer , intent(in) :: JA, JS, JE
166  real(RP), intent(in) :: var(IA,JA)
167  real(RP), intent(in) :: valmin
168  real(RP), intent(in) :: valmax
169  character(len=*), intent(in) :: varname
170  character(len=*), intent(in) :: current_file
171  integer, intent(in) :: current_line
172 
173  logical, intent(in), optional :: mask(IA,JA)
174 
175  logical :: invalid_value
176  integer :: i, j
177  !---------------------------------------------------------------------------
178 
179  call prof_rapstart('Debug', 1)
180 
181  !$acc update host(var) if(acc_is_present(var))
182 
183  invalid_value = .false.
184  if ( present(mask) ) then
185  !$acc update host(mask) if(acc_is_present(mask))
186  outer1:do j = js, je
187  do i = is, ie
188  if ( .not. mask(i,j) ) cycle
189  if ( var(i,j)*0.0_rp /= 0.0_rp &
190  .OR. var(i,j) < valmin &
191  .OR. var(i,j) > valmax ) then
192  invalid_value = .true.
193  exit outer1
194  endif
195  enddo
196  enddo outer1
197  else
198  outer2:do j = js, je
199  do i = is, ie
200  if ( var(i,j)*0.0_rp /= 0.0_rp &
201  .OR. var(i,j) < valmin &
202  .OR. var(i,j) > valmax ) then
203  invalid_value = .true.
204  exit outer2
205  endif
206  enddo
207  enddo outer2
208  end if
209 
210  if ( invalid_value ) then
211  log_error("VALCHECK_2D",*) 'invalid value:', trim(varname), &
212  '(', i, ',', j, ')=', var(i,j)
213  log_error_cont(*) 'in file : ', trim(current_file), ', at line : ', current_line
214  call prc_abort
215  endif
216 
217  call prof_rapend ('Debug', 1)
218 
219  return
220  end subroutine valcheck_2d
221 
222  !-----------------------------------------------------------------------------
224  subroutine valcheck_3d( &
225  KA, KS, KE, &
226  IA, IS, IE, &
227  JA, JS, JE, &
228  var, &
229  valmin, &
230  valmax, &
231  varname, &
232  current_file, &
233  current_line, &
234  mask )
235  use scale_prc, only: &
236  prc_abort, &
237  prc_myrank
238  implicit none
239  integer , intent(in) :: KA, KS, KE
240  integer , intent(in) :: IA, IS, IE
241  integer , intent(in) :: JA, JS, JE
242  real(RP), intent(in) :: var(KA,IA,JA)
243  real(RP), intent(in) :: valmin
244  real(RP), intent(in) :: valmax
245  character(len=*), intent(in) :: varname
246  character(len=*), intent(in) :: current_file
247  integer, intent(in) :: current_line
248 
249  logical, intent(in), optional :: mask(IA,JA)
250 
251  logical :: invalid_value
252  integer :: k, i, j
253  !---------------------------------------------------------------------------
254 
255  call prof_rapstart('Debug', 1)
256 
257  !$acc update host(var) if(acc_is_present(var))
258 
259  invalid_value = .false.
260  if ( present(mask) ) then
261  !$acc update host(mask) if(acc_is_present(mask))
262  outer1:do j = js, je
263  do i = is, ie
264  if ( .not. mask(i,j) ) cycle
265  do k = ks, ke
266  if ( var(k,i,j)*0.0_rp /= 0.0_rp &
267  .OR. var(k,i,j) < valmin &
268  .OR. var(k,i,j) > valmax ) then
269  invalid_value = .true.
270  exit outer1
271  endif
272  enddo
273  enddo
274  enddo outer1
275  else
276  outer2:do j = js, je
277  do i = is, ie
278  do k = ks, ke
279  if ( var(k,i,j)*0.0_rp /= 0.0_rp &
280  .OR. var(k,i,j) < valmin &
281  .OR. var(k,i,j) > valmax ) then
282  invalid_value = .true.
283  exit outer2
284  endif
285  enddo
286  enddo
287  enddo outer2
288  end if
289 
290  if ( invalid_value ) then
291  log_error("VALCHECK_3D",*) 'Invalid value:', trim(varname), &
292  '(', k, ',', i, ',', j, ')=', var(k,i,j)
293  log_error_cont(*) 'in file : ', trim(current_file), ', at line : ', current_line
294  call prc_abort
295  endif
296 
297  call prof_rapend ('Debug', 1)
298 
299  return
300  end subroutine valcheck_3d
301 
302 end module scale_debug
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_prof::prof_rapstart
subroutine, public prof_rapstart(rapname_base, level, disable_barrier)
Start raptime.
Definition: scale_prof.F90:174
scale_debug::valcheck_1d
subroutine valcheck_1d(KA, KS, KE, var, valmin, valmax, varname, current_file, current_line, mask)
Nan & extreme value checker (1D)
Definition: scale_debug.F90:92
scale_prc::prc_myrank
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:91
scale_debug::valcheck_2d
subroutine valcheck_2d(IA, IS, IE, JA, JS, JE, var, valmin, valmax, varname, current_file, current_line, mask)
Nan & extreme value checker (2D)
Definition: scale_debug.F90:160
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_io
module STDIO
Definition: scale_io.F90:10
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_debug::check
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:59
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_debug::valcheck_3d
subroutine valcheck_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, var, valmin, valmax, varname, current_file, current_line, mask)
Nan & extreme value checker (3D)
Definition: scale_debug.F90:235
scale_debug
module DEBUG
Definition: scale_debug.F90:11
scale_prof::prof_rapend
subroutine, public prof_rapend(rapname_base, level, disable_barrier)
Save raptime.
Definition: scale_prof.F90:246
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_debug::debug_domain_num
integer, public debug_domain_num
Definition: scale_debug.F90:41