mms_diffusion_slab_m.f90 Source File


Contents


Source Code

module mms_diffusion_slab_m
    !! MMS Solutions and Sources for diffusion model in slab geometry
    use precision_grillix_m, only : GP
    use params_diffusion_m, only : diffcoeff_perp, diffcoeff_par

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

    public :: mms_sol_diffusion_slab
    public :: mms_source_diffusion_slab

    real(GP), private :: amp = 1.0_GP
    !! Amplitude of MMS solution
    real(GP), private :: kx = 1.0_GP
    !! x-mode number of MMS Solution
    real(GP), private :: ky = 2.0_GP
    !! y-mode number of MMS Solution
    real(GP), private :: phy = 0.2_GP
    !! Poloidal phase shift of MMS Solution
    real(GP), private :: kz = 1.0_GP
    !! Axial mode number of MMS Solution
    real(GP), private :: phz = 1.1_GP
    !! Axial phase shift of MMS Solution
    real(GP), private :: omega = 0.5_GP
    !! Angular time frequency of MMS Solution
    real(GP), private :: pht = 0.4_GP
    !! Temporal phase shift of MMS Solution

contains
    
    real(GP) function mms_sol_diffusion_slab(equi, x, y, z, t)
        !! MMS solution
        class(equilibrium_t) :: equi
        !! Equilibrium
        real(GP), intent(in) :: x
        !! x-coordinate 
        real(GP), intent(in) :: y
        !! y-coordinate
        real(GP), intent(in) :: z
        !! Axial coordinate z (represented by toroidal coordinate)
        real(GP), intent(in) :: t
        !! Time

        real(GP) :: res, xmin, xmax, ymin, ymax
        real(GP), parameter :: twopi = TWO_PI
        include 'mms_sol_slab_head.txt'

        xmin = equi%xmin
        xmax = equi%xmax
        ymin = equi%ymin
        ymax = equi%ymax
                    
        include 'mms_sol_slab.txt'
        mms_sol_diffusion_slab = res

    end function

    real(GP) function mms_source_diffusion_slab(equi, x, y, z, t)
        !! Source corresponding to MMS solution
        class(equilibrium_t), intent(inout) :: equi
        !! Equilibrium
        real(GP), intent(in) :: x
        !! x-coordinate 
        real(GP), intent(in) :: y
        !! y-coordinate
        real(GP), intent(in) :: z
        !! Axial coordinate z (represented by toroidal coordinate)
        real(GP), intent(in) :: t
        !! Time
        
        real(GP), parameter :: twopi = TWO_PI       
        real(GP) :: res, xmin, xmax, ymin, ymax
        real(GP) ::  chiperp, chipar
        include 'mms_source_slab_head.txt'

        xmin = equi%xmin
        xmax = equi%xmax
        ymin = equi%ymin
        ymax = equi%ymax

        chiperp = diffcoeff_perp
        chipar = diffcoeff_par
                
        include 'mms_source_slab.txt'
        mms_source_diffusion_slab = res

    end function
  
end module