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 (KA, KS, KE, var, valmin, valmax, varname, current_file, current_line, mask)
 Nan & extreme value checker (1D) More...
 
subroutine valcheck_2d (IA, IS, IE, JA, JS, JE, var, valmin, valmax, varname, current_file, current_line, mask)
 Nan & extreme value checker (2D) More...
 
subroutine valcheck_3d (KA, KS, KE, IA, IS, IE, JA, JS, JE, var, valmin, valmax, varname, current_file, current_line, mask)
 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 59 of file scale_debug.F90.

59  use scale_prc, only: &
60  prc_abort
61  use scale_const, only: &
63  implicit none
64 
65  integer, intent(in) :: current_line
66  real(RP), intent(in) :: v
67  !---------------------------------------------------------------------------
68 
69  call prof_rapstart('Debug', 1)
70 
71  if ( .NOT. ( abs(v) < abs(const_undef) ) ) then
72  log_error("CHECK",*) 'uninitialized value at line:', current_line
73  call abort
74  end if
75 
76  call prof_rapend ('Debug', 1)
77 
78  return

References scale_const::const_undef, scale_prc::prc_abort(), scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by scale_atmos_dyn_fvm_numfilter::atmos_dyn_fvm_apply_numfilter(), scale_atmos_dyn_fvm_fct::atmos_dyn_fvm_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_cd8::atmos_dyn_fvm_flux_valuew_z_cd8(), 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_ud7::atmos_dyn_fvm_flux_valuew_z_ud7(), 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_cd8::atmos_dyn_fvm_fluxx_uyz_cd8(), 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_ud7::atmos_dyn_fvm_fluxx_uyz_ud7(), 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_cd8::atmos_dyn_fvm_fluxx_xvz_cd8(), 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_ud7::atmos_dyn_fvm_fluxx_xvz_ud7(), 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_cd8::atmos_dyn_fvm_fluxx_xyw_cd8(), 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_ud7::atmos_dyn_fvm_fluxx_xyw_ud7(), 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_cd8::atmos_dyn_fvm_fluxx_xyz_cd8(), 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_ud7::atmos_dyn_fvm_fluxx_xyz_ud7(), 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_cd8::atmos_dyn_fvm_fluxy_uyz_cd8(), 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_ud7::atmos_dyn_fvm_fluxy_uyz_ud7(), 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_cd8::atmos_dyn_fvm_fluxy_xvz_cd8(), 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_ud7::atmos_dyn_fvm_fluxy_xvz_ud7(), 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_cd8::atmos_dyn_fvm_fluxy_xyw_cd8(), 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_ud7::atmos_dyn_fvm_fluxy_xyw_ud7(), 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_cd8::atmos_dyn_fvm_fluxy_xyz_cd8(), 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_ud7::atmos_dyn_fvm_fluxy_xyz_ud7(), 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_cd8::atmos_dyn_fvm_fluxz_uyz_cd8(), 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_ud7::atmos_dyn_fvm_fluxz_uyz_ud7(), 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_cd8::atmos_dyn_fvm_fluxz_xvz_cd8(), 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_ud7::atmos_dyn_fvm_fluxz_xvz_ud7(), 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_cd8::atmos_dyn_fvm_fluxz_xyw_cd8(), 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_ud7::atmos_dyn_fvm_fluxz_xyw_ud7(), 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_cd8::atmos_dyn_fvm_fluxz_xyz_cd8(), 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_fvm_flux_ud7::atmos_dyn_fvm_fluxz_xyz_ud7(), 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(), mod_atmos_vars::atmos_vars_check(), mod_land_vars::land_vars_check(), scale_atmos_dyn_tstep_short_fvm_hivi::make_matrix(), mod_ocean_vars::ocean_vars_check(), scale_atmos_dyn_tstep_short_fvm_hivi::solve_bicgstab(), and mod_urban_vars::urban_vars_check().

Here is the call graph for this function:

◆ valcheck_1d()

subroutine scale_debug::valcheck_1d ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
real(rp), dimension(ka), 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,
logical, dimension(ka), intent(in), optional  mask 
)

Nan & extreme value checker (1D)

Definition at line 92 of file scale_debug.F90.

92  use scale_prc, only: &
93  prc_abort, &
95  implicit none
96  integer , intent(in) :: KA, KS, KE
97  real(RP), intent(in) :: var(KA)
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, intent(in), optional :: mask(KA)
105 
106  logical :: invalid_value
107  integer :: k
108  !---------------------------------------------------------------------------
109 
110  call prof_rapstart('Debug', 1)
111 
112  !$acc update host(var) if(acc_is_present(var))
113 
114  invalid_value = .false.
115  if ( present(mask) ) then
116  !$acc update host(mask) if(acc_is_present(mask))
117  do k = ks, ke
118  if ( .not. mask(k) ) cycle
119  if ( var(k)*0.0_rp /= 0.0_rp &
120  .OR. var(k) < valmin &
121  .OR. var(k) > valmax ) then
122  invalid_value = .true.
123  exit
124  endif
125  enddo
126  else
127  do k = ks, ke
128  if ( var(k)*0.0_rp /= 0.0_rp &
129  .OR. var(k) < valmin &
130  .OR. var(k) > valmax ) then
131  invalid_value = .true.
132  exit
133  endif
134  enddo
135  end if
136 
137  if ( invalid_value ) then
138  log_error("VALCHECK_1D",*) 'invalid value): ', trim(varname), &
139  '(', k, ')=', var(k)
140  log_error_cont(*) 'in file : ', trim(current_file), ', at line : ', current_line
141  call prc_abort
142  endif
143 
144  call prof_rapend ('Debug', 1)
145 
146  return

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

