call exit on errors only in main

This commit is contained in:
Michele Guerini Rocco 2023-09-12 16:29:09 +02:00
parent daf3d500af
commit 7ed76959d4
Signed by: rnhmjoj
GPG Key ID: BFBAF4C975F76450
7 changed files with 122 additions and 83 deletions

View File

@ -4,7 +4,7 @@ module beams
contains contains
subroutine read_beam0(params, unit) subroutine read_beam0(params, err, unit)
! Reads the wave launcher parameters for the simple case ! Reads the wave launcher parameters for the simple case
! where w(z) and 1/R(z) are fixed. ! where w(z) and 1/R(z) are fixed.
! !
@ -34,12 +34,12 @@ contains
implicit none implicit none
! subroutine arguments ! subroutine arguments
type(antenna_parameters), intent(inout) :: params type(antenna_parameters), intent(inout) :: params
integer, intent(in), optional :: unit integer, intent(out) :: err
integer, intent(in), optional :: unit
! local variables ! local variables
integer :: u integer :: u
integer :: err
real(wp_) :: k0, w0(2), z0(2), z_R(2), phi real(wp_) :: k0, w0(2), z0(2), z_R(2), phi
u = get_free_unit(unit) u = get_free_unit(unit)
@ -48,7 +48,7 @@ contains
if (err /= 0) then if (err /= 0) then
call log_error('opening beams file ('//trim(params%filenm)//') failed!', & call log_error('opening beams file ('//trim(params%filenm)//') failed!', &
mod='beams', proc="read_beam0") mod='beams', proc="read_beam0")
call exit(1) return
end if end if
read(u, *) params%fghz ! Wave frequency (GHz) read(u, *) params%fghz ! Wave frequency (GHz)
@ -78,7 +78,7 @@ contains
end subroutine read_beam0 end subroutine read_beam0
subroutine read_beam1(params, unit) subroutine read_beam1(params, err, unit)
! Reads the wave launcher parameters for the case ! Reads the wave launcher parameters for the case
! where w(z, α) and 1/R(z, α) depend on the launcher angle α. ! where w(z, α) and 1/R(z, α) depend on the launcher angle α.
! !
@ -105,7 +105,8 @@ contains
! subroutine arguments ! subroutine arguments
type(antenna_parameters), intent(inout) :: params type(antenna_parameters), intent(inout) :: params
integer, intent(in), optional :: unit integer, intent(out) :: err
integer, optional, intent(in) :: unit
! local variables ! local variables
integer :: u, nisteer, i, k, ii integer :: u, nisteer, i, k, ii
@ -116,7 +117,6 @@ contains
type(spline_simple) :: beta, waist1, waist2, & type(spline_simple) :: beta, waist1, waist2, &
rci1, rci2, phi1, phi2, & rci1, rci2, phi1, phi2, &
x0, y0, z0 x0, y0, z0
integer :: err
u = get_free_unit(unit) u = get_free_unit(unit)
@ -124,7 +124,7 @@ contains
if (err /= 0) then if (err /= 0) then
call log_error('opening beams file ('//trim(params%filenm)//') failed!', & call log_error('opening beams file ('//trim(params%filenm)//') failed!', &
mod='beams', proc="read_beam1") mod='beams', proc="read_beam1")
call exit(1) return
end if end if
read(u,*) params%fghz read(u,*) params%fghz
read(u,*) nisteer read(u,*) nisteer
@ -211,7 +211,7 @@ contains
end subroutine read_beam1 end subroutine read_beam1
subroutine read_beam2(params, beamid, unit) subroutine read_beam2(params, beamid, err, unit)
! Reads the wave launcher parameters for the general case ! Reads the wave launcher parameters for the general case
! where w(z, α, β) and 1/R(z, α, β) depend on the launcher angles α, β. ! where w(z, α, β) and 1/R(z, α, β) depend on the launcher angles α, β.
@ -225,8 +225,9 @@ contains
! subroutine arguments ! subroutine arguments
type(antenna_parameters), intent(inout) :: params type(antenna_parameters), intent(inout) :: params
integer, intent(in) :: beamid integer, intent(in) :: beamid
integer, intent(in), optional :: unit integer, intent(out), optional :: err
integer, intent(in), optional :: unit
! local variables ! local variables
character(len=20) :: beamname character(len=20) :: beamname
@ -252,7 +253,6 @@ contains
real(wp_), dimension(1) :: fi real(wp_), dimension(1) :: fi
integer, parameter :: kspl=1 integer, parameter :: kspl=1
real(wp_), parameter :: sspl=0.01_wp_ real(wp_), parameter :: sspl=0.01_wp_
integer :: err
u = get_free_unit(unit) u = get_free_unit(unit)
@ -260,7 +260,7 @@ contains
if (err /= 0) then if (err /= 0) then
call log_error('opening beams file ('//trim(params%filenm)//') failed!', & call log_error('opening beams file ('//trim(params%filenm)//') failed!', &
mod='beams', proc="read_beam1") mod='beams', proc="read_beam1")
call exit(1) return
end if end if
!======================================================================================= !=======================================================================================

