splitting_m.f90 Source File


Source Code

module splitting_m
    !! Splitting methods
    use precision_m,              only : FP
    use comm_handling_m,          only : is_master
    use csrmat_m,                 only : csrmat_t, csr_times_vec
    use mesh_cart_m,              only : mesh_cart_t
    implicit none

    type, public, abstract :: splitting_t
        !! Splitting solver
        integer, private :: ndim
        !! Dimension of linear equation system
    contains
        procedure(create_abstract), deferred :: create
        procedure(apply_abstract),  deferred :: apply
    end type

    interface

        module subroutine create_abstract(self, mesh)
            !! Creates a splitting solver
            class(splitting_t), intent(inout) :: self
            !! Instance of the type
            type(mesh_cart_t), intent(in), target :: mesh
            !! Mesh
        end subroutine

        module subroutine apply_abstract(self, a, dinv, x, b)
            !! Applies a single iteration of a splitting solver
            class(splitting_t), intent(inout) :: self
            !! Instance of the type
            type(csrmat_t), intent(in) :: a
            !! Matrix
            real(FP), intent(in), dimension(self%ndim) :: dinv
            !! Inverse of the diagonal elements of the matrix "a"
            real(FP), intent(inout), dimension(self%ndim) :: x
            !! Quantity to be smoothed
            real(FP), intent(in), dimension(self%ndim) :: b
            !! Right hand side
        end subroutine

    end interface

    type, public, extends(splitting_t) :: splitting_jacobi_cpu_t
        !! Splitting solver with Jacobi splitting method on CPU
        real(FP), private :: omega
        !! Damping factor
        real(FP), allocatable, private, dimension(:) :: xn
        !! Work array, declared here to be shared(OMP) variable
    contains
        procedure :: create => create_jacobi_cpu
        procedure :: apply => apply_jacobi_cpu
        final     :: destructor_jacobi_cpu
    end type

    interface

        module subroutine create_jacobi_cpu(self, mesh)
            class(splitting_jacobi_cpu_t), intent(inout) :: self
            type(mesh_cart_t), intent(in), target :: mesh
        end subroutine

        module subroutine apply_jacobi_cpu(self, a, dinv, x, b)
            class(splitting_jacobi_cpu_t), intent(inout) :: self
            type(csrmat_t), intent(in) :: a
            real(FP), intent(in), dimension(self%ndim) :: dinv
            real(FP), intent(inout), dimension(self%ndim) :: x
            real(FP), intent(in), dimension(self%ndim) :: b
        end subroutine

        module subroutine destructor_jacobi_cpu(self)
            type(splitting_jacobi_cpu_t), intent(inout) :: self
        end subroutine

    end interface


    type, public, extends(splitting_t) :: splitting_gauss_seidel_cpu_t
        !! Splitting solver with Gauss-Seidel splitting method on CPU
    contains
        procedure :: create => create_gauss_seidel_cpu
        procedure :: apply => apply_gauss_seidel_cpu
        final     :: destructor_gauss_seidel_cpu
    end type

    interface

        module subroutine create_gauss_seidel_cpu(self, mesh)
            class(splitting_gauss_seidel_cpu_t), intent(inout) :: self
            type(mesh_cart_t), intent(in), target :: mesh
        end subroutine

        module subroutine apply_gauss_seidel_cpu(self, a, dinv, x, b)
            class(splitting_gauss_seidel_cpu_t), intent(inout) :: self
            type(csrmat_t), intent(in) :: a
            real(FP), intent(in), dimension(self%ndim) :: dinv
            real(FP), intent(inout), dimension(self%ndim) :: x
            real(FP), intent(in), dimension(self%ndim) :: b
        end subroutine

        module subroutine destructor_gauss_seidel_cpu(self)
            type(splitting_gauss_seidel_cpu_t), intent(inout) :: self
        end subroutine

    end interface


    type, public, extends(splitting_t) :: splitting_gauss_seidel_redblack_cpu_t
        !! Splitting solver with Red-Black Gauss-Seidel splitting method on CPU
        integer, private :: n_points_red
        !! Number of red points
        integer, private :: n_points_black
        !! Number of black points
        integer, private, pointer, dimension(:) :: redblack_indices
        !! Index list storing red/black/ghost points
    contains
        procedure :: create => create_gauss_seidel_redblack_cpu
        procedure :: apply => apply_gauss_seidel_redblack_cpu
        final     :: destructor_gauss_seidel_redblack_cpu
    end type

    interface

        module subroutine create_gauss_seidel_redblack_cpu(self, mesh)
            class(splitting_gauss_seidel_redblack_cpu_t), intent(inout) :: self
            type(mesh_cart_t), intent(in), target :: mesh
        end subroutine

        module subroutine apply_gauss_seidel_redblack_cpu(self, a, dinv, x, b)
            class(splitting_gauss_seidel_redblack_cpu_t), intent(inout) :: self
            type(csrmat_t), intent(in) :: a
            real(FP), intent(in), dimension(self%ndim) :: dinv
            real(FP), intent(inout), dimension(self%ndim) :: x
            real(FP), intent(in), dimension(self%ndim) :: b
        end subroutine

        module subroutine destructor_gauss_seidel_redblack_cpu(self)
            type(splitting_gauss_seidel_redblack_cpu_t), intent(inout) :: self
        end subroutine

    end interface

end module