PyTrilinos::Epetra::FECrsMatrix Class Reference

Inheritance diagram for PyTrilinos::Epetra::FECrsMatrix:

Inheritance graph
[legend]
Collaboration diagram for PyTrilinos::Epetra::FECrsMatrix:

Collaboration graph
[legend]

List of all members.

Public Member Functions

def __init__
def GlobalAssemble
def setIgnoreNonLocalEntries
def __setitem__
def __getitem__
def InsertGlobalValues
def InsertGlobalValue
def __init__
def GlobalAssemble
def setIgnoreNonLocalEntries
def __setitem__
def __getitem__
def InsertGlobalValues
def InsertGlobalValue

Public Attributes

 this

Static Public Attributes

 ROW_MAJOR = _Epetra.FECrsMatrix_ROW_MAJOR
 COLUMN_MAJOR = _Epetra.FECrsMatrix_COLUMN_MAJOR


Detailed Description

Epetra Finite-Element CrsMatrix. This class provides the ability to
input finite-element style sub-matrix data, including sub-matrices
with non-local rows (which could correspond to shared finite-element
nodes for example). This class inherits Epetra_CrsMatrix, and so all
Epetra_CrsMatrix functionality is also available.

It is intended that this class will be used as follows: Construct with
either a map or graph that describes a (non- overlapping) data
distribution.

Input data, including non-local data, using the methods
InsertGlobalValues(), SumIntoGlobalValues() and/or
ReplaceGlobalValues().

Call the method GlobalAssemble(), which gathers all non-local data
onto the owning processors as determined by the map provided at
construction. Users should note that the GlobalAssemble() method has
an optional argument which determines whether GlobalAssemble() in turn
calls FillComplete() after the data-exchange has occurred. If not
explicitly supplied, this argument defaults to true. NOTE***: When
GlobalAssemble() calls FillComplete(), it passes the arguments
'DomainMap()' and 'RangeMap()', which are the map attributes held by
the base-class CrsMatrix and its graph. If a rectangular matrix is
being assembled, the correct domain-map and range-map must be passed
to GlobalAssemble (there are two overloadings of this method) --
otherwise, it has no way of knowing what these maps should really be.

Sub-matrix data, which is assumed to be a rectangular 'table' of
coefficients accompanied by 'scatter-indices', can be provided in
three forms: Fortran-style packed 1-D array.

C-style double-pointer, or list-of-rows.

Epetra_SerialDenseMatrix object.  In all cases, a "format" parameter
specifies whether the data is laid out in row-major or column-major
order (i.e., whether coefficients for a row lie contiguously or
whether coefficients for a column lie contiguously). See the
documentation for the methods SumIntoGlobalValues() and
ReplaceGlobalValues().

Important notes: Since Epetra_FECrsMatrix inherits Epetra_CrsMatrix,
the semantics of the Insert/SumInto/Replace methods are the same as
they are on Epetra_CrsMatrix, which is:  InsertGlobalValues() inserts
values into the matrix only if the graph has not yet been finalized (
FillComplete() has not yet been called). For non-local values, the
call to InsertGlobalValues() may succeed but the GlobalAssemble()
method may then fail because the non-local data is not actually
inserted in the underlying matrix until GlobalAssemble() is called.

SumIntoGlobalValues() and ReplaceGlobalValues() only work for values
that already exist in the matrix. In other words, these methods can
not be used to put new values into the matrix.

C++ includes: Epetra_FECrsMatrix.h 

Member Function Documentation

def PyTrilinos::Epetra::FECrsMatrix::__getitem__ (   self,
  args 
)

__getitem__(self, PyObject args) -> PyObject

Reimplemented from PyTrilinos::Epetra::CrsMatrix.

def PyTrilinos::Epetra::FECrsMatrix::__getitem__ (   self,
  args 
)

__getitem__(self, PyObject args) -> PyObject

Reimplemented from PyTrilinos::Epetra::CrsMatrix.

def PyTrilinos::Epetra::FECrsMatrix::__init__ (   self,
  args 
)

__init__(self, Epetra_DataAccess CV, Map RowMap, int NumEntriesPerRow, 
    bool ignoreNonLocalEntries = False) -> FECrsMatrix
__init__(self, Epetra_DataAccess CV, Map RowMap, int NumEntriesPerRow, 
    bool ignoreNonLocalEntries = False) -> FECrsMatrix
