Skip to content

SPARSE::spiluk

Vinh Dang edited this page Feb 19, 2020 · 7 revisions

KokkosSparse::spiluk()

Header File: KokkosSparse_spiluk.hpp

Usage: KokkosSparse::spiluk_symbolic(handle, k, A_rowmap, A_entries, L_rowmap, L_entries, U_rowmap, U_entries);

Usage: KokkosSparse::spiluk_numeric(handle, k, A_rowmap, A_entries, A_values, L_rowmap, L_entries, L_values, U_rowmap, U_entries, U_values);

Perform incomplete LU factorization ILU(k) on matrices stored in compressed row sparse ("Crs") format.

Interface

  template <typename KernelHandle,
            typename ARowMapType,
            typename AEntriesType,
            typename LRowMapType,
            typename LEntriesType,
            typename URowMapType,
            typename UEntriesType>
  void spiluk_symbolic(KernelHandle *handle,
                       typename KernelHandle::const_nnz_lno_t fill_lev,
                       ARowMapType&  A_rowmap,
                       AEntriesType& A_entries,
                       LRowMapType&  L_rowmap,
                       LEntriesType& L_entries,
                       URowMapType&  U_rowmap,
                       UEntriesType& U_entries);

  template <typename KernelHandle,
            typename ARowMapType,
            typename AEntriesType,
            typename AValuesType,
            typename LRowMapType,
            typename LEntriesType,
            typename LValuesType,
            typename URowMapType,
            typename UEntriesType,
            typename UValuesType>
  void spiluk_numeric(KernelHandle *handle,
                      typename KernelHandle::const_nnz_lno_t fill_lev,
                      ARowMapType&  A_rowmap,
                      AEntriesType& A_entries,
                      AValuesType&  A_values,
                      LRowMapType&  L_rowmap,
                      LEntriesType& L_entries,
                      LValuesType&  L_values,
                      URowMapType&  U_rowmap,
                      UEntriesType& U_entries,
                      UValuesType&  U_values);

Parameters:

spiluk_symbolic
  • KernelHandle
  • InputLocalOrdinalType: Integral-type fill level
  • ARowMapType: a rank-1 Kokkos::View of input matrix row map
  • AEntriesType: a rank-1 Kokkos::View of input matrix entries
  • LRowMapType: a rank-1 Kokkos::View of L matrix row map
  • LEntriesType: a rank-1 Kokkos::View of L matrix entries
  • URowMapType: a rank-1 Kokkos::View of U matrix row map
  • UEntriesType: a rank-1 Kokkos::View of U matrix entries.
spiluk_numeric
  • KernelHandle
  • InputLocalOrdinalType: Integral-type fill level
  • ARowMapType: a rank-1 Kokkos::View of input matrix row map
  • AEntriesType: a rank-1 Kokkos::View of input matrix entries
  • AValuesType: a rank-1 Kokkos::View of input matrix values
  • LRowMapType: a rank-1 Kokkos::View of L matrix row map
  • LEntriesType: a rank-1 Kokkos::View of L matrix entries
  • LValuesType: a rank-1 Kokkos::View of L matrix values
  • URowMapType: a rank-1 Kokkos::View of U matrix row map
  • UEntriesType: a rank-1 Kokkos::View of U matrix entries
  • UValuesType: a rank-1 Kokkos::View of U matrix values.

Requirements:

  • Creation of a KernelHandle
  • All Views must have rank 1
  • ARowMapType::non_const_value_type, LRowMapType::non_const_value_type, and URowMapType::non_const_value_type must match KernelHandle::size_type
  • AEntriesType::non_const_value_type, LEntriesType::non_const_value_type, and UEntriesType::non_const_value_type must match KernelHandle::nnz_lno_t
  • AValuesType::value_type, LValuesType::value_type, and UValuesType::value_type must match KernelHandle::nnz_scalar_t
  • LRowMapType::value_type == LRowMapType::non_const_value_type
  • LEntriesType::value_type == LEntriesType::non_const_value_type
  • LValuesType::value_type == LValuesType::non_const_value_type
  • URowMapType::value_type == URowMapType::non_const_value_type
  • UEntriesType::value_type == UEntriesType::non_const_value_type
  • UValuesType::value_type == UValuesType::non_const_value_type
  • ARowMapType, LRowMapType, and URowMapType, AEntriesType, LEntriesType, and UEntriesType, AValuesType, LValuesType, and UValuesType have same device_types
  • KernelHandle and Views have same execution spaces
  • fill_lev >= 0

