params_equi_slab_m.f90 Source File


This file depends on

sourcefile~~params_equi_slab_m.f90~~EfferentGraph sourcefile~params_equi_slab_m.f90 params_equi_slab_m.f90 sourcefile~error_handling_m.f90 error_handling_m.f90 sourcefile~params_equi_slab_m.f90->sourcefile~error_handling_m.f90 sourcefile~precision_m.f90 precision_m.f90 sourcefile~params_equi_slab_m.f90->sourcefile~precision_m.f90 sourcefile~screen_io_m.f90 screen_io_m.f90 sourcefile~params_equi_slab_m.f90->sourcefile~screen_io_m.f90 sourcefile~status_codes_m.f90 status_codes_m.f90 sourcefile~params_equi_slab_m.f90->sourcefile~status_codes_m.f90 sourcefile~error_handling_m.f90->sourcefile~precision_m.f90 sourcefile~error_handling_m.f90->sourcefile~screen_io_m.f90 sourcefile~error_handling_m.f90->sourcefile~status_codes_m.f90 sourcefile~comm_handling_m.f90 comm_handling_m.f90 sourcefile~error_handling_m.f90->sourcefile~comm_handling_m.f90 sourcefile~screen_io_m.f90->sourcefile~precision_m.f90

Files dependent on this one

