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: &
61  const_undef
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
subroutine, public prc_mpistop
Abort MPI.
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 89 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().

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
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 142 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().

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
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 201 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().

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
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