View File

@ -190,7 +190,7 @@ contains
end function fzeff end function fzeff
subroutine read_profiles(filenm, data, unit) subroutine read_profiles(filenm, data, err, unit)
! Reads the radial plasma profiles from `file` and store them ! Reads the radial plasma profiles from `file` and store them
! into `data`. If given, the file is opened in the `unit` number. ! into `data`. If given, the file is opened in the `unit` number.
! Format notes: ! Format notes:
@ -206,11 +206,11 @@ contains
! subroutine arguments ! subroutine arguments
character(len=*), intent(in) :: filenm character(len=*), intent(in) :: filenm
type(profiles_data), intent(out) :: data type(profiles_data), intent(out) :: data
integer, intent(out) :: err
integer, optional, intent(in) :: unit integer, optional, intent(in) :: unit
! local variables ! local variables
integer :: u, i, nrows integer :: u, i, nrows
integer :: err
! Free the arrays when already allocated ! Free the arrays when already allocated
if (allocated(data%psrad)) deallocate(data%psrad) if (allocated(data%psrad)) deallocate(data%psrad)
@ -225,7 +225,8 @@ contains
if (err /= 0) then if (err /= 0) then
call log_error('opening profiles file ('//trim(filenm)//') failed!', & call log_error('opening profiles file ('//trim(filenm)//') failed!', &
mod='coreprofiles', proc="read_profiles") mod='coreprofiles', proc="read_profiles")
call exit(1) err = 1
return
end if end if
read(u, *) nrows read(u, *) nrows
@ -244,7 +245,7 @@ contains
end subroutine read_profiles end subroutine read_profiles
subroutine read_profiles_an(filenm, data, unit) subroutine read_profiles_an(filenm, data, err, unit)
! Reads the plasma profiles `data` in the analytical format ! Reads the plasma profiles `data` in the analytical format
! from params%filenm. ! from params%filenm.
! If given, the file is opened in the `unit` number. ! If given, the file is opened in the `unit` number.
@ -266,11 +267,11 @@ contains
! subroutine arguments ! subroutine arguments
character(len=*), intent(in) :: filenm character(len=*), intent(in) :: filenm
type(profiles_data), intent(out) :: data type(profiles_data), intent(out) :: data
integer, intent(out) :: err
integer, optional, intent(in) :: unit integer, optional, intent(in) :: unit
! local variables ! local variables
integer :: u integer :: u
integer :: err
u = get_free_unit(unit) u = get_free_unit(unit)
@ -283,7 +284,8 @@ contains
if (err /= 0) then if (err /= 0) then
call log_error('opening profiles file ('//trim(filenm)//') failed!', & call log_error('opening profiles file ('//trim(filenm)//') failed!', &
mod='coreprofiles', proc='read_profiles_an') mod='coreprofiles', proc='read_profiles_an')
call exit(1) err = 1
return
end if end if
read (u,*) data%derad(1:3) ! dens0, n1, n2 read (u,*) data%derad(1:3) ! dens0, n1, n2
@ -338,13 +340,15 @@ contains
end subroutine scale_profiles end subroutine scale_profiles
subroutine set_profiles_spline(params, data, launch_pos) subroutine set_profiles_spline(params, data, err, launch_pos)
! Computes splines for the plasma profiles data and stores them ! Computes splines for the plasma profiles data and stores them
! in their respective global variables, see the top of this file. ! in their respective global variables, see the top of this file.
! !
! When `launch_pos` (cartesian launch coordinates in cm) is present, ! When `launch_pos` (cartesian launch coordinates in cm) is present,
! the subroutine will also check that the wave launcher is strictly ! the subroutine will also check that the wave launcher is strictly
! outside the reconstructed plasma density boundary. ! outside the reconstructed plasma density boundary.
!
! `err` is 1 if I/O errors occured, 2 if other initialisation failed.
use gray_params, only : profiles_parameters, profiles_data use gray_params, only : profiles_parameters, profiles_data
use logger, only : log_debug, log_info, log_warning, log_error use logger, only : log_debug, log_info, log_warning, log_error
@ -353,10 +357,11 @@ contains
! subroutine arguments ! subroutine arguments
type(profiles_parameters), intent(inout) :: params type(profiles_parameters), intent(inout) :: params
type(profiles_data), intent(inout) :: data type(profiles_data), intent(inout) :: data
integer, intent(out) :: err
real(wp_), optional, intent(in) :: launch_pos(3) real(wp_), optional, intent(in) :: launch_pos(3)
! local variables ! local variables
integer :: n, err integer :: n
! for log messages formatting ! for log messages formatting
character(256) :: msg character(256) :: msg
@ -380,7 +385,10 @@ contains
else if (err > 0) then else if (err > 0) then
write (msg, '(a, g0)') 'density fit failed with error ', err write (msg, '(a, g0)') 'density fit failed with error ', err
call log_error(msg, mod='coreprofiles', proc='density') call log_error(msg, mod='coreprofiles', proc='density')
call exit(1) err = 2
return
else
err = 0
end if end if
! Computation of the polynomial tail parameters ! Computation of the polynomial tail parameters
@ -465,7 +473,8 @@ contains
'wave launcher is inside the plasma! ', & 'wave launcher is inside the plasma! ', &
'launcher: ψ=', psi, ' boundary: ψ=', tail%end 'launcher: ψ=', psi, ' boundary: ψ=', tail%end
call log_error(msg, mod='coreprofiles', proc='set_profiles_spline') call log_error(msg, mod='coreprofiles', proc='set_profiles_spline')
call exit(2) err = 2
return
end if end if
end block end block
end if end if

