mms_template_slab_s.f90 Source File


Contents


Source Code

submodule (mms_template_m) mms_template_slab_s
    !! MMS Solutions and Sources for template model in slab geometry
    use precision_grillix_m, only : GP
    use error_handling_grillix_m, only: handle_error
    use status_codes_grillix_m, only : GRILLIX_ERR_OTHER
    use params_template_m, only : path_mms, perptransp_on, partransp_on, pardiff_on, &
        dperp_coeff_dens, dperp_coeff_parmom, dpar_coeff_dens, dpar_coeff_parmom, &
        dperp_coeff_pion, dpar_coeff_pion, ion_heat_cond_coeff
    use static_data_m, only : equi

    use constants_m, only : TWO_PI
    use equilibrium_m, only : equilibrium_t
    use slab_equilibrium_m, only : slab_t
    implicit none

    real(GP) :: bckg_n = 1.1_GP, bckg_pi = 1.2_GP
    !! Constant background of MMS solution
    real(GP) :: amp_n = 1.0_GP, amp_pm = 0.8_GP, amp_pi = 0.5_GP
    !! Amplitude of MMS solution
    real(GP) :: kx_n = 1.0_GP, kx_pm = 2.0_GP, kx_pi = 3.0_GP
    !! x-mode number of MMS Solution
    real(GP) :: ky_n = 2.0_GP, ky_pm = 1.0_GP, ky_pi = 2.0_GP
    !! y-mode number of MMS Solution
    real(GP) :: phy_n = 0.2_GP, phy_pm = 0.4_GP, phy_pi = 0.6_GP
    !! Poloidal phase shift of MMS Solution
    real(GP) :: kz_n = 1.0_GP, kz_pm = 1.0_GP, kz_pi = 1.0_GP
    !! Axial mode number of MMS Solution
    real(GP) :: phz_n = 1.1_GP, phz_pm = 0.2_GP, phz_pi = 0.3_GP
    !! Axial phase shift of MMS Solution
    real(GP) :: omega_n = 4.0_GP, omega_pm = 5.0_GP, omega_pi = 3.0_GP
    !! Angular time frequency of MMS Solution
    real(GP) :: pht_n = 0.4_GP, pht_pm = 0.1_GP, pht_pi = 0.3_GP
    !! Temporal phase shift of MMS Solution

    namelist / mms_slab_template / bckg_n, bckg_pi, &
        amp_n, kx_n, ky_n, phy_n, kz_n, phz_n, omega_n, pht_n, &
        amp_pm, kx_pm, ky_pm, phy_pm, kz_pm, phz_pm, omega_pm, pht_pm, &
        amp_pi, kx_pi, ky_pi, phy_pi, kz_pi, phz_pi, omega_pi, pht_pi

    real(GP) :: switch_perptransp = 1.0_GP
    !! Switch for parallel transport model
    real(GP) :: switch_partransp = 1.0_GP
    !! Switch for parallel transport model
    real(GP) :: switch_pardiff = 1.0_GP
    !! Switch for parallel transport model

