ginkgo-project/ginkgo

How to convert TrilinosWrappers::SparseMatrix and dealii::LinearAlgebra::distributed::Vector for Ginkgo

acse-bn20 opened this issue · 3 comments

Hi,

I am trying to set up the Ginkgo cuda GMRES solver with these inputs. I have read the tutorial and

using mtx = gko::matrix::Csr<double>;
const auto exec = gko::CudaExecutor::create();
auto matrix = gko::share(mtx::create(exec, gko::dim<2>(size)));

Then fill the matrix by looping through and then lend the ownership. I was wondering whether it is possible to do this process with pointers to copy the data TrilinosWrappers::SparseMatrix (LHS) and dealii::LinearAlgebra::distributed::Vector (RHS) for a suitable format for Ginkgo? I am especially curious how to deal with the sparse matrix conversion for ginkgo.

If I understand you correctly, you want to do a single copy instead of looping through the matrix to fill the Ginkgo types. This is very much possible, but it depends on how you can get access to the underlying pointers to the trillinos or dealii data. Assuming you can get the pointers here is what you can do:

// the type of the pointer can be changed, depending on what trillinos stores
int* trillinos_row_ptrs = ...;  // the pointer to the trillinos row pointers/offsets
int* trillinos_col_idxs = ...;  // the pointer to the trillinos column indices
double* trillinos_values = ...;  // the pointer to the trillinos CSR values

// creating Ginkgo arrays and copy assign views of the trillinos data
gko::array<int> row_ptrs{exec};
gko::array<int> col_idxs{exec};
gko::array<double> values{exec};

// this assumes that the trillinos data is stored on the CPU
auto exec_host = exec->get_master();
row_ptrs = gko::make_array_view(exec_host, num_rows + 1, trillinos_row_ptrs);
col_idxs = gko::make_array_view(exec_host, num_nnz, trillinos_col_idxs);
values = gko::make_array_view(exec_host, num_nnz, trillinos_values);

// now create the Ginkgo Csr from the Ginkgo arrays
auto mat = gko::matrix::Csr<double, int>::create(exec, gko::dim<2>{num_rows, num_cols}, std::move(row_ptrs), 
                                                 std::move(col_idxs), std::move(values));

For vectors, the same approach applies.

One complication would be if trillinos or dealii only allow immutable access to their data, i.e. they return const double* or similar. In this case, you would need to do:

const double* values = ...;

gko::array<double> values{exec, gko::make_const_array_view(exec_host, num_nnz, values).copy_to_array()};

One thing to note is that the dealii::LinearAlgebra::distributed::Vector, and maybe also the TrillinosWrappers::SparseMatrix has ghost values, while the closet Ginkgo types do not. So you need to take care how to handle that, until we also provide a Ginkgo types with ghost values.

I see. Thank you very much for the information and the sample.

If there are no follow-up questions, I will close this issue. Please feel free to reopen this if there are still some questions left.