View File

@ -65,7 +65,7 @@ module equilibrium
contains contains
subroutine read_eqdsk(params, data, unit) subroutine read_eqdsk(params, data, err, unit)
! Reads the MHD equilibrium `data` from a G-EQDSK file (params%filenm). ! Reads the MHD equilibrium `data` from a G-EQDSK file (params%filenm).
! If given, the file is opened in the `unit` number. ! If given, the file is opened in the `unit` number.
! For a description of the G-EQDSK, see the GRAY user manual. ! For a description of the G-EQDSK, see the GRAY user manual.
@ -79,6 +79,7 @@ contains
! subroutine arguments ! subroutine arguments
type(equilibrium_parameters), intent(in) :: params type(equilibrium_parameters), intent(in) :: params
type(equilibrium_data), intent(out) :: data type(equilibrium_data), intent(out) :: data
integer, intent(out) :: err
integer, optional, intent(in) :: unit integer, optional, intent(in) :: unit
! local variables ! local variables
@ -86,7 +87,6 @@ contains
character(len=48) :: string character(len=48) :: string
real(wp_) :: dr, dz, dps, rleft, zmid, zleft, psiedge, psiaxis real(wp_) :: dr, dz, dps, rleft, zmid, zleft, psiedge, psiaxis
real(wp_) :: xdum ! dummy variable, used to discard data real(wp_) :: xdum ! dummy variable, used to discard data
integer :: err
u = get_free_unit(unit) u = get_free_unit(unit)
@ -95,7 +95,8 @@ contains
if (err /= 0) then if (err /= 0) then
call log_error('opening eqdsk file ('//trim(params%filenm)//') failed!', & call log_error('opening eqdsk file ('//trim(params%filenm)//') failed!', &
mod='equilibrium', proc='read_eqdsk') mod='equilibrium', proc='read_eqdsk')
call exit(1) err = 1
return
end if end if
! get size of main arrays and allocate them ! get size of main arrays and allocate them
@ -191,7 +192,7 @@ contains
end subroutine read_eqdsk end subroutine read_eqdsk
subroutine read_equil_an(filenm, ipass, data, unit) subroutine read_equil_an(filenm, ipass, data, err, unit)
! Reads the MHD equilibrium `data` in the analytical format ! Reads the MHD equilibrium `data` in the analytical format
! from params%filenm. ! from params%filenm.
! If given, the file is opened in the `unit` number. ! If given, the file is opened in the `unit` number.
@ -207,11 +208,11 @@ contains
character(len=*), intent(in) :: filenm character(len=*), intent(in) :: filenm
integer, intent(in) :: ipass integer, intent(in) :: ipass
type(equilibrium_data), intent(out) :: data type(equilibrium_data), intent(out) :: data
integer, intent(out) :: err
integer, optional, intent(in) :: unit integer, optional, intent(in) :: unit
! local variables ! local variables
integer :: i, u, nlim integer :: i, u, nlim
integer :: err
real(wp_) :: rr0m, zr0m, rpam, b0 real(wp_) :: rr0m, zr0m, rpam, b0
u = get_free_unit(unit) u = get_free_unit(unit)
@ -220,7 +221,8 @@ contains
if (err /= 0) then if (err /= 0) then
call log_error('opening equilibrium file ('//trim(filenm)//') failed!', & call log_error('opening equilibrium file ('//trim(filenm)//') failed!', &
mod='equilibrium', proc='read_equil_an') mod='equilibrium', proc='read_equil_an')
call exit(1) err = 1
return
end if end if
read(u, *) rr0m, zr0m, rpam read(u, *) rr0m, zr0m, rpam
@ -374,7 +376,7 @@ contains
end subroutine scale_equil end subroutine scale_equil
subroutine set_equil_spline(params, data) subroutine set_equil_spline(params, data, err)
! Computes splines for the MHD equilibrium data and stores them ! Computes splines for the MHD equilibrium data and stores them
! in their respective global variables, see the top of this file. ! in their respective global variables, see the top of this file.
use const_and_precisions, only : zero, one use const_and_precisions, only : zero, one
@ -387,8 +389,9 @@ contains
implicit none implicit none
! subroutine arguments ! subroutine arguments
type(equilibrium_parameters), intent(in) :: params type(equilibrium_parameters), intent(in) :: params
type(equilibrium_data), intent(in) :: data type(equilibrium_data), intent(in) :: data
integer, intent(out) :: err
! local variables ! local variables
integer :: nr, nz, nrest, nzest, npsest, nrz, npsi, nbnd, ibinf, ibsup integer :: nr, nz, nrest, nzest, npsest, nrz, npsi, nbnd, ibinf, ibsup
@ -396,7 +399,7 @@ contains
real(wp_) :: rbmin, rbmax, rbinf, rbsup, r1, z1 real(wp_) :: rbmin, rbmax, rbinf, rbsup, r1, z1
real(wp_), dimension(size(data%psinr)) :: rhotn real(wp_), dimension(size(data%psinr)) :: rhotn
real(wp_), dimension(:), allocatable :: rv1d, zv1d, fvpsi, wf real(wp_), dimension(:), allocatable :: rv1d, zv1d, fvpsi, wf
integer :: ier, ixploc, info, i, j, ij integer :: ixploc, info, i, j, ij
character(256) :: msg ! for log messages formatting character(256) :: msg ! for log messages formatting
! compute array sizes ! compute array sizes
@ -468,10 +471,11 @@ contains
rmnm, rmxm, zmnm, zmxm, & rmnm, rmxm, zmnm, zmxm, &
psi_spline%knots_x, psi_spline%nknots_x, & psi_spline%knots_x, psi_spline%nknots_x, &
psi_spline%knots_y, psi_spline%nknots_y, & psi_spline%knots_y, psi_spline%nknots_y, &
psi_spline%coeffs, ier) psi_spline%coeffs, err)
! if failed, re-fit with an interpolating spline (zero tension) ! if failed, re-fit with an interpolating spline (zero tension)
if(ier == -1) then if(err == -1) then
err = 0
tension = 0 tension = 0
psi_spline%nknots_x=nr/4+4 psi_spline%nknots_x=nr/4+4
psi_spline%nknots_y=nz/4+4 psi_spline%nknots_y=nz/4+4
@ -479,7 +483,7 @@ contains
rmnm, rmxm, zmnm, zmxm, & rmnm, rmxm, zmnm, zmxm, &
psi_spline%knots_x, psi_spline%nknots_x, & psi_spline%knots_x, psi_spline%nknots_x, &
psi_spline%knots_y, psi_spline%nknots_y, & psi_spline%knots_y, psi_spline%nknots_y, &
psi_spline%coeffs, ier) psi_spline%coeffs, err)
end if end if
deallocate(rv1d, zv1d, wf, fvpsi) deallocate(rv1d, zv1d, wf, fvpsi)
! reset nrz to the total number of grid points for next allocations ! reset nrz to the total number of grid points for next allocations
@ -494,17 +498,23 @@ contains
! compute spline coefficients ! compute spline coefficients
call psi_spline%init(data%rv, data%zv, fvpsi, nr, nz, & call psi_spline%init(data%rv, data%zv, fvpsi, nr, nz, &
range=[rmnm, rmxm, zmnm, zmxm], & range=[rmnm, rmxm, zmnm, zmxm], &
tension=params%ssplps, err=ier) tension=params%ssplps, err=err)
! if failed, re-fit with an interpolating spline (zero tension) ! if failed, re-fit with an interpolating spline (zero tension)
if(ier == -1) then if(err == -1) then
call psi_spline%init(data%rv, data%zv, fvpsi, nr, nz, & call psi_spline%init(data%rv, data%zv, fvpsi, nr, nz, &
range=[rmnm, rmxm, zmnm, zmxm], & range=[rmnm, rmxm, zmnm, zmxm], &
tension=zero) tension=zero)
err = 0
end if end if
deallocate(fvpsi) deallocate(fvpsi)
end if end if
if (err /= 0) then
err = 2
return
end if
! compute spline coefficients for ψ(R,z) partial derivatives ! compute spline coefficients for ψ(R,z) partial derivatives
call psi_spline%init_deriv(nr, nz, 1, 0) ! ψ/R call psi_spline%init_deriv(nr, nz, 1, 0) ! ψ/R
call psi_spline%init_deriv(nr, nz, 0, 1) ! ψ/z call psi_spline%init_deriv(nr, nz, 0, 1) ! ψ/z