contains

    module subroutine init_mms_slab_template()
        !! Initializes MMS parameters for slab geometry
        integer :: io_unit, io_error
        character(len=256) :: io_errmsg

        ! Read slab MMS parameters
        open(newunit=io_unit, file=path_mms, status='old', action='read', iostat=io_error, iomsg=io_errmsg)
        if (io_error /= 0) call handle_error(io_errmsg, GRILLIX_ERR_NAMELIST, __LINE__, __FILE__)
        read(io_unit, nml=mms_slab_template, iostat=io_error, iomsg=io_errmsg)
        if (io_error /= 0) call handle_error(io_errmsg, GRILLIX_ERR_NAMELIST, __LINE__, __FILE__)
        close(io_unit)

        ! Set procedure pointers
        mms_sol_dens_template => mms_sol_dens_template_slab
        mms_source_dens_template => mms_source_dens_template_slab

        mms_sol_parmom_template => mms_sol_parmom_template_slab
        mms_source_parmom_template => mms_source_parmom_template_slab

        mms_sol_pion_template => mms_sol_pion_template_slab
        mms_source_pion_template => mms_source_pion_template_slab

        ! Set model switches
        switch_perptransp = merge(1.0_GP, 0.0_GP, perptransp_on)
        switch_partransp = merge(1.0_GP, 0.0_GP, partransp_on)
        switch_pardiff = merge(1.0_GP, 0.0_GP, pardiff_on)

    end subroutine

    module subroutine display_mms_params_template_slab()
        if (is_rank_info_writer) then
            write(get_stdout(), nml=mms_slab_template)
        endif

    end subroutine
    
    function mms_sol_dens_template_slab(x, y, z, t) result(res)
        !! MMS solution
        real(GP), intent(in) :: x
        !! x-coordinate 
        real(GP), intent(in) :: y
        !! y-coordinate
        real(GP), intent(in) :: z
        !! Axial coordinate phi (toroidal coordinate)
        real(GP), intent(in) :: t
        !! Time

        real(GP) :: xmin, xmax, ymin, ymax

        include 'generated/mms_sol_Dens_slab_head.inc'
        
        select type(equi)
            type is(slab_t)
                xmin = equi%xmin
                xmax = equi%xmax
                ymin = equi%ymin
                ymax = equi%ymax
            class default
                call handle_error('Equilibrium is not slab', GRILLIX_ERR_OTHER, __LINE__, __FILE__)
        end select

        include 'generated/mms_sol_Dens_slab.inc'

    end function

    function mms_sol_parmom_template_slab(x, y, z, t) result(res)
        !! MMS solution
        real(GP), intent(in) :: x
        !! x-coordinate 
        real(GP), intent(in) :: y
        !! y-coordinate
        real(GP), intent(in) :: z
        !! Axial coordinate phi (toroidal coordinate)
        real(GP), intent(in) :: t
        !! Time

        real(GP) :: xmin, xmax, ymin, ymax

        include 'generated/mms_sol_Parmom_slab_head.inc'
        
        select type(equi)
            type is(slab_t)
                xmin = equi%xmin
                xmax = equi%xmax
                ymin = equi%ymin
                ymax = equi%ymax
            class default
                call handle_error('Equilibrium is not slab', GRILLIX_ERR_OTHER, __LINE__, __FILE__)
        end select

        include 'generated/mms_sol_Parmom_slab.inc'

    end function

    function mms_sol_pion_template_slab(x, y, z, t) result(res)
        !! MMS solution
        real(GP), intent(in) :: x
        !! x-coordinate 
        real(GP), intent(in) :: y
        !! y-coordinate
        real(GP), intent(in) :: z
        !! Axial coordinate phi (toroidal coordinate)
        real(GP), intent(in) :: t
        !! Time

        real(GP) :: xmin, xmax, ymin, ymax

        include 'generated/mms_sol_Pion_slab_head.inc'
        
        select type(equi)
            type is(slab_t)
                xmin = equi%xmin
                xmax = equi%xmax
                ymin = equi%ymin
                ymax = equi%ymax
            class default
                call handle_error('Equilibrium is not slab', GRILLIX_ERR_OTHER, __LINE__, __FILE__)
        end select

        include 'generated/mms_sol_Pion_slab.inc'

    end function

    function mms_source_dens_template_slab(x, y, z, t) result(res)
        !! Source corresponding to MMS solution
        real(GP), intent(in) :: x
        !! x-coordinate 
        real(GP), intent(in) :: y
        !! y-coordinate
        real(GP), intent(in) :: z
        !! Axial coordinate phi (toroidal coordinate)
        real(GP), intent(in) :: t
        !! Time
        
        real(GP) :: xmin, xmax, ymin, ymax
        
        include 'generated/mms_source_Dens_slab_head.inc'

        select type(equi)
            type is(slab_t)
                xmin = equi%xmin
                xmax = equi%xmax
                ymin = equi%ymin
                ymax = equi%ymax
            class default
                call handle_error('Equilibrium is not slab', GRILLIX_ERR_OTHER, __LINE__, __FILE__)
        end select

        include 'generated/mms_source_Dens_slab.inc'

    end function

    function mms_source_parmom_template_slab(x, y, z, t) result(res)
        !! Source corresponding to MMS solution
        real(GP), intent(in) :: x
        !! x-coordinate 
        real(GP), intent(in) :: y
        !! y-coordinate
        real(GP), intent(in) :: z
        !! Axial coordinate phi (toroidal coordinate)
        real(GP), intent(in) :: t
        !! Time
        
        real(GP) :: xmin, xmax, ymin, ymax
        
        include 'generated/mms_source_Parmom_slab_head.inc'

        select type(equi)
            type is(slab_t)
                xmin = equi%xmin
                xmax = equi%xmax
                ymin = equi%ymin
                ymax = equi%ymax
            class default
                call handle_error('Equilibrium is not slab', GRILLIX_ERR_OTHER, __LINE__, __FILE__)
        end select

        include 'generated/mms_source_Parmom_slab.inc'

    end function

    function mms_source_pion_template_slab(x, y, z, t) result(res)
        !! Source corresponding to MMS solution
        real(GP), intent(in) :: x
        !! x-coordinate 
        real(GP), intent(in) :: y
        !! y-coordinate
        real(GP), intent(in) :: z
        !! Axial coordinate phi (toroidal coordinate)
        real(GP), intent(in) :: t
        !! Time
        
        real(GP) :: xmin, xmax, ymin, ymax
        
        include 'generated/mms_source_Pion_slab_head.inc'

        select type(equi)
            type is(slab_t)
                xmin = equi%xmin
                xmax = equi%xmax
                ymin = equi%ymin
                ymax = equi%ymax
            class default
                call handle_error('Equilibrium is not slab', GRILLIX_ERR_OTHER, __LINE__, __FILE__)
        end select

        include 'generated/mms_source_Pion_slab.inc'

    end function

end submodule