__init__(self, Epetra_DataAccess CV, Map RowMap, Map ColMap, int NumEntriesPerRow, 
    bool ignoreNonLocalEntries = False) -> FECrsMatrix
__init__(self, Epetra_DataAccess CV, Map RowMap, Map ColMap, int NumEntriesPerRow, 
    bool ignoreNonLocalEntries = False) -> FECrsMatrix
__init__(self, Epetra_DataAccess CV, CrsGraph Graph, bool ignoreNonLocalEntries = False) -> FECrsMatrix
__init__(self, FECrsMatrix src) -> FECrsMatrix

Epetra_FECrsMatrix::Epetra_FECrsMatrix(const Epetra_FECrsMatrix &src)

Copy Constructor. 

Reimplemented from PyTrilinos::Epetra::CrsMatrix.

def PyTrilinos::Epetra::FECrsMatrix::__init__ (   self,
  args 
)

__init__(self, Epetra_DataAccess CV, Map RowMap, int NumEntriesPerRow, 
    bool ignoreNonLocalEntries = False) -> FECrsMatrix
__init__(self, Epetra_DataAccess CV, Map RowMap, int NumEntriesPerRow, 
    bool ignoreNonLocalEntries = False) -> FECrsMatrix
__init__(self, Epetra_DataAccess CV, Map RowMap, Map ColMap, int NumEntriesPerRow, 
    bool ignoreNonLocalEntries = False) -> FECrsMatrix
__init__(self, Epetra_DataAccess CV, Map RowMap, Map ColMap, int NumEntriesPerRow, 
    bool ignoreNonLocalEntries = False) -> FECrsMatrix
__init__(self, Epetra_DataAccess CV, CrsGraph Graph, bool ignoreNonLocalEntries = False) -> FECrsMatrix
__init__(self, FECrsMatrix src) -> FECrsMatrix

Epetra_FECrsMatrix::Epetra_FECrsMatrix(const Epetra_FECrsMatrix &src)

Copy Constructor. 

Reimplemented from PyTrilinos::Epetra::CrsMatrix.

def PyTrilinos::Epetra::FECrsMatrix::__setitem__ (   self,
  args 
)

__setitem__(self, PyObject args, double val)

Reimplemented from PyTrilinos::Epetra::CrsMatrix.

def PyTrilinos::Epetra::FECrsMatrix::__setitem__ (   self,
  args 
)

__setitem__(self, PyObject args, double val)

Reimplemented from PyTrilinos::Epetra::CrsMatrix.

def PyTrilinos::Epetra::FECrsMatrix::GlobalAssemble (   self,
  args 
)

GlobalAssemble(self, bool callFillComplete = True) -> int
GlobalAssemble(self, Map domain_map, Map range_map, bool callFillComplete = True) -> int

int
Epetra_FECrsMatrix::GlobalAssemble(const Epetra_Map &domain_map, const
Epetra_Map &range_map, bool callFillComplete=true)

Gather any overlapping/shared data into the non-overlapping
partitioning defined by the Map that was passed to this matrix at
construction time. Data imported from other processors is stored on
the owning processor with a "sumInto" or accumulate operation. This
is a collective method -- every processor must enter it before any
will complete it.

NOTE***: When GlobalAssemble() (the other overloading of this method)
calls FillComplete(), it passes the arguments 'DomainMap()' and
'RangeMap()', which are the map attributes already held by the base-
class CrsMatrix and its graph. If a rectangular matrix is being
assembled, the domain-map and range-map must be specified. Otherwise,
GlobalAssemble() has no way of knowing what these maps should really
be.

Parameters:
-----------

domain_map:  user-supplied domain map for this matrix

range_map:  user-supplied range map for this matrix

callFillComplete:  option argument, defaults to true. Determines
whether GlobalAssemble() internally calls the FillComplete() method on
this matrix.

error-code 0 if successful, non-zero if some error occurs 

def PyTrilinos::Epetra::FECrsMatrix::GlobalAssemble (   self,
  args 
)

GlobalAssemble(self, bool callFillComplete = True) -> int
GlobalAssemble(self, Map domain_map, Map range_map, bool callFillComplete = True) -> int

int
Epetra_FECrsMatrix::GlobalAssemble(const Epetra_Map &domain_map, const
Epetra_Map &range_map, bool callFillComplete=true)