View File

@ -1,6 +1,6 @@
subroutine gray_jetto1beam(ijetto, mr, mz, r, z, psin, psia, rax, zax, & subroutine gray_jetto1beam(ijetto, mr, mz, r, z, psin, psia, rax, zax, &
nbnd, rbnd, zbnd, nrho, psrad, fpol, te, dne, zeff, qpsi, ibeam, & nbnd, rbnd, zbnd, nrho, psrad, fpol, te, dne, zeff, qpsi, ibeam, &
p0mw, alphain, betain, dpdv, jcd, pabs, icd, error) p0mw, alphain, betain, dpdv, jcd, pabs, icd, err)
use const_and_precisions, only: wp_ use const_and_precisions, only: wp_
use units, only: set_active_units, close_units use units, only: set_active_units, close_units
use gray_params, only: gray_parameters, gray_data, gray_results use gray_params, only: gray_parameters, gray_data, gray_results
@ -18,7 +18,7 @@ subroutine gray_jetto1beam(ijetto, mr, mz, r, z, psin, psia, rax, zax, &
real(wp_), dimension(nrho), intent(in) :: psrad, fpol, te, dne, zeff, qpsi real(wp_), dimension(nrho), intent(in) :: psrad, fpol, te, dne, zeff, qpsi
real(wp_), dimension(nrho), intent(out) :: dpdv, jcd real(wp_), dimension(nrho), intent(out) :: dpdv, jcd
real(wp_), intent(out) :: pabs, icd real(wp_), intent(out) :: pabs, icd
integer, intent(out) :: error integer, intent(out) :: err
! local variables ! local variables
type(gray_parameters) :: params type(gray_parameters) :: params
@ -43,7 +43,8 @@ subroutine gray_jetto1beam(ijetto, mr, mz, r, z, psin, psia, rax, zax, &
! Read parameters from external file ! Read parameters from external file
init_params: block init_params: block
use gray_params, only : read_gray_params, set_globals use gray_params, only : read_gray_params, set_globals
call read_gray_params('gray.data', params) call read_gray_params('gray.data', params, err)
if (err /= 0) return
! Override some parameters ! Override some parameters
params%misc%rwall = r(1) params%misc%rwall = r(1)
@ -93,7 +94,8 @@ subroutine gray_jetto1beam(ijetto, mr, mz, r, z, psin, psia, rax, zax, &
data%equilibrium%qpsi = qpsi data%equilibrium%qpsi = qpsi
! Compute splines ! Compute splines
call set_equil_spline(params%equilibrium, data%equilibrium) call set_equil_spline(params%equilibrium, data%equilibrium, err)
if (err /= 0) return
end block init_equilibrium end block init_equilibrium
! Set plasma kinetic profiles ! Set plasma kinetic profiles
@ -106,7 +108,8 @@ subroutine gray_jetto1beam(ijetto, mr, mz, r, z, psin, psia, rax, zax, &
data%profiles%zfc = zeff data%profiles%zfc = zeff
! Compute splines ! Compute splines
call set_profiles_spline(params%profiles, data%profiles) call set_profiles_spline(params%profiles, data%profiles, err)
if (err /= 0) return
end block init_profiles end block init_profiles
! Set wave launcher parameters ! Set wave launcher parameters
@ -119,11 +122,12 @@ subroutine gray_jetto1beam(ijetto, mr, mz, r, z, psin, psia, rax, zax, &
params%antenna%power = p0mw params%antenna%power = p0mw
! Read beam description file ! Read beam description file
call read_beam2(params%antenna, beamid=ibeam) call read_beam2(params%antenna, beamid=ibeam, err=err)
if (err /= 0) return
end block init_antenna end block init_antenna
! Call main subroutine for the ibeam-th beam ! Call main subroutine for the ibeam-th beam
call gray_main(params, data, res, error, rhout=sqrt(psrad)) call gray_main(params, data, res, err, rhout=sqrt(psrad))
! Free memory ! Free memory
free_memory: block free_memory: block

View File

@ -243,7 +243,7 @@ contains
end subroutine print_parameters end subroutine print_parameters
function update_parameter(params, name, value) result(error) function update_parameter(params, name, value) result(err)
! Updates the value of a parameter, addressed by a string ! Updates the value of a parameter, addressed by a string
! The return error is: ! The return error is:
! ERR_SUCCESS on success; ! ERR_SUCCESS on success;
@ -258,42 +258,41 @@ contains
! function arguments ! function arguments
type(gray_parameters), intent(inout) :: params type(gray_parameters), intent(inout) :: params
character(*), intent(in) :: name, value character(*), intent(in) :: name, value
integer(kind(ini_error)) :: error integer(kind(ini_error)) :: err
! hope for the best... ! hope for the best...
error = ERR_SUCCESS err = ERR_SUCCESS
select case (name) select case (name)
#include "gray_params.inc" #include "gray_params.inc"
case default case default
! unknown parameter ! unknown parameter
error = ERR_UNKNOWN err = ERR_UNKNOWN
end select end select
end function update_parameter end function update_parameter
subroutine read_gray_config(filename, params) subroutine read_gray_config(filename, params, err)
! Reads the GRAY parameters from the gray.ini configuration file ! Reads the GRAY parameters from the gray.ini configuration file
use ini_parser, only : parse_ini, property_handler, ini_error, ERR_SUCCESS use ini_parser, only : parse_ini, property_handler, ini_error
implicit none implicit none
! subroutine arguments ! subroutine arguments
character(len=*), intent(in) :: filename character(len=*), intent(in) :: filename
type(gray_parameters), intent(out) :: params type(gray_parameters), intent(out) :: params
integer(kind(ini_error)), intent(out) :: err
! local variables ! local variables
procedure(property_handler), pointer :: p procedure(property_handler), pointer :: p
integer :: error
p => ini_handler p => ini_handler
call parse_ini(filename, p, error) call parse_ini(filename, p, err)
if (error /= ERR_SUCCESS) call exit(1)
contains contains
function ini_handler(section, name, value) result(error) function ini_handler(section, name, value) result(err)
! This function handles a single INI property and updates ! This function handles a single INI property and updates
! the `params` structure ! the `params` structure
@ -301,15 +300,15 @@ contains
! function arguments ! function arguments
character(*), intent(in) :: section, name, value character(*), intent(in) :: section, name, value
integer(kind(ini_error)) :: error integer(kind(ini_error)) :: err
error = update_parameter(params, section // "." // name, value) err = update_parameter(params, section // "." // name, value)
end function ini_handler end function ini_handler
end subroutine read_gray_config end subroutine read_gray_config
subroutine read_gray_params(filename, params) subroutine read_gray_params(filename, params, err)
! Reads the GRAY parameters from the legacy gray_params.data file ! Reads the GRAY parameters from the legacy gray_params.data file
use utils, only : get_free_unit use utils, only : get_free_unit
use logger, only : log_error use logger, only : log_error
@ -317,11 +316,12 @@ contains
implicit none implicit none
! subrouting arguments ! subrouting arguments
character(len=*), intent(in) :: filename character(len=*), intent(in) :: filename
type(gray_parameters), intent(out) :: params type(gray_parameters), intent(out) :: params
integer, intent(out) :: err
! local variables ! local variables
integer :: u, err integer :: u
u = get_free_unit() u = get_free_unit()
@ -329,7 +329,7 @@ contains
if (err /= 0) then if (err /= 0) then
call log_error('opening gray_params file ('//filename//') failed!', & call log_error('opening gray_params file ('//filename//') failed!', &
mod='gray_params', proc='read_gray_params') mod='gray_params', proc='read_gray_params')
call exit(1) return
end if end if
read(u, *) params%raytracing%nrayr, params%raytracing%nrayth, & read(u, *) params%raytracing%nrayr, params%raytracing%nrayth, &

View File

@ -23,8 +23,8 @@ for set in $sets; do
for param in $(deref "$set"); do for param in $(deref "$set"); do
cat <<EOF cat <<EOF
case ('$set.$param') case ('$set.$param')
read (value, *, iostat=error) params%$set%$param read (value, *, iostat=err) params%$set%$param
if (error > 0) error = ERR_VALUE if (err > 0) err = ERR_VALUE
EOF EOF
done done
done done

View File

@ -40,12 +40,14 @@ program main
if (allocated(opts%config_file)) then if (allocated(opts%config_file)) then
! Use the newer gray.ini configuration file ! Use the newer gray.ini configuration file
call log_message(level=INFO, mod='main', msg='using GRAY INI config') call log_message(level=INFO, mod='main', msg='using GRAY INI config')
call read_gray_config(opts%config_file, params) call read_gray_config(opts%config_file, params, err)
if (err /= 0) call exit(1)
err = chdir(dirname(opts%config_file)) err = chdir(dirname(opts%config_file))
else else
! Use the legacy gray_params.data file ! Use the legacy gray_params.data file
call log_message(level=INFO, mod='main', msg='using GRAY params file') call log_message(level=INFO, mod='main', msg='using GRAY params file')
call read_gray_params(opts%params_file, params) call read_gray_params(opts%params_file, params, err)
if (err /= 0) call exit(1)
err = chdir(dirname(opts%params_file)) err = chdir(dirname(opts%params_file))
end if end if
@ -62,10 +64,16 @@ program main
! Read the input data and set the global variables ! Read the input data and set the global variables
! of the respective module. Note: order matters. ! of the respective module. Note: order matters.
call init_equilibrium(params, data) call init_antenna(params%antenna, err)
call init_antenna(params%antenna) if (err /= 0) call exit(err)
call init_equilibrium(params, data, err)
if (err /= 0) call exit(err)
call init_profiles(params%profiles, params%equilibrium%factb, & call init_profiles(params%profiles, params%equilibrium%factb, &
params%antenna%pos, data%profiles) params%antenna%pos, data%profiles, err)
if (err /= 0) call exit(err)
call init_misc(params, data) call init_misc(params, data)
! Change the current directory to output files there ! Change the current directory to output files there
@ -188,7 +196,7 @@ program main
contains contains
subroutine init_equilibrium(params, data) subroutine init_equilibrium(params, data, err)
! Reads the MHD equilibrium file (either in the G-EQDSK format ! Reads the MHD equilibrium file (either in the G-EQDSK format
! or an analytical description) and initialises the respective ! or an analytical description) and initialises the respective
! GRAY parameters and data. ! GRAY parameters and data.
@ -201,6 +209,7 @@ contains
! subroutine arguments ! subroutine arguments
type(gray_parameters), intent(inout) :: params type(gray_parameters), intent(inout) :: params
type(gray_data), intent(out) :: data type(gray_data), intent(out) :: data
integer, intent(out) :: err
if (params%equilibrium%iequil < 2) then if (params%equilibrium%iequil < 2) then
! Analytical equilibrium ! Analytical equilibrium
@ -208,7 +217,8 @@ contains
mod='main', proc='init_equilibrium') mod='main', proc='init_equilibrium')
call read_equil_an(params%equilibrium%filenm, & call read_equil_an(params%equilibrium%filenm, &
params%raytracing%ipass, & params%raytracing%ipass, &
data%equilibrium) data%equilibrium, err)
if (err /= 0) return
! Set psia sign to give the correct sign to Iphi ! Set psia sign to give the correct sign to Iphi
! (COCOS=3: psia<0 for Iphi>0) ! (COCOS=3: psia<0 for Iphi>0)
@ -218,7 +228,8 @@ contains
! Numerical equilibrium ! Numerical equilibrium
call log_debug('loading G-EQDK file', & call log_debug('loading G-EQDK file', &
mod='main', proc='init_equilibrium') mod='main', proc='init_equilibrium')
call read_eqdsk(params%equilibrium, data%equilibrium) call read_eqdsk(params%equilibrium, data%equilibrium, err)
if (err /= 0) return
call change_cocos(data%equilibrium, params%equilibrium%icocos, 3) call change_cocos(data%equilibrium, params%equilibrium%icocos, 3)
end if end if
@ -230,7 +241,8 @@ contains
call set_equil_an(data%equilibrium) call set_equil_an(data%equilibrium)
else else
call log_debug('computing splines...', mod='main', proc='init_equilibrium') call log_debug('computing splines...', mod='main', proc='init_equilibrium')
call set_equil_spline(params%equilibrium, data%equilibrium) call set_equil_spline(params%equilibrium, data%equilibrium, err)
if (err /= 0) return
call log_debug('splines computed', mod='main', proc='init_equilibrium') call log_debug('splines computed', mod='main', proc='init_equilibrium')
end if end if
end subroutine init_equilibrium end subroutine init_equilibrium
@ -257,7 +269,7 @@ contains
end subroutine deinit_equilibrium end subroutine deinit_equilibrium
subroutine init_profiles(params, factb, launch_pos, data) subroutine init_profiles(params, factb, launch_pos, data, err)
! Reads the plasma kinetic profiles file (containing the elecron ! Reads the plasma kinetic profiles file (containing the elecron
! temperature, density and plasma effective charge) and initialises ! temperature, density and plasma effective charge) and initialises
! the respective GRAY data structure. ! the respective GRAY data structure.
@ -275,17 +287,20 @@ contains
real(wp_), intent(in) :: factb real(wp_), intent(in) :: factb
real(wp_), intent(in) :: launch_pos(3) real(wp_), intent(in) :: launch_pos(3)
type(profiles_data), intent(out) :: data type(profiles_data), intent(out) :: data
integer, intent(out) :: err
if (params%iprof == 0) then if (params%iprof == 0) then
! Analytical profiles ! Analytical profiles
call log_debug('loading analytical file', & call log_debug('loading analytical file', &
mod='main', proc='init_profiles') mod='main', proc='init_profiles')
call read_profiles_an(params%filenm, data) call read_profiles_an(params%filenm, data, err)
if (err /= 0) return
else else
! Numerical profiles ! Numerical profiles
call log_debug('loading numerical file', & call log_debug('loading numerical file', &
mod='main', proc='init_profiles') mod='main', proc='init_profiles')
call read_profiles(params%filenm, data) call read_profiles(params%filenm, data, err)
if (err /= 0) return
! Convert psrad to ψ ! Convert psrad to ψ
select case (params%irho) select case (params%irho)
@ -307,7 +322,7 @@ contains
else else
! Numerical profiles ! Numerical profiles
call log_debug('computing splines...', mod='main', proc='init_profiles') call log_debug('computing splines...', mod='main', proc='init_profiles')
call set_profiles_spline(params, data, launch_pos) call set_profiles_spline(params, data, err, launch_pos)
call log_debug('splines computed', mod='main', proc='init_profiles') call log_debug('splines computed', mod='main', proc='init_profiles')
end if end if
end subroutine init_profiles end subroutine init_profiles
@ -332,7 +347,7 @@ contains
end subroutine deinit_profiles end subroutine deinit_profiles
subroutine init_antenna(params) subroutine init_antenna(params, err)
! Reads the wave launcher file (containing the wave frequency, launcher ! Reads the wave launcher file (containing the wave frequency, launcher
! position, direction and beam description) and initialises the respective ! position, direction and beam description) and initialises the respective
! GRAY parameters. ! GRAY parameters.
@ -343,6 +358,7 @@ contains
! subroutine arguments ! subroutine arguments
type(antenna_parameters), intent(inout) :: params type(antenna_parameters), intent(inout) :: params
integer, intent(out) :: err
! Note: α, β are loaded from gray_params.data ! Note: α, β are loaded from gray_params.data
select case (params%ibeam) select case (params%ibeam)
@ -350,14 +366,14 @@ contains
! 2 degrees of freedom ! 2 degrees of freedom
! w(z, α, β), 1/R(z, α, β) ! w(z, α, β), 1/R(z, α, β)
! FIXME: 1st beam is always selected, iox read from table ! FIXME: 1st beam is always selected, iox read from table
call read_beam2(params, beamid=1) call read_beam2(params, beamid=1, err=err)
case (1) case (1)
! 1 degree of freedom ! 1 degree of freedom
! w(z, α), 1/R(z, α) ! w(z, α), 1/R(z, α)
call read_beam1(params) call read_beam1(params, err)
case default case default
! fixed w(z), 1/R(z) ! fixed w(z), 1/R(z)
call read_beam0(params) call read_beam0(params, err)
end select end select
end subroutine init_antenna end subroutine init_antenna