splitting_gauss_seidel_cpu_s.f90 Source File


This file depends on

sourcefile~~splitting_gauss_seidel_cpu_s.f90~~EfferentGraph sourcefile~splitting_gauss_seidel_cpu_s.f90 splitting_gauss_seidel_cpu_s.f90 sourcefile~splitting_m.f90 splitting_m.f90 sourcefile~splitting_gauss_seidel_cpu_s.f90->sourcefile~splitting_m.f90 sourcefile~comm_handling_m.f90 comm_handling_m.f90 sourcefile~splitting_m.f90->sourcefile~comm_handling_m.f90 sourcefile~csrmat_m.f90 csrmat_m.f90 sourcefile~splitting_m.f90->sourcefile~csrmat_m.f90 sourcefile~mesh_cart_m.f90 mesh_cart_m.f90 sourcefile~splitting_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~precision_m.f90 precision_m.f90 sourcefile~splitting_m.f90->sourcefile~precision_m.f90 sourcefile~csrmat_m.f90->sourcefile~precision_m.f90 sourcefile~error_handling_m.f90 error_handling_m.f90 sourcefile~csrmat_m.f90->sourcefile~error_handling_m.f90 sourcefile~list_operations_m.f90 list_operations_m.f90 sourcefile~csrmat_m.f90->sourcefile~list_operations_m.f90 sourcefile~screen_io_m.f90 screen_io_m.f90 sourcefile~csrmat_m.f90->sourcefile~screen_io_m.f90 sourcefile~status_codes_m.f90 status_codes_m.f90 sourcefile~csrmat_m.f90->sourcefile~status_codes_m.f90 sourcefile~mesh_cart_m.f90->sourcefile~comm_handling_m.f90 sourcefile~mesh_cart_m.f90->sourcefile~precision_m.f90 sourcefile~descriptors_m.f90 descriptors_m.f90 sourcefile~mesh_cart_m.f90->sourcefile~descriptors_m.f90 sourcefile~equilibrium_m.f90 equilibrium_m.f90 sourcefile~mesh_cart_m.f90->sourcefile~equilibrium_m.f90 sourcefile~mesh_cart_m.f90->sourcefile~error_handling_m.f90 sourcefile~slab_equilibrium.f90 slab_equilibrium.f90 sourcefile~mesh_cart_m.f90->sourcefile~slab_equilibrium.f90 sourcefile~mesh_cart_m.f90->sourcefile~status_codes_m.f90 sourcefile~descriptors_m.f90->sourcefile~screen_io_m.f90 sourcefile~equilibrium_m.f90->sourcefile~precision_m.f90 sourcefile~error_handling_m.f90->sourcefile~comm_handling_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~list_operations_m.f90->sourcefile~precision_m.f90 sourcefile~list_operations_m.f90->sourcefile~screen_io_m.f90 sourcefile~screen_io_m.f90->sourcefile~precision_m.f90 sourcefile~slab_equilibrium.f90->sourcefile~precision_m.f90 sourcefile~slab_equilibrium.f90->sourcefile~descriptors_m.f90 sourcefile~slab_equilibrium.f90->sourcefile~equilibrium_m.f90 sourcefile~slab_equilibrium.f90->sourcefile~screen_io_m.f90 sourcefile~params_equi_slab_m.f90 params_equi_slab_m.f90 sourcefile~slab_equilibrium.f90->sourcefile~params_equi_slab_m.f90 sourcefile~params_equi_slab_m.f90->sourcefile~precision_m.f90 sourcefile~params_equi_slab_m.f90->sourcefile~error_handling_m.f90 sourcefile~params_equi_slab_m.f90->sourcefile~screen_io_m.f90 sourcefile~params_equi_slab_m.f90->sourcefile~status_codes_m.f90

Source Code

submodule(splitting_m) splitting_gauss_seidel_cpu_s
    !! Splitting solver with Gauss-Seidel splitting method on CPU
    implicit none

contains

    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

        ! Store dimension
        self%ndim = mesh%get_n_points()

    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

        integer :: m, n
        real(FP) :: sg

        !$OMP MASTER
        do m = 1, a%ndim
            sg = 0.0_FP
            do n = a%i(m), a%i(m+1)-1
                if (a%j(n) /= m) then
                    sg = sg + a%val(n)*x(a%j(n))
                endif
            enddo
            x(m) = (b(m) - sg) * dinv(m)
        enddo
        !$OMP END MASTER
        !$OMP BARRIER

    end subroutine

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

    end subroutine

end submodule