FEI  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Classes | Public Member Functions | List of all members
fei::MatrixGraph Class Referenceabstract

#include <fei_MatrixGraph.hpp>

Inheritance diagram for fei::MatrixGraph:
Inheritance graph
[legend]

Classes

class  Factory
 

Public Member Functions

virtual ~MatrixGraph ()
 
virtual void setParameters (const fei::ParameterSet &params)=0
 
virtual void setRowSpace (fei::SharedPtr< fei::VectorSpace > rowSpace)=0
 
virtual fei::SharedPtr
< fei::VectorSpace
getRowSpace ()=0
 
virtual void setColumnSpace (fei::SharedPtr< fei::VectorSpace > columnSpace)=0
 
virtual fei::SharedPtr
< fei::VectorSpace
getColSpace ()=0
 
virtual int definePattern (int numIDs, int idType)=0
 
virtual int definePattern (int numIDs, int idType, int fieldID)=0
 
virtual int definePattern (int numIDs, int idType, const int *numFieldsPerID, const int *fieldIDs)=0
 
virtual int definePattern (int numIDs, const int *idTypes, const int *numFieldsPerID, const int *fieldIDs)=0
 
virtual int initConnectivityBlock (int blockID, int numConnectivityLists, int patternID, bool diagonal=false)=0
 
virtual int initConnectivityBlock (int numConnectivityLists, int patternID, bool diagonal=false)=0
 
virtual int initConnectivityBlock (int blockID, int numConnectivityLists, int rowPatternID, int colPatternID)=0
 
virtual int initConnectivity (int blockID, int connectivityID, const int *connectedIdentifiers)=0
 
virtual int initConnectivity (int blockID, int connectivityID, const int *rowConnectedIdentifiers, const int *colConnectedIdentifiers)=0
 
virtual int initConnectivity (int patternID, const int *connectedIdentifiers)=0
 
virtual int initConnectivity (int rowPatternID, const int *rowConnectedIdentifiers, int colPatternID, const int *colConnectedIdentifiers)=0
 
virtual int initConnectivity (int idType, int numRows, const int *rowIDs, const int *rowOffsets, const int *packedColumnIDs)=0
 
virtual int initConnectivity (int idType, int fieldID, int numRows, const int *rowIDs, const int *rowOffsets, const int *packedColumnIDs)=0
 
virtual int initConnectivity (int idType, int numRows, const int *rowIDs, const int *rowLengths, const int *const *columnIDs)=0
 
virtual int initLagrangeConstraint (int constraintID, int constraintIDType, int numIDs, const int *idTypes, const int *IDs, const int *fieldIDs)=0
 
virtual int initPenaltyConstraint (int constraintID, int constraintIDType, int numIDs, const int *idTypes, const int *IDs, const int *fieldIDs)=0
 
virtual int initSlaveConstraint (int numIDs, const int *idTypes, const int *IDs, const int *fieldIDs, int offsetOfSlave, int offsetIntoSlaveField, const double *weights, double rhsValue)=0
 
virtual bool hasSlaveDof (int ID, int idType)=0
 
virtual int initComplete ()=0
 
virtual fei::SharedPtr
< fei::SparseRowGraph
createGraph (bool blockEntryGraph, bool localRowGraph_includeSharedRows=false)=0
 
virtual int compareStructure (const fei::MatrixGraph &matrixGraph, bool &equivalent) const =0
 
virtual int getNumConnectivityBlocks () const =0
 
virtual std::map< int,
fei::ConnectivityBlock * > & 
getConnectivityBlocks ()=0
 
virtual int getConnectivityBlockIDs (std::vector< int > &blockIDs) const =0
 
virtual int getNumIDsPerConnectivityList (int blockID) const =0
 
virtual int getConnectivityNumIndices (int blockID) const =0
 
virtual int getConnectivityNumIndices (int blockID, int &numRowIndices, int &numColIndices)=0
 
virtual int getConnectivityIndices (int blockID, int connectivityID, int indicesAllocLen, int *indices, int &numIndices)=0
 
virtual int getConnectivityIndices (int blockID, int connectivityID, int rowIndicesAllocLen, int *rowIndices, int &numRowIndices, int colIndicesAllocLen, int *colIndices, int &numColIndices)=0
 
virtual int getPatternNumIndices (int patternID, int &numIndices)=0
 
virtual int getPatternIndices (int patternID, const int *IDs, std::vector< int > &indices)=0
 
virtual int getLocalNumLagrangeConstraints () const =0
 