sourcefile~~params_equi_slab_m.f90~~AfferentGraph sourcefile~params_equi_slab_m.f90 params_equi_slab_m.f90 sourcefile~slab_equilibrium.f90 slab_equilibrium.f90 sourcefile~slab_equilibrium.f90->sourcefile~params_equi_slab_m.f90 sourcefile~coords_polar_m.f90 coords_polar_m.f90 sourcefile~coords_polar_m.f90->sourcefile~slab_equilibrium.f90 sourcefile~equilibrium_factory_m.f90 equilibrium_factory_m.f90 sourcefile~equilibrium_factory_m.f90->sourcefile~slab_equilibrium.f90 sourcefile~mesh_cart_m.f90 mesh_cart_m.f90 sourcefile~mesh_cart_m.f90->sourcefile~slab_equilibrium.f90 sourcefile~polar_grid_m.f90 polar_grid_m.f90 sourcefile~polar_grid_m.f90->sourcefile~slab_equilibrium.f90 sourcefile~polar_grid_m.f90->sourcefile~coords_polar_m.f90 sourcefile~testfunctions_m.f90 testfunctions_m.f90 sourcefile~testfunctions_m.f90->sourcefile~slab_equilibrium.f90 sourcefile~benchmark_helmholtz_solvers.f90 benchmark_helmholtz_solvers.f90 sourcefile~benchmark_helmholtz_solvers.f90->sourcefile~equilibrium_factory_m.f90 sourcefile~benchmark_helmholtz_solvers.f90->sourcefile~mesh_cart_m.f90 sourcefile~benchmark_helmholtz_solvers.f90->sourcefile~testfunctions_m.f90 sourcefile~helmholtz_netcdfio_m.f90 helmholtz_netcdfio_m.f90 sourcefile~benchmark_helmholtz_solvers.f90->sourcefile~helmholtz_netcdfio_m.f90 sourcefile~multigrid_m.f90 multigrid_m.f90 sourcefile~benchmark_helmholtz_solvers.f90->sourcefile~multigrid_m.f90 sourcefile~helmholtz_solver_factory_m.f90 helmholtz_solver_factory_m.f90 sourcefile~benchmark_helmholtz_solvers.f90->sourcefile~helmholtz_solver_factory_m.f90 sourcefile~boundaries_perp.f90 boundaries_perp.f90 sourcefile~boundaries_perp.f90->sourcefile~mesh_cart_m.f90 sourcefile~diagnose_poincare.f90 diagnose_poincare.f90 sourcefile~diagnose_poincare.f90->sourcefile~coords_polar_m.f90 sourcefile~diagnose_poincare.f90->sourcefile~equilibrium_factory_m.f90 sourcefile~helmholtz_netcdfio_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~helmholtz_operator_m.f90 helmholtz_operator_m.f90 sourcefile~helmholtz_netcdfio_m.f90->sourcefile~helmholtz_operator_m.f90 sourcefile~helmholtz_operator_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~helmholtz_operator_m.f90->sourcefile~boundaries_perp.f90 sourcefile~immersed_factory_m.f90 immersed_factory_m.f90 sourcefile~immersed_factory_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~immersed_m.f90 immersed_m.f90 sourcefile~immersed_factory_m.f90->sourcefile~immersed_m.f90 sourcefile~immersed_rho_m.f90 immersed_rho_m.f90 sourcefile~immersed_factory_m.f90->sourcefile~immersed_rho_m.f90 sourcefile~immersed_trace_m.f90 immersed_trace_m.f90 sourcefile~immersed_factory_m.f90->sourcefile~immersed_trace_m.f90 sourcefile~immersed_vessel_m.f90 immersed_vessel_m.f90 sourcefile~immersed_factory_m.f90->sourcefile~immersed_vessel_m.f90 sourcefile~immersed_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~immersed_rho_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~immersed_rho_m.f90->sourcefile~immersed_m.f90 sourcefile~immersed_trace_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~immersed_trace_m.f90->sourcefile~immersed_m.f90 sourcefile~immersed_vessel_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~immersed_vessel_m.f90->sourcefile~immersed_m.f90 sourcefile~map_factory_m.f90 map_factory_m.f90 sourcefile~map_factory_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~mesh_cart_build_s.f90 mesh_cart_build_s.f90 sourcefile~mesh_cart_build_s.f90->sourcefile~mesh_cart_m.f90 sourcefile~mesh_cart_communicate_s.f90 mesh_cart_communicate_s.f90 sourcefile~mesh_cart_communicate_s.f90->sourcefile~mesh_cart_m.f90 sourcefile~mesh_cart_netcdfio_s.f90 mesh_cart_netcdfio_s.f90 sourcefile~mesh_cart_netcdfio_s.f90->sourcefile~mesh_cart_m.f90 sourcefile~mesh_cart_reorder_s.f90 mesh_cart_reorder_s.f90 sourcefile~mesh_cart_reorder_s.f90->sourcefile~mesh_cart_m.f90 sourcefile~mesh_cart_s.f90 mesh_cart_s.f90 sourcefile~mesh_cart_s.f90->sourcefile~mesh_cart_m.f90 sourcefile~multigrid_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~multigrid_m.f90->sourcefile~boundaries_perp.f90 sourcefile~multigrid_solver_m.f90 multigrid_solver_m.f90 sourcefile~multigrid_solver_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~multigrid_solver_m.f90->sourcefile~boundaries_perp.f90 sourcefile~multigrid_solver_m.f90->sourcefile~helmholtz_operator_m.f90 sourcefile~multigrid_solver_m.f90->sourcefile~multigrid_m.f90 sourcefile~splitting_m.f90 splitting_m.f90 sourcefile~multigrid_solver_m.f90->sourcefile~splitting_m.f90 sourcefile~parbnd_taylor_m.f90 parbnd_taylor_m.f90 sourcefile~parbnd_taylor_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~polar_grid_s.f90 polar_grid_s.f90 sourcefile~polar_grid_s.f90->sourcefile~polar_grid_m.f90 sourcefile~polar_map_factory_m.f90 polar_map_factory_m.f90 sourcefile~polar_map_factory_m.f90->sourcefile~coords_polar_m.f90 sourcefile~polar_map_factory_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~polar_map_factory_m.f90->sourcefile~polar_grid_m.f90 sourcefile~safety_factor_m.f90 safety_factor_m.f90 sourcefile~safety_factor_m.f90->sourcefile~coords_polar_m.f90 sourcefile~splitting_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~test_diffusion.f90 test_diffusion.f90 sourcefile~test_diffusion.f90->sourcefile~equilibrium_factory_m.f90 sourcefile~test_diffusion.f90->sourcefile~mesh_cart_m.f90 sourcefile~test_diffusion.f90->sourcefile~map_factory_m.f90 sourcefile~vis_vtk3d_m.f90 vis_vtk3d_m.f90 sourcefile~test_diffusion.f90->sourcefile~vis_vtk3d_m.f90 sourcefile~vis_vtk3d_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~zonal_averages_factory_m.f90 zonal_averages_factory_m.f90 sourcefile~zonal_averages_factory_m.f90->sourcefile~coords_polar_m.f90 sourcefile~zonal_averages_factory_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~zonal_averages_factory_m.f90->sourcefile~polar_grid_m.f90 sourcefile~helmholtz_solver_factory_m.f90->sourcefile~multigrid_m.f90 sourcefile~helmholtz_solver_factory_m.f90->sourcefile~splitting_m.f90 sourcefile~immersed_netcdf_s.f90 immersed_netcdf_s.f90 sourcefile~immersed_netcdf_s.f90->sourcefile~immersed_m.f90 sourcefile~map_factory_s.f90 map_factory_s.f90 sourcefile~map_factory_s.f90->sourcefile~map_factory_m.f90 sourcefile~multigrid_s.f90 multigrid_s.f90 sourcefile~multigrid_s.f90->sourcefile~multigrid_m.f90 sourcefile~multigrid_solver_s.f90 multigrid_solver_s.f90 sourcefile~multigrid_solver_s.f90->sourcefile~boundaries_perp.f90 sourcefile~multigrid_solver_s.f90->sourcefile~multigrid_solver_m.f90 sourcefile~parbnd_taylor_netcdf_s.f90 parbnd_taylor_netcdf_s.f90 sourcefile~parbnd_taylor_netcdf_s.f90->sourcefile~parbnd_taylor_m.f90 sourcefile~splitting_gauss_seidel_cpu_s.f90 splitting_gauss_seidel_cpu_s.f90 sourcefile~splitting_gauss_seidel_cpu_s.f90->sourcefile~splitting_m.f90 sourcefile~splitting_gauss_seidel_redblack_cpu_s.f90 splitting_gauss_seidel_redblack_cpu_s.f90 sourcefile~splitting_gauss_seidel_redblack_cpu_s.f90->sourcefile~splitting_m.f90 sourcefile~splitting_jacobi_cpu_s.f90 splitting_jacobi_cpu_s.f90 sourcefile~splitting_jacobi_cpu_s.f90->sourcefile~splitting_m.f90

