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  var, &
82  valmin, &
83  valmax, &
84  varname, &
85  current_file, &
86  current_line )
87  use scale_prc, only: &
88  prc_abort, &
90  implicit none
91 
92  real(RP), intent(in) :: var(:)
93  real(RP), intent(in) :: valmin
94  real(RP), intent(in) :: valmax
95  character(len=*), intent(in) :: varname
96  character(len=*), intent(in) :: current_file
97  integer, intent(in) :: current_line
98 
99  logical :: invalid_value
100  integer :: k, kstr, kend
101  !---------------------------------------------------------------------------
102 
103  call prof_rapstart('Debug', 1)
104 
105  kstr = lbound( var(:), 1 )
106  kend = ubound( var(:), 1 )
107 
108  invalid_value = .false.
109  do k = kstr, kend
110  if ( var(k)*0.0_rp /= 0.0_rp &
111  .OR. var(k) < valmin &
112  .OR. var(k) > valmax ) then
113  invalid_value = .true.
114  exit
115  endif
116  enddo
117 
118  if ( invalid_value ) then
119  log_error("VALCHECK_1D",*) 'invalid value:', trim(varname), &
120  '(', prc_myrank, ',', k, ')=', var(k)
121  log_error_cont(*) 'in file : ', trim(current_file), ', at line : ', current_line
122  log_error_cont(*) 'in domain : ', debug_domain_num
123  call prc_abort
124  endif
125 
126  call prof_rapend ('Debug', 1)
127 
128  return
129  end subroutine valcheck_1d
130 
131  !-----------------------------------------------------------------------------
133  subroutine valcheck_2d( &
134  var, &
135  valmin, &
136  valmax, &
137  varname, &
138  current_file, &
139  current_line )
140  use scale_prc, only: &
141  prc_abort, &
142  prc_myrank
143  implicit none
144 
145  real(RP), intent(in) :: var(:,:)
146  real(RP), intent(in) :: valmin
147  real(RP), intent(in) :: valmax
148  character(len=*), intent(in) :: varname
149  character(len=*), intent(in) :: current_file
150  integer, intent(in) :: current_line
151 
152  logical :: invalid_value
153  integer :: k, kstr, kend
154  integer :: i, istr, iend
155  !---------------------------------------------------------------------------
156 
157  call prof_rapstart('Debug', 1)
158 
159  kstr = lbound( var(:,:), 1 )
160  kend = ubound( var(:,:), 1 )
161 
162  istr = lbound( var(:,:), 2 )
163  iend = ubound( var(:,:), 2 )
164 
165  invalid_value = .false.
166  outer:do i = istr, iend
167  do k = kstr, kend
168  if ( var(k,i)*0.0_rp /= 0.0_rp &
169  .OR. var(k,i) < valmin &
170  .OR. var(k,i) > valmax ) then
171  invalid_value = .true.
172  exit outer
173  endif
174  enddo
175  enddo outer
176 
177  if ( invalid_value ) then
178  log_error("VALCHECK_2D",*) 'invalid value:', trim(varname), &
179  '(', prc_myrank, ',', k, ',', i, ')=', var(k,i)
180  log_error_cont(*) 'in file : ', trim(current_file), ', at line : ', current_line
181  log_error_cont(*) 'in domain : ', debug_domain_num
182  call prc_abort
183  endif
184 
185  call prof_rapend ('Debug', 1)
186 
187  return
188  end subroutine valcheck_2d
189 
190  !-----------------------------------------------------------------------------
192  subroutine valcheck_3d( &
193  var, &
194  valmin, &
195  valmax, &
196  varname, &
197  current_file, &
198  current_line )
199  use scale_prc, only: &
200  prc_abort, &
201  prc_myrank
202  implicit none
203 
204  real(RP), intent(in) :: var(:,:,:)
205  real(RP), intent(in) :: valmin
206  real(RP), intent(in) :: valmax
207  character(len=*), intent(in) :: varname
208  character(len=*), intent(in) :: current_file
209  integer, intent(in) :: current_line
210 
211  logical :: invalid_value
212  integer :: k, kstr, kend
213  integer :: i, istr, iend
214  integer :: j, jstr, jend
215  !---------------------------------------------------------------------------
216 
217  call prof_rapstart('Debug', 1)
218 
219  kstr = lbound( var(:,:,:), 1 )
220  kend = ubound( var(:,:,:), 1 )
221 
222  istr = lbound( var(:,:,:), 2 )
223  iend = ubound( var(:,:,:), 2 )
224 
225  jstr = lbound( var(:,:,:), 3 )
226  jend = ubound( var(:,:,:), 3 )
227 
228  invalid_value = .false.
229  outer:do j = jstr, jend
230  do i = istr, iend
231  do k = kstr, kend
232  if ( var(k,i,j)*0.0_rp /= 0.0_rp &
233  .OR. var(k,i,j) < valmin &
234  .OR. var(k,i,j) > valmax ) then
235  invalid_value = .true.
236  exit outer
237  endif
238  enddo
239  enddo
240  enddo outer
241 
242  if ( invalid_value ) then
243  log_error("VALCHECK_3D",*) 'Invalid value:', trim(varname), &
244  '(', prc_myrank, ',', k, ',', i, ',', j, ')=', var(k,i,j)
245  log_error_cont(*) 'in file : ', trim(current_file), ', at line : ', current_line
246  log_error_cont(*) 'in domain : ', debug_domain_num
247  call prc_abort
248  endif
249 
250  call prof_rapend ('Debug', 1)
251 
252  return
253  end subroutine valcheck_3d
254 
255 end module scale_debug
module DEBUG
Definition: scale_debug.F90:11
subroutine valcheck_1d(var, valmin, valmax, varname, current_file, current_line)
Nan & extreme value checker (1D)
Definition: scale_debug.F90:87
subroutine valcheck_2d(var, valmin, valmax, varname, current_file, current_line)
Nan & extreme value checker (2D)
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:56
real(rp), public const_undef
Definition: scale_const.F90:41
module PROCESS
Definition: scale_prc.F90:11
subroutine valcheck_3d(var, valmin, valmax, varname, current_file, current_line)
Nan & extreme value checker (3D)
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:89
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:157
module profiler
Definition: scale_prof.F90:11
module PRECISION
module STDIO
Definition: scale_io.F90:10
integer, public debug_domain_num
Definition: scale_debug.F90:38
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:210