virtual int getGlobalNumSlaveConstraints () const =0
 
virtual ConstraintTypegetLagrangeConstraint (int constraintID)=0
 
virtual std::map< int,
ConstraintType * > & 
getLagrangeConstraints ()=0
 
virtual ConstraintTypegetPenaltyConstraint (int constraintID)=0
 
virtual ConstraintTypegetSlaveConstraint (int constraintID)=0
 
virtual int getConstraintConnectivityIndices (ConstraintType *cr, std::vector< int > &globalIndices)=0
 
virtual const
fei::ConnectivityBlock
getConnectivityBlock (int blockID) const =0
 
virtual fei::ConnectivityBlockgetConnectivityBlock (int blockID)=0
 
virtual void setIndicesMode (int mode)=0
 
virtual fei::SharedPtr
< fei::FillableMat > 
getSlaveDependencyMatrix ()=0
 
virtual fei::PatterngetPattern (int patternID)=0
 
virtual int createSlaveMatrices ()=0
 
virtual fei::SharedPtr
< fei::Reducer > 
getReducer ()=0
 
virtual fei::SharedPtr
< fei::SparseRowGraph
getRemotelyOwnedGraphRows ()=0
 
virtual void getConstrainedIndices (std::vector< int > &crindices) const =0
 

Detailed Description

A container for the data that defines connectivity, and which will ultimately be used to generate a matrix graph.

Definition at line 33 of file fei_MatrixGraph.hpp.

Constructor & Destructor Documentation

virtual fei::MatrixGraph::~MatrixGraph ( )
inlinevirtual

Destructor.

Definition at line 57 of file fei_MatrixGraph.hpp.

Member Function Documentation

virtual void fei::MatrixGraph::setParameters ( const fei::ParameterSet params)
pure virtual

Set parameters from a ParameterSet object. Currently two parameters are recognized: "debugOutput 'path'" where 'path' is the path to the location where debug-log files will be produced.
"name 'string'" where 'string' is an identifier that will be used in debug-log file-names.

Implemented in fei::MatrixGraph_Impl2.

virtual void fei::MatrixGraph::setRowSpace ( fei::SharedPtr< fei::VectorSpace rowSpace)
pure virtual

Provide a VectorSpace to be used for looking up indices, field-masks, etc., for the row-space. If no column-VectorSpace is provided, it will be assumed that the column-space equals the row-space.

Returns
error-code 0 if successful

Implemented in fei::MatrixGraph_Impl2.

virtual fei::SharedPtr<fei::VectorSpace> fei::MatrixGraph::getRowSpace ( )
pure virtual

Obtain the VectorSpace that corresponds to the row-space for this MatrixGraph object.

Implemented in fei::MatrixGraph_Impl2.

virtual void fei::MatrixGraph::setColumnSpace ( fei::SharedPtr< fei::VectorSpace columnSpace)
pure virtual

Provide a VectorSpace to be used for looking up indices, field-masks, etc., for the column-space. If no column-VectorSpace is provided, it will be assumed that the column-space equals the row-space.

Returns
error-code 0 if successful

Implemented in fei::MatrixGraph_Impl2.

virtual fei::SharedPtr<fei::VectorSpace> fei::MatrixGraph::getColSpace ( )
pure virtual

Obtain the VectorSpace that corresponds to the column-space for this MatrixGraph object.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::definePattern ( int  numIDs,
int  idType 
)
pure virtual

Define a pattern to use for subsequent blocked-contributions. Examples include element-contributions. Return an int patternID that can be used to reference this pattern in calls to initConnectivityBlock, etc

This is the simplest of the pattern-definition methods. IMPORTANT NOTE: this method does not associate a field with the identifiers. Only use this method for problems where you explicitly don't want or need to associate fields with identifiers. Examples would include problems where only a single scalar field exists across the entire mesh and thus doesn't need to be explicitly referenced. Other cases where this might be used is for non finite-element problems that don't have identifier/field pairs.

Parameters
numIDsInput. number of identifiers per pattern 'instance'.
idTypeInput. Specifies which type of identifiers are associated with instances of this pattern. Must be one of the idTypes defined for a VectorSpace that is associated with this MatrixGraph. idTypes are defined via the method VectorSpace::defineIDTypes().
Returns
patternID

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::definePattern ( int  numIDs,
int  idType,
int  fieldID 
)
pure virtual

Define a pattern to use for subsequent blocked-contributions. Examples include element-contributions. Return an int patternID that can be used to reference this pattern in calls to initConnectivityBlock, etc

