FEI Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Public Member Functions | Protected Attributes | List of all members
fei::LinearSystem Class Referenceabstract

#include <fei_LinearSystem.hpp>

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

Classes

class  Factory
 

Public Member Functions

 LinearSystem (fei::SharedPtr< fei::MatrixGraph > &matrixGraph)
 
virtual ~LinearSystem ()
 
virtual int parameters (int numParams, const char *const *paramStrings)=0
 
virtual int parameters (const fei::ParameterSet &params)=0
 
virtual void setMatrix (fei::SharedPtr< fei::Matrix > &matrix)
 
virtual fei::SharedPtr
< fei::Matrix
getMatrix ()
 
virtual fei::SharedPtr< const
fei::Matrix
getMatrix () const
 
virtual void setRHS (fei::SharedPtr< fei::Vector > &rhs)
 
virtual fei::SharedPtr
< fei::Vector
getRHS ()
 
virtual fei::SharedPtr< const
fei::Vector
getRHS () const
 
virtual void setSolutionVector (fei::SharedPtr< fei::Vector > &soln)
 
virtual fei::SharedPtr
< fei::Vector
getSolutionVector ()
 
virtual fei::SharedPtr< const
fei::Vector
getSolutionVector () const
 
virtual int putAttribute (const char *name, void *attribute)
 
virtual int getAttribute (const char *name, void *&attribute)
 
virtual int loadEssentialBCs (int numIDs, const int *IDs, int idType, int fieldID, int offsetIntoField, const double *prescribedValues)
 
virtual int loadEssentialBCs (int numIDs, const int *IDs, int idType, int fieldID, const int *offsetsIntoField, const double *prescribedValues)
 
virtual int loadLagrangeConstraint (int constraintID, const double *weights, double rhsValue)=0
 
virtual int loadPenaltyConstraint (int constraintID, const double *weights, double penaltyValue, double rhsValue)=0
 
virtual int loadComplete (bool applyBCs=true, bool globalAssemble=true)=0
 
virtual int setBCValuesOnVector (fei::Vector *vector)=0
 
virtual bool eqnIsEssentialBC (int globalEqnIndex) const =0
 
virtual void getEssentialBCs (std::vector< int > &bcEqns, std::vector< double > &bcVals) const =0
 
virtual void getConstrainedEqns (std::vector< int > &crEqns) const =0
 

Protected Attributes

fei::SharedPtr< fei::Matrixmatrix_
 
fei::SharedPtr< fei::Vectorsoln_
 
fei::SharedPtr< fei::Vectorrhs_
 
fei::SharedPtr< fei::MatrixGraphmatrixGraph_
 
fei::DirichletBCManagerdbcManager_
 
std::vector< char * > attributeNames_
 
std::vector< void * > attributes_
 

Detailed Description

A simple container to bind a matrix and two vectors together as the matrix, rhs and solution of a linear system.

Definition at line 26 of file fei_LinearSystem.hpp.

Constructor & Destructor Documentation

fei::LinearSystem::LinearSystem ( fei::SharedPtr< fei::MatrixGraph > &  matrixGraph)

Constructor

Definition at line 17 of file fei_LinearSystem.cpp.

fei::LinearSystem::~LinearSystem ( )
virtual

Destructor

Definition at line 27 of file fei_LinearSystem.cpp.

Member Function Documentation

virtual int fei::LinearSystem::parameters ( int  numParams,
const char *const *  paramStrings 
)
pure virtual
  Set parameters on this 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 snl_fei::LinearSystem_General, and snl_fei::LinearSystem_FEData.

Referenced by beam_main(), main(), poisson3_main(), snl_fei_tester::setParameter(), test_LinearSystem::test3(), and snl_fei_tester::testLoading().

virtual int fei::LinearSystem::parameters ( const fei::ParameterSet params)
pure virtual

Set parameters on this object.

Implemented in snl_fei::LinearSystem_General, and snl_fei::LinearSystem_FEData.

void fei::LinearSystem::setMatrix ( fei::SharedPtr< fei::Matrix > &  matrix)
virtual
virtual fei::SharedPtr<fei::Matrix> fei::LinearSystem::getMatrix ( )
inlinevirtual

Get the matrix for this linear system.

Definition at line 62 of file fei_LinearSystem.hpp.

References matrix_.

Referenced by fei_Solver_solve().

virtual fei::SharedPtr<const fei::Matrix> fei::LinearSystem::getMatrix ( ) const
inlinevirtual

