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