This is the simplest of the 3 pattern-definition methods that associate fields with identifiers (there is one pattern-definition method above that allows for specifying a pattern of identifiers that don't have associated fields). This method defines patterns for contributions where a single field is associated with each identifier in a list of identifiers, and all the identifiers in the list are of the same type.

Parameters
numIDsInput. number of identifiers per pattern 'instance'.
idTypeInput. Specifies which type of identifiers are associated with instances of this pattern. Must be one of the idTypes defined for a VectorSpace that is associated with this MatrixGraph. idTypes are defined via the method VectorSpace::defineIDTypes().
fieldIDInput. field-identifier for the single field that is to reside at each identifier.
Returns
patternID

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::definePattern ( int  numIDs,
int  idType,
const int *  numFieldsPerID,
const int *  fieldIDs 
)
pure virtual

Define a pattern to use for subsequent blocked-contributions. Examples include element-contributions. Return an int patternID that can be used to reference this pattern in calls to initConnectivityBlock, etc.

This is the 'middle' of the pattern-definition methods, in terms of the complexity of pattern that can be defined. This method defines patterns for contributions where the identifiers are all of the same type, but an arbitrary list of fields can be associated with each identifier.

Parameters
numIDsInput. number of identifiers per pattern 'instance'.
idTypeInput. Specifies which type of identifiers are associated with instances of this pattern. Must be one of the idTypes defined for a VectorSpace that is associated with this MatrixGraph. idTypes are defined via the method VectorSpace::defineIDTypes().
numFieldsPerIDInput. List of length numIDs. i-th entry ives the number of fields to be associated with the i-th identifier in a contribution.
fieldIDsInput. Packed list of length sum(numFieldsPerID[i]). Contains the fieldIDs to be associated with the identifiers for a contribution.
Returns
patternID Input. Identifier to be used later when referring to this pattern.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::definePattern ( int  numIDs,
const int *  idTypes,
const int *  numFieldsPerID,
const int *  fieldIDs 
)
pure virtual

Define a pattern to use for subsequent blocked-contributions. Examples include element-contributions. Return an int patternID that can be used to reference this pattern in calls to initConnectivityBlock, etc.

This is the most general of the pattern-definition methods. This method defines a pattern consisting of a mixture of identifier-types, with each identifier having an arbitrary list of associated fields.

Parameters
numIDsInput. number of identifiers per pattern 'instance'.
idTypesInput. List of length numIDs. Specifies the type of each identifier to be contributed for instances of this pattern. Each of the idTypes must be one of the idTypes defined for a VectorSpace that is associated with this MatrixGraph. idTypes are defined via the method VectorSpace::defineIDTypes().
numFieldsPerIDInput. List of length numIDs. i-th entry gives the number of fields to be associated with the i-th identifier in a contribution.
fieldIDsInput. Packed list of length sum(numFieldsPerID[i]). Contains the fieldIDs to be associated with the identifiers for a contribution.
Returns
patternID Input. Identifier to be used later when referring to this pattern.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::initConnectivityBlock ( int  blockID,
int  numConnectivityLists,
int  patternID,
bool  diagonal = false 
)
pure virtual

Initialize a block of connectivity contributions. An example is a block of elements which share a common layout of nodes/fields per element.
This method accepts only one pattern-id, implying that connectivities in this block describe a symmetric structure. See the other overloading of this method for the non-symmetric case.

Parameters
blockIDInput. User-specified identifier for this block. Will generally be required to be non-negative.
numConnectivityListsInput. Number of connectivity-lists that will be supplied for this block.
patternIDInput. Descriptor for the connectivities to be provided. Must be a pattern that was previously defined via definePattern().
diagonalOptional argument, defaults to false. If specified as true, each connectivity list will only contribute diagonal entries to the graph. This is used if the connectivity-block represents a collection of lumped- mass submatrix contributions, or something similar.
Returns
error-code 0 if successful

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::initConnectivityBlock ( int  numConnectivityLists,
int  patternID,
bool  diagonal = false 
)
pure virtual

Initialize a block of connectivity contributions. An example is a block of elements which share a common layout of nodes/fields per element.
This method accepts only one pattern-id, implying that connectivities in this block describe a symmetric structure. See the other overloading of this method for the non-symmetric case.

Parameters
blockIDInput. User-specified identifier for this block. Will generally be required to be non-negative.
numConnectivityListsInput. Number of connectivity-lists that will be supplied for this block.
patternIDInput. Descriptor for the connectivities to be provided. Must be a pattern that was previously defined via definePattern().
diagonalOptional argument, defaults to false. If specified as true, each connectivity list will only contribute diagonal entries to the graph. This is used if the connectivity-block represents a collection of lumped- mass submatrix contributions, or something similar.
Returns
identifier for the new connectivity-block.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::initConnectivityBlock ( int  blockID,
int  numConnectivityLists,
int  rowPatternID,
int  colPatternID 
)
pure virtual

Initialize a block of connectivity contributions. An example is a block of elements which share a common layout of nodes/fields per element.
This method accepts two pattern-ids, implying that connectivities in this block describe a non-symmetric structure. See the other overloading of this method for the symmetric case.

Parameters
blockIDInput. User-specified identifier for this block. Will generally be required to be non-negative.
numConnectivityListsInput. Number of connectivity-lists that will be supplied for this block.
rowPatternIDInput. Descriptor for the row-connectivities to be provided. Must be a pattern that was previously defined via definePattern().
colPatternIDInput. Descriptor for the column-connectivities to be provided. Must be a pattern that was previously defined via definePattern().
Returns
error-code 0 if successful

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::initConnectivity ( int  blockID,
int  connectivityID,
const int *  connectedIdentifiers 
)
pure virtual

Make a contribution to the MatrixGraph's connectivity. Examples would include element-node connectivity lists, etc.

Parameters
blockIDInput. Must correspond to a blockID that was previously used in a call to initConnectivityBlock().
connectivityIDInput. Identifier for this connectivity list. May be an element-identifier, etc.
connectedIdentifiersInput. List of the identifiers that form this connectivity list.
Returns
error-code 0 if successful

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::initConnectivity ( int  blockID,
int  connectivityID,
const int *  rowConnectedIdentifiers,
const int *  colConnectedIdentifiers 
)
pure virtual

Make a contribution to the MatrixGraph's connectivity. This overloading of initConnectivity() provides for structurally non-symmetric entries.

Parameters
blockIDInput. Must correspond to a blockID that was previously used in a call to initConnectivityBlock().
connectivityIDInput. Identifier for this connectivity list. May be an element-identifier, etc.
rowConnectedIdentifiersInput. List of the identifiers that form the connectivity list for the row-space.
colConnectedIdentifiersInput. List of the identifiers that form the connectivity list for the column-space.
Returns
error-code 0 if successful

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::initConnectivity ( int  patternID,
const int *  connectedIdentifiers 
)
pure virtual

Make a contribution to the MatrixGraph's connectivity. This overloading of initConnectivity() assumes structurally symmetric entries.

Parameters
patternIDInput. Must correspond to a Pattern ID that was previously used in a call to definePattern().
connectedIdentifiersInput. List of the identifiers that form the connectivity list for the row-space (and the column-space, since this is a structurally symmetric contribution).
Returns
error-code 0 if successful

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::initConnectivity ( int  rowPatternID,
const int *  rowConnectedIdentifiers,
int  colPatternID,
const int *  colConnectedIdentifiers 
)
pure virtual

Make a contribution to the MatrixGraph's connectivity. This overloading of initConnectivity() provides for structurally non-symmetric entries.

Parameters
rowPatternIDInput. Must correspond to a Pattern ID that was previously used in a call to definePattern().
rowConnectedIdentifiersInput. List of the identifiers that form the connectivity list for the row-space.
colPatternIDInput. Must correspond to a Pattern ID that was previously used in a call to definePattern().
colConnectedIdentifiersInput. List of the identifiers that form the connectivity list for the column-space.
Returns
error-code 0 if successful

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::initConnectivity ( int  idType,
int  numRows,
const int *  rowIDs,
const int *  rowOffsets,
const int *  packedColumnIDs 
)
pure virtual

Initialize a set of arbitrary positions in the graph by providing data in a "raw" or "purely algebraic" format similar to what might be used with a standard sparse CSR (compressed sparse row) matrix.

Parameters
idTypeidentifier-type
numRowsNumber of rows, length of the following 'rowIDs' list.
rowIDsList of length 'numRows', specifying identifiers in the row-space.
rowOffsetsList of length numRows+1, giving offsets into the 'packedColumnIDs' list at which each row begins. i.e., the column IDs for rowIDs[i] are packedColumnIDs[rowOffsets[i]...rowOffsets[i+1]-1].
packedColumnIDsPacked list of length rowOffsets[numRows], containing the column IDs.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::initConnectivity ( int  idType,
int  fieldID,
int  numRows,
const int *  rowIDs,
const int *  rowOffsets,
const int *  packedColumnIDs 
)
pure virtual

Initialize a set of arbitrary positions in the graph by providing data in a "raw" or "purely algebraic" format similar to what might be used with a standard sparse CSR (compressed sparse row) matrix. Also specify a fieldID to be associated with these graph positions.

Parameters
idTypeidentifier-type
fieldIDfield-identifier
numRowsNumber of rows, length of the following 'rowIDs' list.
rowIDsList of length 'numRows', specifying identifiers in the row-space.
rowOffsetsList of length numRows+1, giving offsets into the 'packedColumnIDs' list at which each row begins. i.e., the column IDs for rowIDs[i] are packedColumnIDs[rowOffsets[i]...rowOffsets[i+1]-1].
packedColumnIDsPacked list of length rowOffsets[numRows], containing the column IDs.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::initConnectivity ( int  idType,
int  numRows,
const int *  rowIDs,
const int *  rowLengths,
const int *const *  columnIDs 
)
pure virtual

Initialize a set of arbitrary positions in the graph by providing data in a "raw" or "purely algebraic" format similar to what might be used with a standard sparse CSR (compressed sparse row) matrix.

Parameters
idTypeidentifier-type
numRowsNumber of rows, length of the following 'rowIDs' list.
rowIDsList of length 'numRows', specifying identifiers in the row-space.
rowLengthsList of length numRows, giving the number of column IDs for each row ID.
columnIDsC-style table (list of lists) containing the column IDs. Number of rows is numRows, length of i-th row is rowLengths[i].

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::initLagrangeConstraint ( int  constraintID,
int  constraintIDType,
int  numIDs,
const int *  idTypes,
const int *  IDs,
const int *  fieldIDs 
)
pure virtual

Initialize a lagrange-multiplier constraint.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::initPenaltyConstraint ( int  constraintID,
int  constraintIDType,
int  numIDs,
const int *  idTypes,
const int *  IDs,
const int *  fieldIDs 
)
pure virtual

Initialize a penalty constraint.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::initSlaveConstraint ( int  numIDs,
const int *  idTypes,
const int *  IDs,
const int *  fieldIDs,
int  offsetOfSlave,
int  offsetIntoSlaveField,
const double *  weights,
double  rhsValue 
)
pure virtual

Initialize a slave constraint. (Note to self: document the parameters.)

Implemented in fei::MatrixGraph_Impl2.

virtual bool fei::MatrixGraph::hasSlaveDof ( int  ID,
int  idType 
)
pure virtual

Query whether a given mesh object has one or more slave DOFs.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::initComplete ( )
pure virtual

Signal the MatrixGraph object that initialization is complete. At this point the MatrixGraph implementation performs internal synchronizations etc. This is a collective method.

Implemented in fei::MatrixGraph_Impl2.

virtual fei::SharedPtr<fei::SparseRowGraph> fei::MatrixGraph::createGraph ( bool  blockEntryGraph,
bool  localRowGraph_includeSharedRows = false 
)
pure virtual

Generate a sparse row-based graph from structural data that has been accumulated. Don't use this until after initComplete() has been called.

Parameters
locallyOwnedRowsThose rows of a matrix that would be owned by the local processor.
blockEntryGraphSpecifies whether the graph should be constructed on a block-entry or point-entry basis. If there is only 1 scalar DOF at each mesh-object, then a block-entry graph is the same as a point-entry graph.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::compareStructure ( const fei::MatrixGraph matrixGraph,
bool &  equivalent 
) const
pure virtual

Query whether the specified MatrixGraph is structurally equivalent to this MatrixGraph.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::getNumConnectivityBlocks ( ) const
pure virtual

Query how many connectivity blocks have been initialized.

Implemented in fei::MatrixGraph_Impl2.

virtual std::map<int,fei::ConnectivityBlock*>& fei::MatrixGraph::getConnectivityBlocks ( )
pure virtual

Query for the container of connectivity-blocks.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::getConnectivityBlockIDs ( std::vector< int > &  blockIDs) const
pure virtual

Query for the list of connectivity-block-IDs.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::getNumIDsPerConnectivityList ( int  blockID) const
pure virtual

Query how many IDs are in each connectivity list in the specified connectivity block.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::getConnectivityNumIndices ( int  blockID) const
pure virtual

Query how many scatter-indices are associated with each connectivity list for a given connectivity-block.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::getConnectivityNumIndices ( int  blockID,
int &  numRowIndices,
int &  numColIndices 
)
pure virtual

Query how many scatter-indices are associated with each connectivity list for a given connectivity-block, in both the row-dimension and the column-dimension.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::getConnectivityIndices ( int  blockID,
int  connectivityID,
int  indicesAllocLen,
int *  indices,
int &  numIndices 
)
pure virtual

Obtain the scatter-indices associated with a connectivity list.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::getConnectivityIndices ( int  blockID,
int  connectivityID,
int  rowIndicesAllocLen,
int *  rowIndices,
int &  numRowIndices,
int  colIndicesAllocLen,
int *  colIndices,
int &  numColIndices 
)
pure virtual

Obtain the scatter-indices for both the row- and column-dimension, associated with a connectivity list.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::getPatternNumIndices ( int  patternID,
int &  numIndices 
)
pure virtual

Query associated with Pattern rather than connectivity-block.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::getPatternIndices ( int  patternID,
const int *  IDs,
std::vector< int > &  indices 
)
pure virtual

Query associated with Pattern rather than connectivity-block.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::getLocalNumLagrangeConstraints ( ) const
pure virtual

Query number of local lagrange constraints

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::getGlobalNumSlaveConstraints ( ) const
pure virtual

Query number of slave-constraints

Implemented in fei::MatrixGraph_Impl2.

virtual ConstraintType* fei::MatrixGraph::getLagrangeConstraint ( int  constraintID)
pure virtual

Won't typically be of interest to application users of fei:: methods.

Implemented in fei::MatrixGraph_Impl2.

virtual std::map<int, ConstraintType* >& fei::MatrixGraph::getLagrangeConstraints ( )
pure virtual

Won't typically be of interest to application users of fei:: methods.

Implemented in fei::MatrixGraph_Impl2.

virtual ConstraintType* fei::MatrixGraph::getPenaltyConstraint ( int  constraintID)
pure virtual

Won't typically be of interest to application users of fei:: methods.

Implemented in fei::MatrixGraph_Impl2.

virtual ConstraintType* fei::MatrixGraph::getSlaveConstraint ( int  constraintID)
pure virtual

Won't typically be of interest to application users of fei:: methods.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::getConstraintConnectivityIndices ( ConstraintType cr,
std::vector< int > &  globalIndices 
)
pure virtual