Get the matrix for this linear system.

Definition at line 66 of file fei_LinearSystem.hpp.

References matrix_.

virtual void fei::LinearSystem::setRHS ( fei::SharedPtr< fei::Vector > &  rhs)
inlinevirtual

Set the right-hand-side for this linear system.

Definition at line 71 of file fei_LinearSystem.hpp.

References rhs_.

Referenced by beam_main(), main(), poisson3_main(), test_LinearSystem::test2(), test_LinearSystem::test3(), test_LinearSystem::test4(), and snl_fei_tester::testLoading().

virtual fei::SharedPtr<fei::Vector> fei::LinearSystem::getRHS ( )
inlinevirtual

Get the right-hand-side for this linear system.

Definition at line 75 of file fei_LinearSystem.hpp.

References rhs_.

virtual fei::SharedPtr<const fei::Vector> fei::LinearSystem::getRHS ( ) const
inlinevirtual

Get the right-hand-side for this linear system.

Definition at line 79 of file fei_LinearSystem.hpp.

References rhs_.

virtual void fei::LinearSystem::setSolutionVector ( fei::SharedPtr< fei::Vector > &  soln)
inlinevirtual

Set the solution for this linear system.

Definition at line 84 of file fei_LinearSystem.hpp.

References soln_.

Referenced by beam_main(), main(), poisson3_main(), test_LinearSystem::test2(), test_LinearSystem::test3(), test_LinearSystem::test4(), and snl_fei_tester::testLoading().

virtual fei::SharedPtr<fei::Vector> fei::LinearSystem::getSolutionVector ( )
inlinevirtual

Get the solution for this linear system.

Definition at line 88 of file fei_LinearSystem.hpp.

References soln_.

virtual fei::SharedPtr<const fei::Vector> fei::LinearSystem::getSolutionVector ( ) const
inlinevirtual

Get the solution for this linear system.

Definition at line 92 of file fei_LinearSystem.hpp.

References soln_.

int fei::LinearSystem::putAttribute ( const char *  name,
void *  attribute 
)
virtual

Store an attribute on this LinearSystem which can be retrieved later by name.

Definition at line 53 of file fei_LinearSystem.cpp.

References snl_fei::storeNamedAttribute().

int fei::LinearSystem::getAttribute ( const char *  name,
void *&  attribute 
)
virtual

Retrieve an attribute which was previously stored.

Definition at line 62 of file fei_LinearSystem.cpp.

References snl_fei::retrieveNamedAttribute().

int fei::LinearSystem::loadEssentialBCs ( int  numIDs,
const int *  IDs,
int  idType,
int  fieldID,
int  offsetIntoField,
const double *  prescribedValues 
)
virtual

Essential boundary-condition function that simply accepts a list of prescribed values, rather than the 'old' FEI's confusing approach of accepting arrays of alpha, beta and gamma values that nobody every really understood.

For each specified ID, a value is being prescribed for a specified fieldID and a specified offset into that field.

Parameters
numIDs
IDs
idType
fieldID
offsetIntoField
prescribedValuesInput. List of values. Has length numIDs.

Reimplemented in snl_fei::LinearSystem_General.

Definition at line 70 of file fei_LinearSystem.cpp.

References fei::DirichletBCManager::addBCRecords(), fei::console_out(), and FEI_ENDL.

Referenced by load_BC_data(), HexBeam_Functions::load_BC_data(), snl_fei::LinearSystem_General::loadEssentialBCs(), and snl_fei_tester::testLoading().

int fei::LinearSystem::loadEssentialBCs ( int  numIDs,
const int *  IDs,
int  idType,
int  fieldID,
const int *  offsetsIntoField,
const double *  prescribedValues 
)
virtual

Essential boundary-condition function that simply accepts a list of prescribed values, rather than the 'old' FEI's confusing approach of accepting arrays of alpha, beta and gamma values that nobody every really understood.

For each specified ID, a value is being prescribed for a specified fieldID and a specified offset into that field. The offset into the field can be different for each prescribed value.

Parameters
numIDs
IDs
idType
fieldID
offsetsIntoFieldInput. List of values, length numIDs.
prescribedValuesInput. List of values. Has length numIDs.

Reimplemented in snl_fei::LinearSystem_General.

Definition at line 94 of file fei_LinearSystem.cpp.

References fei::DirichletBCManager::addBCRecords(), fei::console_out(), and FEI_ENDL.

