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  write(*,*) 'xxx uninitialized value at line:', current_line
72  call abort
73  end if
74 
75  call prof_rapend ('Debug', 1)
76 
77  return
78  end subroutine check
79 
80  !-----------------------------------------------------------------------------
82  subroutine valcheck_1d( &
83  var, &
84  valmin, &
85  valmax, &
86  varname, &
87  current_file, &
88  current_line )
89  use scale_process, only: &
90  prc_mpistop, &
92  implicit none
93 
94  real(RP), intent(in) :: var(:)
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 :: invalid_value
102  integer :: k, kstr, kend
103  !---------------------------------------------------------------------------
104 
105  call prof_rapstart('Debug', 1)
106 
107  kstr = lbound( var(:), 1 )
108  kend = ubound( var(:), 1 )
109 
110  invalid_value = .false.
111  do k = kstr, kend
112  if ( var(k)*0.0_rp /= 0.0_rp &
113  .OR. var(k) < valmin &
114  .OR. var(k) > valmax ) then
115  invalid_value = .true.
116  exit
117  endif
118  enddo
119 
120  if ( invalid_value ) then
121  write(*,*) 'xxx [VALCHECK_1D] invalid value:', trim(varname), &
122  '(', prc_myrank, ',', k, ')=', var(k)
123  write(*,*) 'xxx in file : ', trim(current_file), ', at line : ', current_line
124  write(*,*) 'xxx in domain : ', debug_domain_num
125  call prc_mpistop
126  endif
127 
128  call prof_rapend ('Debug', 1)
129 
130  return
131  end subroutine valcheck_1d
132 
133  !-----------------------------------------------------------------------------
135  subroutine valcheck_2d( &
136  var, &
137  valmin, &
138  valmax, &
139  varname, &
140  current_file, &
141  current_line )
142  use scale_process, only: &
143  prc_mpistop, &
144  prc_myrank
145  implicit none
146 
147  real(RP), intent(in) :: var(:,:)
148  real(RP), intent(in) :: valmin
149  real(RP), intent(in) :: valmax
150  character(len=*), intent(in) :: varname
151  character(len=*), intent(in) :: current_file
152  integer, intent(in) :: current_line
153 
154  logical :: invalid_value
155  integer :: k, kstr, kend
156  integer :: i, istr, iend
157  !---------------------------------------------------------------------------
158 
159  call prof_rapstart('Debug', 1)
160 
161  kstr = lbound( var(:,:), 1 )
162  kend = ubound( var(:,:), 1 )
163 
164  istr = lbound( var(:,:), 2 )
165  iend = ubound( var(:,:), 2 )
166 
167  invalid_value = .false.
168  outer:do i = istr, iend
169  do k = kstr, kend
170  if ( var(k,i)*0.0_rp /= 0.0_rp &
171  .OR. var(k,i) < valmin &
172  .OR. var(k,i) > valmax ) then
173  invalid_value = .true.
174  exit outer
175  endif
176  enddo
177  enddo outer
178 
179  if ( invalid_value ) then
180  write(*,*) 'xxx [VALCHECK_2D] invalid value:', trim(varname), &
181  '(', prc_myrank, ',', k, ',', i, ')=', var(k,i)
182  write(*,*) 'xxx in file : ', trim(current_file), ', at line : ', current_line
183  write(*,*) 'xxx in domain : ', debug_domain_num
184  call prc_mpistop
185  endif
186 
187  call prof_rapend ('Debug', 1)
188 
189  return
190  end subroutine valcheck_2d
191 
192  !-----------------------------------------------------------------------------
194  subroutine valcheck_3d( &
195  var, &
196  valmin, &
197  valmax, &
198  varname, &
199  current_file, &
200  current_line )
201  use scale_process, only: &
202  prc_mpistop, &
203  prc_myrank
204  implicit none
205 
206  real(RP), intent(in) :: var(:,:,:)
207  real(RP), intent(in) :: valmin
208  real(RP), intent(in) :: valmax
209  character(len=*), intent(in) :: varname
210  character(len=*), intent(in) :: current_file
211  integer, intent(in) :: current_line
212 
213  logical :: invalid_value
214  integer :: k, kstr, kend
215  integer :: i, istr, iend
216  integer :: j, jstr, jend
217  !---------------------------------------------------------------------------
218 
219  call prof_rapstart('Debug', 1)
220 
221  kstr = lbound( var(:,:,:), 1 )
222  kend = ubound( var(:,:,:), 1 )
223 
224  istr = lbound( var(:,:,:), 2 )
225  iend = ubound( var(:,:,:), 2 )
226 
227  jstr = lbound( var(:,:,:), 3 )
228  jend = ubound( var(:,:,:), 3 )
229 
230  invalid_value = .false.
231  outer:do j = jstr, jend
232  do i = istr, iend
233  do k = kstr, kend
234  if ( var(k,i,j)*0.0_rp /= 0.0_rp &
235  .OR. var(k,i,j) < valmin &
236  .OR. var(k,i,j) > valmax ) then
237  invalid_value = .true.
238  exit outer
239  endif
240  enddo
241  enddo
242  enddo outer
243 
244  if ( invalid_value ) then
245  write(*,*) 'xxx [VALCHECK_3D] Invalid value:', trim(varname), &
246  '(', prc_myrank, ',', k, ',', i, ',', j, ')=', var(k,i,j)
247  write(*,*) 'xxx in file : ', trim(current_file), ', at line : ', current_line
248  write(*,*) 'xxx in domain : ', debug_domain_num
249  call prc_mpistop
250  endif
251 
252  call prof_rapend ('Debug', 1)
253 
254  return
255  end subroutine valcheck_3d
256 
257 end module scale_debug
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:89
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:156
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:204