Here is the call graph for this function:

◆ valcheck_2d()

subroutine scale_debug::valcheck_2d ( integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension(ia,ja), 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,
logical, dimension(ia,ja), intent(in), optional  mask 
)

Nan & extreme value checker (2D)

Definition at line 160 of file scale_debug.F90.

160  use scale_prc, only: &
161  prc_abort, &
162  prc_myrank
163  implicit none
164  integer , intent(in) :: IA, IS, IE
165  integer , intent(in) :: JA, JS, JE
166  real(RP), intent(in) :: var(IA,JA)
167  real(RP), intent(in) :: valmin
168  real(RP), intent(in) :: valmax
169  character(len=*), intent(in) :: varname
170  character(len=*), intent(in) :: current_file
171  integer, intent(in) :: current_line
172 
173  logical, intent(in), optional :: mask(IA,JA)
174 
175  logical :: invalid_value
176  integer :: i, j
177  !---------------------------------------------------------------------------
178 
179  call prof_rapstart('Debug', 1)
180 
181  !$acc update host(var) if(acc_is_present(var))
182 
183  invalid_value = .false.
184  if ( present(mask) ) then
185  !$acc update host(mask) if(acc_is_present(mask))
186  outer1:do j = js, je
187  do i = is, ie
188  if ( .not. mask(i,j) ) cycle
189  if ( var(i,j)*0.0_rp /= 0.0_rp &
190  .OR. var(i,j) < valmin &
191  .OR. var(i,j) > valmax ) then
192  invalid_value = .true.
193  exit outer1
194  endif
195  enddo
196  enddo outer1
197  else
198  outer2:do j = js, je
199  do i = is, ie
200  if ( var(i,j)*0.0_rp /= 0.0_rp &
201  .OR. var(i,j) < valmin &
202  .OR. var(i,j) > valmax ) then
203  invalid_value = .true.
204  exit outer2
205  endif
206  enddo
207  enddo outer2
208  end if
209 
210  if ( invalid_value ) then
211  log_error("VALCHECK_2D",*) 'invalid value:', trim(varname), &
212  '(', i, ',', j, ')=', var(i,j)
213  log_error_cont(*) 'in file : ', trim(current_file), ', at line : ', current_line
214  call prc_abort
215  endif
216 
217  call prof_rapend ('Debug', 1)
218 
219  return

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

Here is the call graph for this function:

◆ valcheck_3d()

subroutine scale_debug::valcheck_3d ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension(ka,ia,ja), 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,
logical, dimension(ia,ja), intent(in), optional  mask 
)

Nan & extreme value checker (3D)

Definition at line 235 of file scale_debug.F90.

235  use scale_prc, only: &
236  prc_abort, &
237  prc_myrank
238  implicit none
239  integer , intent(in) :: KA, KS, KE
240  integer , intent(in) :: IA, IS, IE
241  integer , intent(in) :: JA, JS, JE
242  real(RP), intent(in) :: var(KA,IA,JA)
243  real(RP), intent(in) :: valmin
244  real(RP), intent(in) :: valmax
245  character(len=*), intent(in) :: varname
246  character(len=*), intent(in) :: current_file
247  integer, intent(in) :: current_line
248 
249  logical, intent(in), optional :: mask(IA,JA)
250 
251  logical :: invalid_value
252  integer :: k, i, j
253  !---------------------------------------------------------------------------
254 
255  call prof_rapstart('Debug', 1)
256 
257  !$acc update host(var) if(acc_is_present(var))
258 
259  invalid_value = .false.
260  if ( present(mask) ) then
261  !$acc update host(mask) if(acc_is_present(mask))
262  outer1:do j = js, je
263  do i = is, ie
264  if ( .not. mask(i,j) ) cycle
265  do k = ks, ke
266  if ( var(k,i,j)*0.0_rp /= 0.0_rp &
267  .OR. var(k,i,j) < valmin &
268  .OR. var(k,i,j) > valmax ) then
269  invalid_value = .true.
270  exit outer1
271  endif
272  enddo
273  enddo
274  enddo outer1
275  else
276  outer2:do j = js, je
277  do i = is, ie
278  do k = ks, ke
279  if ( var(k,i,j)*0.0_rp /= 0.0_rp &
280  .OR. var(k,i,j) < valmin &
281  .OR. var(k,i,j) > valmax ) then
282  invalid_value = .true.
283  exit outer2
284  endif
285  enddo
286  enddo
287  enddo outer2
288  end if
289 
290  if ( invalid_value ) then
291  log_error("VALCHECK_3D",*) 'Invalid value:', trim(varname), &
292  '(', k, ',', i, ',', j, ')=', var(k,i,j)
293  log_error_cont(*) 'in file : ', trim(current_file), ', at line : ', current_line
294  call prc_abort
295  endif
296 
297  call prof_rapend ('Debug', 1)
298 
299  return

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

Here is the call graph for this function:

Variable Documentation

◆ debug_domain_num

integer, public scale_debug::debug_domain_num = 0

Definition at line 41 of file scale_debug.F90.

41  integer, public :: DEBUG_DOMAIN_NUM = 0

Referenced by scale_comm_cartesc_nest::comm_cartesc_nest_setup().

scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_prc::prc_myrank
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:91
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56