params_equi_salpha_m.f90 Source File


Source Code

module params_equi_salpha_m
    !! Module for parameter reading for the salpha equilibrium type
    use, intrinsic :: iso_fortran_env, only: IOSTAT_END, IOSTAT_EOR
    use screen_io_m, only: get_stdout
    use precision_m, only: FP
    use error_handling_m, only: handle_error, error_info_t
    use status_codes_m, only: PARALLAX_ERR_PARAMETERS

    implicit none
    private

    real(FP) :: rhomin = 0.1_FP
    !! Minimum rho in the simulation domain
    real(FP) :: rhomax = 0.7_FP
    !! Maximum rho in the simulation domain
    real(FP) :: minor_r = 0.365_FP
    !! Minor radius of the tokamak normalized to L_ref
    real(FP) :: B_ref = 1.0_FP
    !! Reference magnetic field
    real(FP) :: L_ref = 1.0_FP
    !! Reference length scale
    real(FP) :: q_ref = 0.854_FP
    !! Safety factor value at the core
    real(FP) :: shear = 2.184_FP
    !! Magnetic shear
    real(FP) :: dtheta_limiter = 0.4_FP
    !! Poloidal width of the limiter. If smaller equal than zero, no limiter
    !! is used.
    real(FP) :: rho_limiter = 1.0_FP
    !! Radial position (start) of the limiter
    real(FP) :: theta_limiter = 3.1415_FP
    !! Poloidal position of the limiter

    public :: read_params_salpha
    public :: write_params_salpha

    public :: get_salpha_rhomin
    public :: get_salpha_rhomax
    public :: get_salpha_minor_r
    public :: get_salpha_B_ref
    public :: get_salpha_L_ref
    public :: get_salpha_q_ref
    public :: get_salpha_shear
    public :: get_salpha_dtheta_limiter
    public :: get_salpha_rho_limiter
    public :: get_salpha_theta_limiter

    namelist / equi_salpha_params / &
        rhomin, rhomax, q_ref, shear, minor_r, L_ref, B_ref, &
        dtheta_limiter, rho_limiter, theta_limiter

contains

    pure real(FP) function get_salpha_rhomin()
        get_salpha_rhomin = rhomin
    end function

    pure real(FP) function get_salpha_rhomax()
        get_salpha_rhomax = rhomax
    end function

    pure real(FP) function get_salpha_minor_r()
        get_salpha_minor_r = minor_r
    end function

    pure real(FP) function get_salpha_B_ref()
        get_salpha_B_ref = B_ref
    end function

    pure real(FP) function get_salpha_L_ref()
        get_salpha_L_ref = L_ref
    end function

    pure real(FP) function get_salpha_q_ref()
        get_salpha_q_ref = q_ref
    end function

    pure real(FP) function get_salpha_shear()
        get_salpha_shear = shear
    end function

    pure real(FP) function get_salpha_dtheta_limiter()
        get_salpha_dtheta_limiter = dtheta_limiter
    end function

    pure real(FP) function get_salpha_rho_limiter()
        get_salpha_rho_limiter = rho_limiter
    end function

    pure real(FP) function get_salpha_theta_limiter()
        get_salpha_theta_limiter = theta_limiter
    end function

    subroutine read_params_salpha(filename)
        !! Reads the equi salpha namelist from the given filename
        character(len=*), intent(in) :: filename
        !! Filename to read from

        integer :: funit, io_error

        open(newunit=funit, file=filename, status='old', action='read', &
             iostat=io_error)
        if (io_error /= 0) then
            call handle_error("Opening parameter file "//filename//" failed!", &
                              PARALLAX_ERR_PARAMETERS, __LINE__, __FILE__, &
                              additional_info=&
                                error_info_t("Iostat:", [io_error]))
        end if

        read(funit, nml=equi_salpha_params, iostat=io_error)
        ! END and EOR are returned if namelist can not be found. We treat
        ! these cases as non fatal and write a notification.
        if (io_error == IOSTAT_END .or. io_error == IOSTAT_EOR) then
            write(get_stdout(), *) &
                "Info: Namelist equi_salpha_params not present."
        elseif (io_error /= 0) then
            call handle_error("Reading salpha equi parameters failed!", &
                              PARALLAX_ERR_PARAMETERS, __LINE__, __FILE__, &
                              additional_info=&
                                error_info_t("Iostat:", [io_error]))
        end if

        close(funit)
    end subroutine

    subroutine write_params_salpha(filename)
        !! Writes the equi salpha namelist to the given filename
        character(len=*), intent(in) :: filename
        !! Filename to read from

        integer :: funit, io_error

        open(newunit=funit, file=filename, action='write', position='append', &
             iostat=io_error)
        if (io_error /= 0) then
            call handle_error("Opening parameter file "//filename//" failed!", &
                              PARALLAX_ERR_PARAMETERS, __LINE__, __FILE__, &
                              additional_info=&
                                error_info_t("Iostat:", [io_error]))
        end if

        write(funit, nml=equi_salpha_params, iostat=io_error)
        if (io_error /= 0) then
            call handle_error("Writing salpha equi parameters failed!", &
                              PARALLAX_ERR_PARAMETERS, __LINE__, __FILE__, &
                              additional_info=&
                                error_info_t("Iostat:", [io_error]))
        end if

        close(funit)
    end subroutine

end module