[Trilinos-Users] Epetra_CrsMatrix Basics

Baker, Christopher G. bakercg at ornl.gov
Mon Nov 16 09:31:38 MST 2009


Tim,

In general, you do not want to specify a column map for the sparse matrix; this will be constructed for you based on the non-zeros entries of the matrix, and is not necessary or desired in the local, non-distributed case.

If you remove the declaration of column_map and call the Epetra_CrsMatrix constructor that accepts row map and NNZ information, that should address the problem.

Note, if your matrix is rectangular,  you will need to pass maps for the domain and range space into FillComplete(); otherwise, don't worry about this.

Lastly, if your VTK matrix can easily give you information regarding the number of non-zeros per row, it would good to pass that info to the CrsMatrix constructor so that it can more efficiently allocate the data. The first step would be to generate a guess at the number of non-zeros per row and use that, along with a dynamic allocation (StaticProfile==false). The second would be to specify the number of non-zeros the rows and use a static allocation (StaticProfile==0).

Chris


On 11/16/09 11:16 AM, "Timothy M. Shead" <tshead at sandia.gov> wrote:

Folks:

I'm working on the integration between Trilinos and VTK, and running
into an area where I'm not sure what I'm doing wrong.  I have a toy
example where I'm converting a VTK sparse matrix to Epetra_CrsMatrix,
then multiplying it with an Epetra_MultiVector.  At this point, all I'm
trying to do is serial, local computation:

   Epetra_LocalMap row_map(vtk_matrix->RowCount, 0, communicator);
   Epetra_LocalMap column_map(vtk_matrix->ColumnCount, 0, communicator);

   Epetra_CrsMatrix* const epetra_matrix = new Epetra_CrsMatrix(Copy,
row_map, column_map, 0, false);

   for(vtkIdType n = 0; n != vtk_matrix->NonZeroCount; ++n)
     {
     int row = vtk_matrix->GetRowN(n);
     int column = vtk_matrix->GetColumnN(n);
     double value = vtk_matrix->GetValueN(n);
     epetra_matrix->InsertGlobalValues(row, 1, &value, &column);
     }

   epetra_matrix->FillComplete();

However, when I try to perform the multiplication, I'm getting some
errors that I don't understand:

   epetra_matrix->Multiply(false, epetra_input_multivec,
epetra_output_multivec);

   Epetra ERROR -3,
/Users/tshead/src/trilinos/packages/epetra/src/Epetra_DistObject.cpp,
line 110
   Epetra ERROR -3,
/Users/tshead/src/trilinos/packages/epetra/src/Epetra_CrsMatrix.cpp,
line 2461

These errors have something to do with an "importer", which leads me to
believe that I'm not setting-up my matrix correctly for the serial,
undistributed use-case ... any thoughts?

Thanks in advance ...
Tim
_______________________________________________
Trilinos-Users mailing list
Trilinos-Users at software.sandia.gov
http://software.sandia.gov/mailman/listinfo/trilinos-users






More information about the Trilinos-Users mailing list