[Trilinos-Users] Rectangular Epetra_CrsMatrix sparse matrices

Baker, Christopher G. bakercg at ornl.gov
Wed Jun 15 09:58:50 MDT 2011


Hello Alberto,

The shape of a sparse matrix, of course, isn't necessarily defined by the
entries. For Epetra, the shape of the matrix is defined by the range and
domain maps, I.e., the maps for the linear subspaces between which the
matrix is an operator.

What this means is that defining a rectangular Epetra_CrsMatrix consists
of specifying distinct domain and range maps the matrix. These are passed
to FillComplete(). In the case of your example, you can pass the row map
used to construct the matrix as the range map of the matrix, but you
should pass a 3-dimensional map as the domain map. Documentation for this
version of FillComplete() is here:
http://trilinos.sandia.gov/packages/docs/r10.6/packages/epetra/doc/html/cla
ssEpetra__CrsMatrix.html#aac3c13275e8dbe80826e820abbb8e49c


An example of this usage is available here:
Trilinos/packages/epetra/test/CrsRectMatrix/cxx_main.cpp

Chris

On 6/15/11 11:37 AM, "Alberto F. Martín-Huertas" <amartin at cimne.upc.edu>
wrote:

>
> Dear all,
>
> the Epetra_CrsMatrix class reference explicitly mentions
> that Epetra_CrsMatrix matrices can be rectangular. I have thoroghly
> read Epetra_CrsMatrix class reference I did not manage to see
> any further reference on how to create rectangular Epetra_CrsMatrix.
> I tried with  piece of code below, which tries to create a 2x3
> Epetra_CrsMatrix without success (in particular
> it throws an uncatched exception on matrix->FillComplete()).
>
> I would greatly appreciate any hint on how to deal with
> rectangular.
>
> Thanks.
>
> Best regards,
>  Alberto.
>
>
>
> RCP<Epetra_CrsMatrix> build2x3( double a00, double a01, double a02,
>                                 double a10, double a11, double a12,
> Epetra_Comm & Comm)
> {
>    RCP<Epetra_Map> map = rcp(new Epetra_Map(2,0,Comm));
>    RCP<Epetra_CrsMatrix> matrix = rcp(new
> Epetra_CrsMatrix(Copy,*map,3));
>
>    int MyPID = Comm.MyPID();
>    int indicies[3] = { 0, 1, 2 };
>    double values[3];
>
>    int num_my_elements,   i ;
>    int * my_global_elements ;
>
>    num_my_elements    = map->NumMyElements();
>    my_global_elements = map->MyGlobalElements();
>
>    for ( i = 0 ; i < num_my_elements ; i++)
>    {
>         switch ( my_global_elements[i] )
>         {
>           case 0:
>               values[0] = a00;
>               values[1] = a01;
>               values[2] = a02;
>           break;
>           case 1:
>               values[0] = a10;
>               values[1] = a11;
>               values[2] = a12;
>           break;
>         }
>         matrix->InsertGlobalValues(my_global_elements[i], 3,
> values,indicies);
>    }
>
>    // pack it up
>    matrix->FillComplete();
>
>    return matrix;
> }
>
>
>
>
>
>-- 
> Alberto F. Martín-Huertas
> Centre Internacional de Mètodes Numèrics a l'Enginyeria (CIMNE)
> Parc Mediterrani de la Tecnologia, UPC
> Esteve Terrades 5, Building C3, Office 207,
> 08860 Castelldefels (Barcelona, Spain)
> e-mail: amartin at cimne . upc . edu
>
> ________________
> IMPORTANT NOTICE
> All personal data contained on this mail will be processed
> confidentially and registered in a file property of CIMNE in
> order to manage corporate communications. You may exercise the rights
> of access, rectification, erasure and object by
> letter sent to Ed. C1 Campus Norte UPC. Gran Capitán s/n Barcelona.
>
>
>_______________________________________________
>Trilinos-Users mailing list
>Trilinos-Users at software . sandia . gov
>hxxp://software.sandia.gov/mailman/listinfo/trilinos-users





More information about the Trilinos-Users mailing list