SCALE-RM
Data Types | Functions/Subroutines | Variables
scale_debug Module Reference

module DEBUG More...

Functions/Subroutines

subroutine, public check (current_line, v)
 Undefined value checker. More...
 
subroutine valcheck_1d (var, valmin, valmax, varname, current_file, current_line)
 Nan & extreme value checker (1D) More...
 
subroutine valcheck_2d (var, valmin, valmax, varname, current_file, current_line)
 Nan & extreme value checker (2D) More...
 
subroutine valcheck_3d (var, valmin, valmax, varname, current_file, current_line)
 Nan & extreme value checker (3D) More...
 

Variables

integer, public debug_domain_num = 0
 

Detailed Description

module DEBUG

Description
Debug & Value check tools
Author
Team SCALE
History
  • 2014-01-15 (H.Yashiro) Merge scale_debug & scale_misc(valcheck)

Function/Subroutine Documentation

◆ check()

subroutine, public scale_debug::check ( integer, intent(in)  current_line,
real(rp), intent(in)  v 
)

Undefined value checker.

Definition at line 58 of file scale_debug.F90.

References scale_const::const_undef, scale_process::prc_mpistop(), scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by scale_atmos_dyn_common::atmos_dyn_fct(), scale_atmos_dyn_fvm_flux_cd2::atmos_dyn_fvm_flux_valuew_z_cd2(), scale_atmos_dyn_fvm_flux_cd4::atmos_dyn_fvm_flux_valuew_z_cd4(), scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_flux_valuew_z_cd6(), scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_flux_valuew_z_ud1(), scale_atmos_dyn_fvm_flux_ud3::atmos_dyn_fvm_flux_valuew_z_ud3(), scale_atmos_dyn_fvm_flux_ud3koren1993::atmos_dyn_fvm_flux_valuew_z_ud3koren1993(), scale_atmos_dyn_fvm_flux_ud5::atmos_dyn_fvm_flux_valuew_z_ud5(), scale_atmos_dyn_fvm_flux_cd2::atmos_dyn_fvm_fluxx_uyz_cd2(), scale_atmos_dyn_fvm_flux_cd4::atmos_dyn_fvm_fluxx_uyz_cd4(), scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxx_uyz_cd6(), scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxx_uyz_ud1(), scale_atmos_dyn_fvm_flux_ud3::atmos_dyn_fvm_fluxx_uyz_ud3(), scale_atmos_dyn_fvm_flux_ud3koren1993::atmos_dyn_fvm_fluxx_uyz_ud3koren1993(), scale_atmos_dyn_fvm_flux_ud5::atmos_dyn_fvm_fluxx_uyz_ud5(), scale_atmos_dyn_fvm_flux_cd2::atmos_dyn_fvm_fluxx_xvz_cd2(), scale_atmos_dyn_fvm_flux_cd4::atmos_dyn_fvm_fluxx_xvz_cd4(), scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxx_xvz_cd6(), scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxx_xvz_ud1(), scale_atmos_dyn_fvm_flux_ud3::atmos_dyn_fvm_fluxx_xvz_ud3(), scale_atmos_dyn_fvm_flux_ud3koren1993::atmos_dyn_fvm_fluxx_xvz_ud3koren1993(), scale_atmos_dyn_fvm_flux_ud5::atmos_dyn_fvm_fluxx_xvz_ud5(), scale_atmos_dyn_fvm_flux_cd2::atmos_dyn_fvm_fluxx_xyw_cd2(), scale_atmos_dyn_fvm_flux_cd4::atmos_dyn_fvm_fluxx_xyw_cd4(), scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxx_xyw_cd6(), scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxx_xyw_ud1(), scale_atmos_dyn_fvm_flux_ud3::atmos_dyn_fvm_fluxx_xyw_ud3(), scale_atmos_dyn_fvm_flux_ud3koren1993::atmos_dyn_fvm_fluxx_xyw_ud3koren1993(), scale_atmos_dyn_fvm_flux_ud5::atmos_dyn_fvm_fluxx_xyw_ud5(), scale_atmos_dyn_fvm_flux_cd2::atmos_dyn_fvm_fluxx_xyz_cd2(), scale_atmos_dyn_fvm_flux_cd4::atmos_dyn_fvm_fluxx_xyz_cd4(), scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxx_xyz_cd6(), scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxx_xyz_ud1(), scale_atmos_dyn_fvm_flux_ud3::atmos_dyn_fvm_fluxx_xyz_ud3(), scale_atmos_dyn_fvm_flux_ud3koren1993::atmos_dyn_fvm_fluxx_xyz_ud3koren1993(), scale_atmos_dyn_fvm_flux_ud5::atmos_dyn_fvm_fluxx_xyz_ud5(), scale_atmos_dyn_fvm_flux_cd2::atmos_dyn_fvm_fluxy_uyz_cd2(), scale_atmos_dyn_fvm_flux_cd4::atmos_dyn_fvm_fluxy_uyz_cd4(), scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxy_uyz_cd6(), scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxy_uyz_ud1(), scale_atmos_dyn_fvm_flux_ud3::atmos_dyn_fvm_fluxy_uyz_ud3(), scale_atmos_dyn_fvm_flux_ud3koren1993::atmos_dyn_fvm_fluxy_uyz_ud3koren1993(), scale_atmos_dyn_fvm_flux_ud5::atmos_dyn_fvm_fluxy_uyz_ud5(), scale_atmos_dyn_fvm_flux_cd2::atmos_dyn_fvm_fluxy_xvz_cd2(), scale_atmos_dyn_fvm_flux_cd4::atmos_dyn_fvm_fluxy_xvz_cd4(), scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxy_xvz_cd6(), scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxy_xvz_ud1(), scale_atmos_dyn_fvm_flux_ud3::atmos_dyn_fvm_fluxy_xvz_ud3(), scale_atmos_dyn_fvm_flux_ud3koren1993::atmos_dyn_fvm_fluxy_xvz_ud3koren1993(), scale_atmos_dyn_fvm_flux_ud5::atmos_dyn_fvm_fluxy_xvz_ud5(), scale_atmos_dyn_fvm_flux_cd2::atmos_dyn_fvm_fluxy_xyw_cd2(), scale_atmos_dyn_fvm_flux_cd4::atmos_dyn_fvm_fluxy_xyw_cd4(), scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxy_xyw_cd6(), scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxy_xyw_ud1(), scale_atmos_dyn_fvm_flux_ud3::atmos_dyn_fvm_fluxy_xyw_ud3(), scale_atmos_dyn_fvm_flux_ud3koren1993::atmos_dyn_fvm_fluxy_xyw_ud3koren1993(), scale_atmos_dyn_fvm_flux_ud5::atmos_dyn_fvm_fluxy_xyw_ud5(), scale_atmos_dyn_fvm_flux_cd2::atmos_dyn_fvm_fluxy_xyz_cd2(), scale_atmos_dyn_fvm_flux_cd4::atmos_dyn_fvm_fluxy_xyz_cd4(), scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxy_xyz_cd6(), scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxy_xyz_ud1(), scale_atmos_dyn_fvm_flux_ud3::atmos_dyn_fvm_fluxy_xyz_ud3(), scale_atmos_dyn_fvm_flux_ud3koren1993::atmos_dyn_fvm_fluxy_xyz_ud3koren1993(), scale_atmos_dyn_fvm_flux_ud5::atmos_dyn_fvm_fluxy_xyz_ud5(), scale_atmos_dyn_fvm_flux_cd2::atmos_dyn_fvm_fluxz_uyz_cd2(), scale_atmos_dyn_fvm_flux_cd4::atmos_dyn_fvm_fluxz_uyz_cd4(), scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxz_uyz_cd6(), scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxz_uyz_ud1(), scale_atmos_dyn_fvm_flux_ud3::atmos_dyn_fvm_fluxz_uyz_ud3(), scale_atmos_dyn_fvm_flux_ud3koren1993::atmos_dyn_fvm_fluxz_uyz_ud3koren1993(), scale_atmos_dyn_fvm_flux_ud5::atmos_dyn_fvm_fluxz_uyz_ud5(), scale_atmos_dyn_fvm_flux_cd2::atmos_dyn_fvm_fluxz_xvz_cd2(), scale_atmos_dyn_fvm_flux_cd4::atmos_dyn_fvm_fluxz_xvz_cd4(), scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxz_xvz_cd6(), scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxz_xvz_ud1(), scale_atmos_dyn_fvm_flux_ud3::atmos_dyn_fvm_fluxz_xvz_ud3(), scale_atmos_dyn_fvm_flux_ud3koren1993::atmos_dyn_fvm_fluxz_xvz_ud3koren1993(), scale_atmos_dyn_fvm_flux_ud5::atmos_dyn_fvm_fluxz_xvz_ud5(), scale_atmos_dyn_fvm_flux_cd2::atmos_dyn_fvm_fluxz_xyw_cd2(), scale_atmos_dyn_fvm_flux_cd4::atmos_dyn_fvm_fluxz_xyw_cd4(), scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxz_xyw_cd6(), scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxz_xyw_ud1(), scale_atmos_dyn_fvm_flux_ud3::atmos_dyn_fvm_fluxz_xyw_ud3(), scale_atmos_dyn_fvm_flux_ud3koren1993::atmos_dyn_fvm_fluxz_xyw_ud3koren1993(), scale_atmos_dyn_fvm_flux_ud5::atmos_dyn_fvm_fluxz_xyw_ud5(), scale_atmos_dyn_fvm_flux_cd2::atmos_dyn_fvm_fluxz_xyz_cd2(), scale_atmos_dyn_fvm_flux_cd4::atmos_dyn_fvm_fluxz_xyz_cd4(), scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxz_xyz_cd6(), scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxz_xyz_ud1(), scale_atmos_dyn_fvm_flux_ud3::atmos_dyn_fvm_fluxz_xyz_ud3(), scale_atmos_dyn_fvm_flux_ud3koren1993::atmos_dyn_fvm_fluxz_xyz_ud3koren1993(), scale_atmos_dyn_fvm_flux_ud5::atmos_dyn_fvm_fluxz_xyz_ud5(), scale_atmos_dyn_tstep_short_fvm_heve::atmos_dyn_tstep_short_fvm_heve(), scale_atmos_dyn_tstep_short_fvm_hevi::atmos_dyn_tstep_short_fvm_hevi(), scale_atmos_dyn_tstep_short_fvm_hivi::atmos_dyn_tstep_short_fvm_hivi(), scale_atmos_phy_tb_common::atmos_phy_tb_calc_flux_phi(), scale_atmos_phy_tb_common::atmos_phy_tb_calc_strain_tensor(), scale_atmos_phy_tb_d1980::atmos_phy_tb_d1980(), scale_atmos_phy_tb_smg::atmos_phy_tb_smg(), scale_atmos_dyn_common::calc_diff3(), scale_atmos_dyn_tstep_short_fvm_hivi::make_matrix(), and scale_atmos_dyn_tstep_short_fvm_hivi::solve_bicgstab().

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
subroutine, public prc_mpistop
Abort MPI.
real(rp), public const_undef
Definition: scale_const.F90:43
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
Here is the call graph for this function:

