SCALE-RM
scale_debug.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14  !-----------------------------------------------------------------------------
15  !
16  !++ used modules
17  !
18  use scale_precision
19  use scale_stdio
20  use scale_prof
21  !-----------------------------------------------------------------------------
22  implicit none
23  private
24  !-----------------------------------------------------------------------------
25  !
26  !++ Public procedure
27  !
28  public :: check
29  public :: valcheck
30 
31  interface valcheck
32  module procedure valcheck_1d
33  module procedure valcheck_2d
34  module procedure valcheck_3d
35  end interface valcheck
36 
37  !-----------------------------------------------------------------------------
38  !
39  !++ Public parameters & variables
40  integer, public :: debug_domain_num = 0
41 
42  !
43  !-----------------------------------------------------------------------------
44  !
45  !++ Private procedure
46  !
47  !-----------------------------------------------------------------------------
48  !
49  !++ Private parameters & variables
50  !
51  !-----------------------------------------------------------------------------
52 contains
53  !-----------------------------------------------------------------------------
55  subroutine check( &
56  current_line, &
57  v )
58  use scale_process, only: &
60  use scale_const, only: &
62  implicit none
63 
64  integer, intent(in) :: current_line
65  real(RP), intent(in) :: v
66  !---------------------------------------------------------------------------
67 
68  call prof_rapstart('Debug', 1)
69 
70  if ( .NOT. ( abs(v) < abs(const_undef) ) ) then
71 ! if( IO_L ) write(IO_FID_LOG,,*) 'xxx uninitialized value at line:', current_line
72 
73  write(*,*) 'xxx uninitialized value at line:', current_line
74 
75  call abort
76  end if
77 
78  call prof_rapend ('Debug', 1)
79 
80  return
81  end subroutine check
82 
83  !-----------------------------------------------------------------------------
85  subroutine valcheck_1d( &
86  var, &
87  valmin, &
88  valmax, &
89  varname, &
90  current_file, &
91  current_line )
92  use scale_process, only: &
93  prc_mpistop, &
95  implicit none
96 
97  real(RP), intent(in) :: var(:)
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 :: invalid_value
105  integer :: k, kstr, kend
106  !---------------------------------------------------------------------------
107 
108  call prof_rapstart('Debug', 1)
109 
110  kstr = lbound( var(:), 1 )
111  kend = ubound( var(:), 1 )
112 
113  invalid_value = .false.
114  do k = kstr, kend
115  if ( var(k)*0.0_rp /= 0.0_rp &
116  .OR. var(k) < valmin &
117  .OR. var(k) > valmax ) then
118  invalid_value = .true.
119  exit
120  endif
121  enddo
122 
123  if ( invalid_value ) then
124 ! if( IO_L ) write(IO_FID_LOG,,*) 'xxx invalid value:', trim(varname), '(', k, ')=', var(k)
125 ! if( IO_L ) write(IO_FID_LOG,,*) 'xxx in file:', trim(current_file), ', at line:', current_line
126 
127  write(*,*) 'xxx invalid value:', trim(varname), &
128  '(', prc_myrank, ',', k, ')=', var(k)
129  write(*,*) 'xxx in file:', trim(current_file), ', at line:', current_line
130  write(*,*) 'xxx in domain:', debug_domain_num
131 
132  call prc_mpistop
133  endif
134 
135  call prof_rapend ('Debug', 1)
136 
137  return
138  end subroutine valcheck_1d
139 
140  !-----------------------------------------------------------------------------
142  subroutine valcheck_2d( &
143  var, &
144  valmin, &
145  valmax, &
146  varname, &
147  current_file, &
148  current_line )
149  use scale_process, only: &
150  prc_mpistop, &
151  prc_myrank
152  implicit none
153 
154  real(RP), intent(in) :: var(:,:)
155  real(RP), intent(in) :: valmin
156  real(RP), intent(in) :: valmax
157  character(len=*), intent(in) :: varname
158  character(len=*), intent(in) :: current_file
159  integer, intent(in) :: current_line
160 
161  logical :: invalid_value
162  integer :: k, kstr, kend
163  integer :: i, istr, iend
164  !---------------------------------------------------------------------------
165 
166  call prof_rapstart('Debug', 1)
167 
168  kstr = lbound( var(:,:), 1 )
169  kend = ubound( var(:,:), 1 )
170 
171  istr = lbound( var(:,:), 2 )
172  iend = ubound( var(:,:), 2 )
173 
174  invalid_value = .false.
175  outer:do i = istr, iend
176  do k = kstr, kend
177  if ( var(k,i)*0.0_rp /= 0.0_rp &
178  .OR. var(k,i) < valmin &
179  .OR. var(k,i) > valmax ) then
180  invalid_value = .true.
181  exit outer
182  endif
183  enddo
184  enddo outer
185 
186  if ( invalid_value ) then
187 ! if( IO_L ) write(IO_FID_LOG,,*) 'xxx invalid value:', trim(varname), '(', k, ',', i, ')=', var(k,i)
188 ! if( IO_L ) write(IO_FID_LOG,,*) 'xxx in file:', trim(current_file), ', at line:', current_line
189 
190  write(*,*) 'xxx invalid value:', trim(varname), &
191  '(', prc_myrank, ',', k, ',', i, ')=', var(k,i)
192  write(*,*) 'xxx in file:', trim(current_file), ', at line:', current_line
193  write(*,*) 'xxx in domain:', debug_domain_num
194 
195  call prc_mpistop
196  endif
197 
198  call prof_rapend ('Debug', 1)
199 
200  return
201  end subroutine valcheck_2d
202 
203  !-----------------------------------------------------------------------------
205  subroutine valcheck_3d( &
206  var, &
207  valmin, &
208  valmax, &
209  varname, &
210  current_file, &
211  current_line )
212  use scale_process, only: &
213  prc_mpistop, &
214  prc_myrank
215  implicit none
216 
217  real(RP), intent(in) :: var(:,:,:)
218  real(RP), intent(in) :: valmin
219  real(RP), intent(in) :: valmax
220  character(len=*), intent(in) :: varname
221  character(len=*), intent(in) :: current_file
222  integer, intent(in) :: current_line
223 
224  logical :: invalid_value
225  integer :: k, kstr, kend
226  integer :: i, istr, iend
227  integer :: j, jstr, jend
228  !---------------------------------------------------------------------------
229 
230  call prof_rapstart('Debug', 1)
231 
232  kstr = lbound( var(:,:,:), 1 )
233  kend = ubound( var(:,:,:), 1 )
234 
235  istr = lbound( var(:,:,:), 2 )
236  iend = ubound( var(:,:,:), 2 )
237 
238  jstr = lbound( var(:,:,:), 3 )
239  jend = ubound( var(:,:,:), 3 )
240 
241  invalid_value = .false.
242  outer:do j = jstr, jend
243  do i = istr, iend
244  do k = kstr, kend
245  if ( var(k,i,j)*0.0_rp /= 0.0_rp &
246  .OR. var(k,i,j) < valmin &
247  .OR. var(k,i,j) > valmax ) then
248  invalid_value = .true.
249  exit outer
250  endif
251  enddo
252  enddo
253  enddo outer
254 
255  if ( invalid_value ) then
256 ! if( IO_L ) write(IO_FID_LOG,,*) 'xxx invalid value:', trim(varname), '(', k, ',', i, ',', j, ')=', var(k,i,j)
257 ! if( IO_L ) write(IO_FID_LOG,,*) 'xxx in file:', trim(current_file), ', at line:', current_line
258 
259  write(*,*) 'xxx invalid value:', trim(varname), &
260  '(', prc_myrank, ',', k, ',', i, ',', j, ')=', var(k,i,j)
261  write(*,*) 'xxx in file:', trim(current_file), ', at line:', current_line
262  write(*,*) 'xxx in domain:', debug_domain_num
263 
264  call prc_mpistop
265  endif
266 
267  call prof_rapend ('Debug', 1)
268 
269  return
270  end subroutine valcheck_3d
271 
272 end module scale_debug
273 !-------------------------------------------------------------------------------
module DEBUG
Definition: scale_debug.F90:13
subroutine valcheck_1d(var, valmin, valmax, varname, current_file, current_line)
Nan & extreme value checker (1D)
Definition: scale_debug.F90:92
subroutine, public prc_mpistop
Abort MPI.
subroutine valcheck_2d(var, valmin, valmax, varname, current_file, current_line)
Nan & extreme value checker (2D)
module STDIO
Definition: scale_stdio.F90:12
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:58
real(rp), public const_undef
Definition: scale_const.F90:43
subroutine valcheck_3d(var, valmin, valmax, varname, current_file, current_line)
Nan & extreme value checker (3D)
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, public prc_myrank
process num in local communicator
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
module profiler
Definition: scale_prof.F90:10
module PRECISION
integer, public debug_domain_num
Definition: scale_debug.F90:40
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178