Example

#include <Kokkos_Core.hpp>
#include <KokkosSparse_CrsMatrix.hpp>
#include <KokkosSparse_spiluk.hpp>
#include <KokkosKernels_IOUtils.hpp>

int main(int argc, char* argv[]) {
   Kokkos::initialize();

   using scalar_t  = double;
   using lno_t     = int;
   using size_type = int;
   using crsMat_t  = typename KokkosSparse::CrsMatrix<scalar_t, lno_t, Kokkos::DefaultExecutionSpace, void, size_type>;

   using graph_t         = typename crsmat_t::StaticCrsGraphType;
   using lno_view_t      = typename graph_t::row_map_type::non_const_type;
   using lno_nnz_view_t  = typename graph_t::entries_type::non_const_type;
   using scalar_view_t   = typename crsmat_t::values_type::non_const_type;

   using ViewVectorType  = Kokkos::View<scalar_t*>;
   using execution_space = typename ViewVectorType::device_type::execution_space;
   using memory_space    = typename ViewVectorType::device_type::memory_space;

   using KernelHandle    = KokkosKernels::Experimental::KokkosKernelsHandle <size_type, lno_t, scalar_t, execution_space, memory_space, memory_space> ;

   // Read and fill matrix
   crsmat_t A        = KokkosKernels::Impl::read_kokkos_crst_matrix<crsmat_t>("mtx filename");
   graph_t  graph    = A.graph;
   const size_type N = graph.numRows();
   typename KernelHandle::const_nnz_lno_t fill_lev = lno_t(2) ;
   const size_type nnzA = A.graph.entries.extent(0);

   // Create KokkosKernelHandle with an spiluk algorithm, limited by configuration at compile-time and set via the handle
   // Some options: {SEQLVLSCHD_RP, SEQLVLSCHD_TP1}
   KernelHandle kh;

   //kh.create_spiluk_handle(KokkosSparse::Experimental::SPILUKAlgorithm::SEQLVLSCHD_RP, N, EXPAND_FACT*nnzA*(fill_lev+1), EXPAND_FACT*nnzA*(fill_lev+1));
   kh.create_spiluk_handle(KokkosSparse::Experimental::SPILUKAlgorithm::SEQLVLSCHD_TP1, N, EXPAND_FACT*nnzA*(fill_lev+1), EXPAND_FACT*nnzA*(fill_lev+1));

   auto spiluk_handle = kh.get_spiluk_handle();
  
   lno_view_t     L_row_map("L_row_map", N + 1);
   lno_nnz_view_t L_entries("L_entries", spiluk_handle->get_nnzL());
   scalar_view_t  L_values ("L_values",  spiluk_handle->get_nnzL());
   lno_view_t     U_row_map("U_row_map", N + 1);
   lno_nnz_view_t U_entries("U_entries", spiluk_handle->get_nnzU());
   scalar_view_t  U_values ("U_values",  spiluk_handle->get_nnzU());
  
   KokkosSparse::Experimental::spiluk_symbolic( &kh, fill_lev, A.graph.row_map, A.graph.entries, 
                                                               L_row_map, L_entries, U_row_map, U_entries );

   Kokkos::resize(L_entries, spiluk_handle->get_nnzL());
   Kokkos::resize(L_values,  spiluk_handle->get_nnzL());
   Kokkos::resize(U_entries, spiluk_handle->get_nnzU());
   Kokkos::resize(U_values,  spiluk_handle->get_nnzU());

   spiluk_handle->set_team_size(16);
	  
   KokkosSparse::Experimental::spiluk_numeric( &kh, fill_lev, 
                                               A.graph.row_map, A.graph.entries, A.values, 
                                               L_row_map, L_entries, L_values, U_row_map, U_entries, U_values );

   kh.destroy_spiluk_handle();

   Kokkos::finalize();
}
Clone this wiki locally