218 lines
7.8 KiB
Fortran
218 lines
7.8 KiB
Fortran
! This module contains the error codes handled by the gray_main subroutine
|
||
!
|
||
! An error is actually a 32bit integer bitmask that can encode several
|
||
! error cases. The individual errors can be raised using the `raise_error`
|
||
! function and combined simply with the intrinsic `ior` function.
|
||
module gray_errors
|
||
|
||
use logger, only : log_warning
|
||
use, intrinsic :: iso_fortran_env, only : int32
|
||
|
||
implicit none
|
||
|
||
! The type of a GRAY error is `integer(kind=gray_error)`
|
||
integer, parameter :: gray_error = int32
|
||
|
||
! Specification of a GRAY error
|
||
type error_spec
|
||
integer :: offset ! bitmask offset
|
||
integer :: subcases ! number of error subcases (max: 10)
|
||
character(32) :: mod ! module raising the error
|
||
character(32) :: proc ! procedure raising the error
|
||
character(64) :: msg(0:9) ! error message (one for each subcase)
|
||
end type
|
||
|
||
! macros used for defining errors
|
||
#define after(x) x%offset + x%subcases
|
||
#define str character(64)
|
||
#define pad1 '', '', '', '', '', '', '', '', ''
|
||
#define pad2 '', '', '', '', '', '', '', ''
|
||
#define pad3 '', '', '', '', '', '', ''
|
||
|
||
! All GRAY errors
|
||
|
||
type(error_spec), parameter :: unstable_beam = &
|
||
error_spec(offset=0, subcases=1, &
|
||
mod='gray_core', proc='gray_main', &
|
||
msg=[str :: 'beamtracing may be unstable', pad1])
|
||
|
||
type(error_spec), parameter :: dielectric_tensor = &
|
||
error_spec(offset=after(unstable_beam), subcases=2, &
|
||
mod='gray_core', proc='gray_main', &
|
||
msg=[str :: 'ε tensor, overflow in `fsup`', &
|
||
'ε tensor, integration error in `hermitian_2`', &
|
||
pad2])
|
||
|
||
type(error_spec), parameter :: warmdisp_convergence = &
|
||
error_spec(offset=after(dielectric_tensor), subcases=2, &
|
||
mod='dispersion', proc='warmdisp', &
|
||
msg=[str :: 'failed to converge, returned fallback value', &
|
||
'failed to converge, returned last value', &
|
||
pad2])
|
||
|
||
type(error_spec), parameter :: warmdisp_result = &
|
||
error_spec(offset=after(warmdisp_convergence), subcases=2, &
|
||
mod='dispersion', proc='warmdisp', &
|
||
msg=[str :: 'final N⊥² is NaN or ±Infinity', &
|
||
'final N⊥² in 3rd quadrant', &
|
||
pad2])
|
||
|
||
type(error_spec), parameter :: negative_absorption = &
|
||
error_spec(offset=after(warmdisp_result), subcases=1, &
|
||
mod='gray_core', proc='alpha_effj', &
|
||
msg=[str :: 'negative absorption coeff.', pad1])
|
||
|
||
type(error_spec), parameter :: fpp_integration = &
|
||
error_spec(offset=after(negative_absorption), subcases=1, &
|
||
mod='eccd', proc='eccdeff', &
|
||
msg=[str :: 'fpp integration error', pad1])
|
||
|
||
type(error_spec), parameter :: fcur_integration = &
|
||
error_spec(offset=after(fpp_integration), subcases=3, &
|
||
mod='eccd', proc='eccdeff', &
|
||
msg=[str :: 'fcur integration error (no trapping)', &
|
||
'fcur integration error (1st trapping region)', &
|
||
'fcur integration error (2st trapping region)', &
|
||
pad3])
|
||
|
||
! Errors occuring during raytracing
|
||
type(error_spec), parameter :: raytracing_errors(*) = [unstable_beam]
|
||
|
||
! Errors occuring during absorption and current drive computations
|
||
type(error_spec), parameter :: ecrh_cd_errors(*) = &
|
||
[dielectric_tensor, &
|
||
warmdisp_convergence, &
|
||
warmdisp_result, &
|
||
negative_absorption, &
|
||
fpp_integration, &
|
||
fcur_integration]
|
||
|
||
contains
|
||
|
||
pure function is_critical(error)
|
||
! Checks whether critical errors have occurred
|
||
|
||
! subroutines arguments
|
||
integer(kind=gray_error), intent(in) :: error
|
||
logical :: is_critical
|
||
|
||
is_critical = has_error(error, negative_absorption)
|
||
end function is_critical
|
||
|
||
|
||
pure function has_error(error, spec)
|
||
! Checks whether the `error` bitmask contains the error given by `spec`
|
||
|
||
! function arguments
|
||
integer(kind=gray_error), intent(in) :: error
|
||
type(error_spec), intent(in) :: spec
|
||
logical :: has_error
|
||
|
||
has_error = ibits(error, spec%offset, spec%subcases) == 1
|
||
end function has_error
|
||
|
||
|
||
pure function has_new_errors(error0, error1)
|
||
! Checks whether the `error1` bitmask contains errors
|
||
! not already present in `error0`
|
||
|
||
! function arguments
|
||
integer(kind=gray_error), intent(in) :: error0, error1
|
||
logical :: has_new_errors
|
||
|
||
! check if any bit flipped from 0→1
|
||
has_new_errors = iand(not(error0), error1) /= 0
|
||
end function has_new_errors
|
||
|
||
|
||
pure function raise_error(error, spec, subcase)
|
||
! Raise the bits of error `spec` (with optional `subcase` number)
|
||
! in the `error` bitmask.
|
||
|
||
! function arguments
|
||
integer(kind=gray_error), intent(in) :: error
|
||
type(error_spec), intent(in) :: spec
|
||
integer, intent(in), optional :: subcase
|
||
integer(kind=gray_error) :: raise_error
|
||
|
||
raise_error = ibset(error, spec%offset &
|
||
+ merge(subcase, 0, present(subcase)))
|
||
end function raise_error
|
||
|
||
|
||
subroutine print_err_raytracing(error, step, ray)
|
||
! Pretty prints raytracing errors
|
||
!
|
||
! The error and some context (integration step, ray)
|
||
! is logged using the logger module.
|
||
|
||
! subroutines arguments
|
||
integer, intent(in) :: error, step, ray
|
||
|
||
! local variables
|
||
integer :: slice ! a slice of the bitmask
|
||
character(256) :: line ! formatted log line
|
||
type(error_spec) :: spec
|
||
integer :: i, j
|
||
|
||
! format specifier of the log line
|
||
character(*), parameter :: fmt = &
|
||
'(a,": ","error=",g0," step=",g0," ray=",g0)'
|
||
|
||
! iterate on the known errors
|
||
do i = 1, size(raytracing_errors)
|
||
spec = raytracing_errors(i)
|
||
slice = ibits(error, spec%offset, spec%subcases)
|
||
if (slice == 0) cycle
|
||
|
||
! iterate on the subcases
|
||
do j = 0, spec%subcases - 1
|
||
if (ibits(slice, j, 1) == 0) cycle
|
||
write(line, fmt) trim(spec%msg(j)), slice * 2**j, step, ray
|
||
call log_warning(line, mod=spec%mod, proc=spec%proc)
|
||
end do
|
||
end do
|
||
end subroutine print_err_raytracing
|
||
|
||
|
||
subroutine print_err_ecrh_cd(error, step, ray, Npr, alpha)
|
||
! Pretty prints ECRH & CD errors
|
||
!
|
||
! The error and some context (integration step, ray, N⊥, α)
|
||
! is logged using the logger module.
|
||
|
||
use const_and_precisions, only : wp_
|
||
|
||
! subroutines arguments
|
||
integer, intent(in) :: error, step, ray
|
||
complex(wp_), intent(in) :: Npr
|
||
real(wp_), intent(in) :: alpha
|
||
|
||
! local variables
|
||
integer :: slice ! a slice of the bitmask
|
||
character(256) :: line ! formatted log line
|
||
type(error_spec) :: spec
|
||
integer :: i, j
|
||
|
||
! format specifier of the log line
|
||
character(*), parameter :: fmt = &
|
||
'(a,": ", "error=",g0, " N⊥=",f0.0,sp,f0.0,"i", " α=",g0, " step=",g0," ray=",g0)'
|
||
|
||
! iterate on the known errors
|
||
do i = 1, size(ecrh_cd_errors)
|
||
spec = ecrh_cd_errors(i)
|
||
slice = ibits(error, spec%offset, spec%subcases)
|
||
if (slice == 0) cycle
|
||
|
||
! iterate on the subcases
|
||
do j = 0, spec%subcases - 1
|
||
if (ibits(slice, j, 1) == 0) cycle
|
||
write(line, fmt) trim(spec%msg(j)), slice * 2**j, &
|
||
Npr%re, Npr%im, alpha, step, ray
|
||
call log_warning(line, mod=spec%mod, proc=spec%proc)
|
||
end do
|
||
end do
|
||
end subroutine print_err_ecrh_cd
|
||
|
||
end module gray_errors
|