Source Code

module params_equi_slab_m
    !! Module for parameter reading for the slab 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) :: boxsize = 64.0_FP
    !! Size of the domain in x and y
    logical :: sol = .true.
    !! Indicates if slab is open or closed field lines
    ! NOTE: We choose open field lines as the default to enable the
    !       testing of boundary conditions for open field lines with
    !       unit tests.
    logical :: yperiodic = .true.
    !! Indicates if box is periodic in y direction

    public :: read_params_slab
    public :: write_params_slab

    public :: get_slab_boxsize
    public :: get_slab_sol
    public :: get_slab_yperiodic

    namelist / equi_slab_params / boxsize, sol, yperiodic

contains

    pure real(FP) function get_slab_boxsize()
        get_slab_boxsize = boxsize
    end function

    pure logical function get_slab_sol()
        get_slab_sol = sol
    end function

    pure logical function get_slab_yperiodic()
        get_slab_yperiodic = yperiodic
    end function

    subroutine read_params_slab(filename)
        !! Reads the equi slab 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_slab_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_slab_params not present."
        elseif (io_error /= 0) then
            call handle_error("Reading slab 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_slab(filename)
        !! Writes the equi slab 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_slab_params, iostat=io_error)
        if (io_error /= 0) then
            call handle_error("Writing slab equi parameters failed!", &
                              PARALLAX_ERR_PARAMETERS, __LINE__, __FILE__, &
                              additional_info=&
                                error_info_t("Iostat:", [io_error]))
        end if

        close(funit)
    end subroutine

end module