◆ valcheck_1d()

subroutine scale_debug::valcheck_1d ( real(rp), dimension(:), intent(in)  var,
real(rp), intent(in)  valmin,
real(rp), intent(in)  valmax,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  current_file,
integer, intent(in)  current_line 
)

Nan & extreme value checker (1D)

Definition at line 92 of file scale_debug.F90.

References debug_domain_num, scale_process::prc_mpistop(), scale_process::prc_myrank, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

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
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
integer, public prc_myrank
process num in local communicator
Here is the call graph for this function:

◆ valcheck_2d()

subroutine scale_debug::valcheck_2d ( real(rp), dimension(:,:), intent(in)  var,
real(rp), intent(in)  valmin,
real(rp), intent(in)  valmax,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  current_file,
integer, intent(in)  current_line 
)

Nan & extreme value checker (2D)

Definition at line 149 of file scale_debug.F90.

References debug_domain_num, scale_process::prc_mpistop(), scale_process::prc_myrank, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

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
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
integer, public prc_myrank
process num in local communicator
Here is the call graph for this function:

◆ valcheck_3d()

subroutine scale_debug::valcheck_3d ( real(rp), dimension(:,:,:), intent(in)  var,
real(rp), intent(in)  valmin,
real(rp), intent(in)  valmax,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  current_file,
integer, intent(in)  current_line 
)

Nan & extreme value checker (3D)

Definition at line 212 of file scale_debug.F90.

References debug_domain_num, scale_process::prc_mpistop(), scale_process::prc_myrank, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

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
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
integer, public prc_myrank
process num in local communicator
Here is the call graph for this function:

Variable Documentation

◆ debug_domain_num

integer, public scale_debug::debug_domain_num = 0

Definition at line 40 of file scale_debug.F90.

Referenced by scale_grid_nest::nest_setup(), valcheck_1d(), valcheck_2d(), and valcheck_3d().

40  integer, public :: debug_domain_num = 0