Won't typically be of interest to application users of fei:: methods.

Implemented in fei::MatrixGraph_Impl2.

virtual const fei::ConnectivityBlock* fei::MatrixGraph::getConnectivityBlock ( int  blockID) const
pure virtual

Won't typically be of interest to application users of fei:: methods.

Implemented in fei::MatrixGraph_Impl2.

virtual fei::ConnectivityBlock* fei::MatrixGraph::getConnectivityBlock ( int  blockID)
pure virtual

Won't typically be of interest to application users of fei:: methods.

Implemented in fei::MatrixGraph_Impl2.

virtual void fei::MatrixGraph::setIndicesMode ( int  mode)
pure virtual

Utility method.

Implemented in fei::MatrixGraph_Impl2.

virtual fei::SharedPtr<fei::FillableMat> fei::MatrixGraph::getSlaveDependencyMatrix ( )
pure virtual

Utility method.

Implemented in fei::MatrixGraph_Impl2.

virtual fei::Pattern* fei::MatrixGraph::getPattern ( int  patternID)
pure virtual

Retrieve pointer to specified Pattern object. If specified pattern is not found, return NULL.

Implemented in fei::MatrixGraph_Impl2.

virtual int fei::MatrixGraph::createSlaveMatrices ( )
pure virtual

power-users only

Implemented in fei::MatrixGraph_Impl2.

virtual fei::SharedPtr<fei::Reducer> fei::MatrixGraph::getReducer ( )
pure virtual

Query for the object that manages the equation-space reductions associated with removing slave constraints from the linear system.

Implemented in fei::MatrixGraph_Impl2.

virtual fei::SharedPtr<fei::SparseRowGraph> fei::MatrixGraph::getRemotelyOwnedGraphRows ( )
pure virtual

query for shared-but-not-owned graph rows

Implemented in fei::MatrixGraph_Impl2.

virtual void fei::MatrixGraph::getConstrainedIndices ( std::vector< int > &  crindices) const
pure virtual

fill a vector with eqn-numbers of constrained ids

Implemented in fei::MatrixGraph_Impl2.


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