[Trilinos-Users] Simple way to convert Epetra CrsMatrix to Tpetra CrsMatrix
Hoemmen, Mark
mhoemme at sandia.gov
Mon Jan 20 14:03:37 MST 2014
On Jan 20, 2014, at 12:00 PM, <trilinos-users-request at software.sandia.gov>
<trilinos-users-request at software.sandia.gov> wrote:
> Date: Mon, 20 Jan 2014 14:02:08 +0000
> From: Timo Betcke <t.betcke at ucl.ac.uk>
> Subject: [Trilinos-Users] Simple way to convert Epetra CrsMatrix to
> Tpetra CrsMatrix
> To: trilinos-users at software.sandia.gov
>
> Dear Trilinos Users,
>
> I am starting to move some code from Epetra to Tpetra. To begin I'd like to
> create Tpetra Crs matrices given Epetra Crs matrices as input. Is there a
> more simple way to do this rather than iterating through the rows of the
> Epetra matrix and copy elements over to a Tpetra matrix?
Hi Timo --
You could use the setAllValues() and expertStaticFillComplete() methods of Tpetra::CrsMatrix to speed up the process. I won't say that this is "simpler" but it may perform better. Here is how it would work:
1. Create the Epetra_CrsMatrix A_e.
2. Call FillComplete on A_e.
3. Extract raw pointers to the data in A_e, and their lengths.
int* ptr;
int* ind;
double* val;
int info = A_e.ExtractCrsDataPointers (&ptr, &ind, &val);
if (info != 0) {
// … report error and exit …
}
const int numRows = A_e.Graph ().NumMyRows ();
const int nnz = A_e.Graph ().NumMyEntries ();
4. Copy data into new Teuchos::ArrayRCP arrays. (Note that ptr2 and ptr store data of different types.)
Teuchos::ArrayRCP<size_t> ptr2 (numRows+1);
Teuchos::ArrayRCP<int> ind2 (nnz);
Teuchos::ArrayRCP<double> val2 (nnz);
std::copy (ptr, ptr + numRows + 1, ptr2.begin ());
std::copy (ind, ind + nnz, ind2.begin ());
std::copy (val, val + nnz, val2.begin ());
5. Create a Tpetra::CrsMatrix A_t with the appropriate row (and column) Map(s), just as you created A_e above.
6. Give A_t the new data arrays:
A_t.setAllValues (ptr2, ind2, val2);
7. Call fillComplete on A_t. You may also use expertStaticFillComplete if you already have Import (and Export) objects constructed.
If there is sufficient interest, we may also provide a Tpetra::RowMatrix interface to an Epetra_CrsMatrix object. This would let you use Ifpack2 with Epetra objects.
mfh
More information about the Trilinos-Users
mailing list