virtual int fei::LinearSystem::loadLagrangeConstraint ( int  constraintID,
const double *  weights,
double  rhsValue 
)
pure virtual
  Lagrange constraint coefficient loading function.
Parameters
constraintIDInput. Must be an identifier of a lagrange constraint that was initialized on the fei::MatrixGraph object which was used to construct the matrix for this linear system.
weightsInput. List, with length given by the sum of the sizes of the constrained fields.
rhsValue

Implemented in snl_fei::LinearSystem_General, and snl_fei::LinearSystem_FEData.

Referenced by HexBeam_Functions::load_constraints(), snl_fei_tester::loadConstraints(), test_LinearSystem::test2(), and test_LinearSystem::test4().

virtual int fei::LinearSystem::loadPenaltyConstraint ( int  constraintID,
const double *  weights,
double  penaltyValue,
double  rhsValue 
)
pure virtual
  Penalty constraint coefficient loading function.
Parameters
constraintIDInput. Must be an identifier of a lagrange constraint that was initialized on the fei::MatrixGraph object which was used to construct the matrix for this linear system.
weightsInput. List, with length given by the sum of the sizes of the constrained fields.
penaltyValue
rhsValue

Implemented in snl_fei::LinearSystem_General, and snl_fei::LinearSystem_FEData.

Referenced by snl_fei_tester::loadConstraints(), and test_LinearSystem::test3().

virtual int fei::LinearSystem::loadComplete ( bool  applyBCs = true,
bool  globalAssemble = true 
)
pure virtual
  Signal that all boundary-conditions and constraint coefficients have been

loaded, and they may now be applied to the linear system.

Implemented in snl_fei::LinearSystem_General, and snl_fei::LinearSystem_FEData.

Referenced by beam_main(), main(), poisson3_main(), test_LinearSystem::test2(), test_LinearSystem::test3(), test_LinearSystem::test4(), and snl_fei_tester::testLoading().

virtual int fei::LinearSystem::setBCValuesOnVector ( fei::Vector vector)
pure virtual
  Request that any boundary-condition values that have been provided via

loadEssentialBCs() be set in the specified vector.

Implemented in snl_fei::LinearSystem_General, and snl_fei::LinearSystem_FEData.

virtual bool fei::LinearSystem::eqnIsEssentialBC ( int  globalEqnIndex) const
pure virtual
  Query whether a specified equation-index has a prescribed

essential boundary-condition.

Implemented in snl_fei::LinearSystem_FEData, and snl_fei::LinearSystem_General.

virtual void fei::LinearSystem::getEssentialBCs ( std::vector< int > &  bcEqns,
std::vector< double > &  bcVals 
) const
pure virtual
  Fill caller-supplied vectors with the global equation-indices (which

reside on the local processor) that have essential boundary-conditions prescribed, and fill a second vector with the prescribed values.

Implemented in snl_fei::LinearSystem_FEData, and snl_fei::LinearSystem_General.

virtual void fei::LinearSystem::getConstrainedEqns ( std::vector< int > &  crEqns) const
pure virtual

Fill a caller-supplied vector with the global equation-indices (which reside on the local processor) that are involved in constraint-relations.

Implemented in snl_fei::LinearSystem_FEData, and snl_fei::LinearSystem_General.

Referenced by test_LinearSystem::test2().

Member Data Documentation

fei::SharedPtr<fei::Matrix> fei::LinearSystem::matrix_
protected

Definition at line 205 of file fei_LinearSystem.hpp.

Referenced by getMatrix().

fei::SharedPtr<fei::Vector> fei::LinearSystem::soln_
protected

Definition at line 206 of file fei_LinearSystem.hpp.

Referenced by getSolutionVector(), and setSolutionVector().

fei::SharedPtr<fei::Vector> fei::LinearSystem::rhs_
protected

Definition at line 207 of file fei_LinearSystem.hpp.

Referenced by getRHS(), and setRHS().

fei::SharedPtr<fei::MatrixGraph> fei::LinearSystem::matrixGraph_
protected

Definition at line 209 of file fei_LinearSystem.hpp.

fei::DirichletBCManager* fei::LinearSystem::dbcManager_
protected

Definition at line 210 of file fei_LinearSystem.hpp.

std::vector<char*> fei::LinearSystem::attributeNames_
protected

Definition at line 212 of file fei_LinearSystem.hpp.

std::vector<void*> fei::LinearSystem::attributes_
protected

Definition at line 213 of file fei_LinearSystem.hpp.


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