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

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 56 of file scale_debug.F90.

References scale_const::const_undef, scale_prc::prc_abort(), 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().

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
real(rp), public const_undef
Definition: scale_const.F90:41
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
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 87 of file scale_debug.F90.

References debug_domain_num, scale_prc::prc_abort(), scale_prc::prc_myrank, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

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
module PROCESS
Definition: scale_prc.F90:11
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
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 140 of file scale_debug.F90.

References debug_domain_num, scale_prc::prc_abort(), scale_prc::prc_myrank, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

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
module PROCESS
Definition: scale_prc.F90:11
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
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 199 of file scale_debug.F90.

References debug_domain_num, scale_prc::prc_abort(), scale_prc::prc_myrank, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

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
module PROCESS
Definition: scale_prc.F90:11
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
Here is the call graph for this function:

Variable Documentation

◆ debug_domain_num

integer, public scale_debug::debug_domain_num = 0

Definition at line 38 of file scale_debug.F90.

Referenced by scale_comm_cartesc_nest::comm_cartesc_nest_setup(), valcheck_1d(), valcheck_2d(), and valcheck_3d().

38  integer, public :: debug_domain_num = 0