Gather any overlapping/shared data into the non-overlapping
partitioning defined by the Map that was passed to this matrix at
construction time. Data imported from other processors is stored on
the owning processor with a "sumInto" or accumulate operation. This
is a collective method -- every processor must enter it before any
will complete it.

NOTE***: When GlobalAssemble() (the other overloading of this method)
calls FillComplete(), it passes the arguments 'DomainMap()' and
'RangeMap()', which are the map attributes already held by the base-
class CrsMatrix and its graph. If a rectangular matrix is being
assembled, the domain-map and range-map must be specified. Otherwise,
GlobalAssemble() has no way of knowing what these maps should really
be.

Parameters:
-----------

domain_map:  user-supplied domain map for this matrix

range_map:  user-supplied range map for this matrix

callFillComplete:  option argument, defaults to true. Determines
whether GlobalAssemble() internally calls the FillComplete() method on
this matrix.

error-code 0 if successful, non-zero if some error occurs 

def PyTrilinos::Epetra::FECrsMatrix::InsertGlobalValue (   self,
  args 
)

InsertGlobalValue(self, int i, int j, double val) -> int

def PyTrilinos::Epetra::FECrsMatrix::InsertGlobalValue (   self,
  args 
)

InsertGlobalValue(self, int i, int j, double val) -> int

def PyTrilinos::Epetra::FECrsMatrix::InsertGlobalValues (   self,
  args 
)

InsertGlobalValues(self, int Row, int Size, Epetra_SerialDenseVector Values, 
    Epetra_IntSerialDenseVector Entries) -> int

int
Epetra_FECrsMatrix::InsertGlobalValues(const
Epetra_IntSerialDenseVector &rows, const Epetra_IntSerialDenseVector
&cols, const Epetra_SerialDenseMatrix &values, int
format=Epetra_FECrsMatrix::COLUMN_MAJOR)

Insert a general sub-matrix into the global matrix. For square
structurally-symmetric sub-matrices, see the other overloading of this
method.

Parameters:
-----------

rows:  List of row-indices. rows.Length() must be the same as
values.M().

cols:  List of column-indices. cols.Length() must be the same as
values.N().

values:  Sub-matrix of coefficients.

format:  Optional format specifier, defaults to COLUMN_MAJOR. 

Reimplemented from PyTrilinos::Epetra::CrsMatrix.

def PyTrilinos::Epetra::FECrsMatrix::InsertGlobalValues (   self,
  args 
)

InsertGlobalValues(self, int Row, int Size, Epetra_SerialDenseVector Values, 
    Epetra_IntSerialDenseVector Entries) -> int

int
Epetra_FECrsMatrix::InsertGlobalValues(const
Epetra_IntSerialDenseVector &rows, const Epetra_IntSerialDenseVector
&cols, const Epetra_SerialDenseMatrix &values, int
format=Epetra_FECrsMatrix::COLUMN_MAJOR)

Insert a general sub-matrix into the global matrix. For square
structurally-symmetric sub-matrices, see the other overloading of this
method.

Parameters:
-----------

rows:  List of row-indices. rows.Length() must be the same as
values.M().

cols:  List of column-indices. cols.Length() must be the same as
values.N().

values:  Sub-matrix of coefficients.

format:  Optional format specifier, defaults to COLUMN_MAJOR. 

Reimplemented from PyTrilinos::Epetra::CrsMatrix.

def PyTrilinos::Epetra::FECrsMatrix::setIgnoreNonLocalEntries (   self,
  args 
)

setIgnoreNonLocalEntries(self, bool flag)

void Epetra_FECrsMatrix::setIgnoreNonLocalEntries(bool flag)

Set whether or not non-local data values should be ignored. By
default, non-local data values are NOT ignored. 

def PyTrilinos::Epetra::FECrsMatrix::setIgnoreNonLocalEntries (   self,
  args 
)

setIgnoreNonLocalEntries(self, bool flag)

void Epetra_FECrsMatrix::setIgnoreNonLocalEntries(bool flag)

Set whether or not non-local data values should be ignored. By
default, non-local data values are NOT ignored. 


The documentation for this class was generated from the following files:

Generated on Thu Dec 17 11:00:21 2009 for PyTrilinos by  doxygen 1.5.9