helpers_mesh_sample_m.f90 Source File


Contents


Source Code

module helpers_mesh_sample_m
    !! Creates test mehses that can be used further on for unit tests
    !! Mesh itself is tested within fci-core
    use MPI
    use precision_grillix_m, only : GP 
    
    ! From FCI-CORE
    use equilibrium_m,          only : equilibrium_t
    use equilibrium_factory_m,  only : create_equilibrium, CIRCULAR, CERFONS
    use mesh_cart_m,            only : mesh_cart_t
    use multigrid_m,            only : multigrid_t
    implicit none

    public :: mesh_sample_circular
    public :: mesh_sample_cerfons

contains

    subroutine mesh_sample_circular(equi, phi, mesh, multigrid)
        ! Creates a testmesh for circular geometry
        class(equilibrium_t), allocatable, intent(out) :: equi
        !! Equilibrium
        real(GP), intent(in) :: phi
        !! Toroidal angle
        type(mesh_cart_t), intent(out) :: mesh
        !! Mesh within poloidal plane
        type(multigrid_t), intent(out), optional :: multigrid
        !! Multigrid 

        real(GP) :: spacing_f

        spacing_f = 8.0E-3_GP
        
        call create_equilibrium(equi, CIRCULAR)

        if (present(multigrid)) then
            call multigrid%init(equi, phi, 2, spacing_f, 2, 2, mesh, [8,8], .false., 0)
        else
            call mesh%init(equi, phi, 1, spacing_f, &
                           size_neighbor=2, size_ghost_layer=2, &
                           extend_beyond_wall=.false., dbgout=0)
            call mesh%reorder(8, dbgout=0)
        endif       
        
    end


    subroutine mesh_sample_cerfons(equi, phi, mesh, multigrid)
        ! Creates a testmesh for circular geometry
        class(equilibrium_t), allocatable, intent(out) :: equi
        !! Equilibrium
        real(GP), intent(in) :: phi
        !! Toroidal angle
        type(mesh_cart_t), intent(out) :: mesh
        !! Mesh within poloidal plane  
        type(multigrid_t), intent(out), optional :: multigrid
        !! Multigrid 
        
        real(GP) :: spacing_f
        
        call create_equilibrium(equi, CERFONS)

        spacing_f = 1.0E-2_GP
        
        if (present(multigrid)) then
            call multigrid%init(equi, phi, 2, spacing_f, 2, 2, mesh, [8,8], .false., 0)
        else
            call mesh%init(equi, phi, 1, spacing_f, &
                        size_neighbor=2, size_ghost_layer=2, &
                        extend_beyond_wall=.false., dbgout=0)
            call mesh%reorder(8, dbgout=0)
        endif
    end

end module