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

56  use scale_prc, only: &
57  prc_abort
58  use scale_const, only: &
59  const_undef
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

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

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

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

226  use scale_prc, only: &
227  prc_abort, &
228  prc_myrank
229  implicit none
230  integer , intent(in) :: KA, KS, KE
231  integer , intent(in) :: IA, IS, IE
232  integer , intent(in) :: JA, JS, JE
233  real(RP), intent(in) :: var(KA,IA,JA)
234  real(RP), intent(in) :: valmin
235  real(RP), intent(in) :: valmax
236  character(len=*), intent(in) :: varname
237  character(len=*), intent(in) :: current_file
238  integer, intent(in) :: current_line
239 
240  logical, intent(in), optional :: mask(IA,JA)
241 
242  logical :: invalid_value
243  integer :: k, i, j
244  !---------------------------------------------------------------------------
245 
246  call prof_rapstart('Debug', 1)
247 
248  invalid_value = .false.
249  if ( present(mask) ) then
250  outer1:do j = js, je
251  do i = is, ie
252  if ( .not. mask(i,j) ) cycle
253  do k = ks, ke
254  if ( var(k,i,j)*0.0_rp /= 0.0_rp &
255  .OR. var(k,i,j) < valmin &
256  .OR. var(k,i,j) > valmax ) then
257  invalid_value = .true.
258  exit outer1
259  endif
260  enddo
261  enddo
262  enddo outer1
263  else
264  outer2:do j = js, je
265  do i = is, ie
266  do k = ks, ke
267  if ( var(k,i,j)*0.0_rp /= 0.0_rp &
268  .OR. var(k,i,j) < valmin &
269  .OR. var(k,i,j) > valmax ) then
270  invalid_value = .true.
271  exit outer2
272  endif
273  enddo
274  enddo
275  enddo outer2
276  end if
277 
278  if ( invalid_value ) then
279  log_error("VALCHECK_3D",*) 'Invalid value:', trim(varname), &
280  '(', k, ',', i, ',', j, ')=', var(k,i,j)
281  log_error_cont(*) 'in file : ', trim(current_file), ', at line : ', current_line
282  call prc_abort
283  endif
284 
285  call prof_rapend ('Debug', 1)
286 
287  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 38 of file scale_debug.F90.

38  integer, public :: DEBUG_DOMAIN_NUM = 0

Referenced by scale_comm_cartesc_nest::comm_cartesc_nest_setup().

scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_prc::prc_myrank
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:90
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_const
module CONSTANT
Definition: scale_const.F90:11