OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
Base class for all matrices stored in sparse format. More...
#include <sparsemtrx.h>
Public Types | |
typedef long | SparseMtrxVersionType |
Public Member Functions | |
SparseMtrx (int n, int m) | |
Constructor, creates (n,m) sparse matrix. More... | |
SparseMtrx () | |
Constructor. More... | |
virtual | ~SparseMtrx () |
Destructor. More... | |
SparseMtrxVersionType | giveVersion () |
Return receiver version. More... | |
void | checkBounds (int i, int j) const |
Checks size of receiver towards requested bounds. More... | |
int | giveNumberOfRows () const |
Returns number of rows of receiver. More... | |
int | giveNumberOfColumns () const |
Returns number of columns of receiver. More... | |
bool | isSquare () const |
Returns nonzero if receiver is square matrix. More... | |
bool | isNotEmpty () const |
Tests for empty matrix. More... | |
virtual SparseMtrx * | GiveCopy () const |
Returns a newly allocated copy of receiver. More... | |
virtual void | times (const FloatArray &x, FloatArray &answer) const |
Evaluates . More... | |
virtual void | timesT (const FloatArray &x, FloatArray &answer) const |
Evaluates . More... | |
virtual void | times (const FloatMatrix &B, FloatMatrix &answer) const |
Evaluates . More... | |
virtual void | timesT (const FloatMatrix &B, FloatMatrix &answer) const |
Evaluates . More... | |
virtual void | times (double x) |
Multiplies receiver by scalar value. More... | |
virtual void | add (double x, SparseMtrx &m) |
Adds x * m. More... | |
virtual void | addDiagonal (double x, FloatArray &m) |
Adds x * m (treats m as a diagonal matrix, stored as an array) More... | |
virtual int | buildInternalStructure (EngngModel *eModel, int n, int m, const IntArray &I, const IntArray &J) |
Builds internal structure of receiver based on I and J. More... | |
virtual int | buildInternalStructure (EngngModel *eModel, int di, const UnknownNumberingScheme &s)=0 |
Builds internal structure of receiver. More... | |
virtual int | buildInternalStructure (EngngModel *eModel, int di, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s) |
Build internal structure of receiver. More... | |
virtual int | assemble (const IntArray &loc, const FloatMatrix &mat)=0 |
Assembles sparse matrix from contribution of local elements. More... | |
virtual int | assemble (const IntArray &rloc, const IntArray &cloc, const FloatMatrix &mat)=0 |
Assembles sparse matrix from contribution of local elements. More... | |
virtual int | assembleBegin () |
Starts assembling the elements. More... | |
virtual int | assembleEnd () |
Returns when assemble is completed. More... | |
virtual bool | canBeFactorized () const =0 |
Determines, whether receiver can be factorized. More... | |
virtual SparseMtrx * | factorized () |
Returns the receiver factorized. More... | |
virtual FloatArray * | backSubstitutionWith (FloatArray &y) const |
Computes the solution of linear system where A is receiver. More... | |
virtual void | zero ()=0 |
Zeroes the receiver. More... | |
virtual double | computeNorm () const |
Returns the norm of receiver. More... | |
virtual SparseMtrx * | giveSubMatrix (const IntArray &rows, const IntArray &cols) |
virtual double & | at (int i, int j)=0 |
Returns coefficient at position (i,j). More... | |
virtual double | at (int i, int j) const =0 |
Returns coefficient at position (i,j). More... | |
virtual bool | isAllocatedAt (int i, int j) const |
Checks whether memory is allocated at position (i,j). More... | |
virtual void | toFloatMatrix (FloatMatrix &answer) const |
Converts receiving sparse matrix to a dense float matrix. More... | |
virtual void | printStatistics () const |
Prints the receiver statistics (one-line) to stdout. More... | |
virtual void | printYourself () const |
Prints receiver to stdout. Works only for relatively small matrices. More... | |
virtual void | writeToFile (const char *fname) const |
Helpful for debugging, writes the matrix to given file. More... | |
virtual SparseMtrxType | giveType () const =0 |
Sparse matrix type identification. More... | |
virtual bool | isAsymmetric () const =0 |
Returns true if asymmetric. More... | |
virtual const char * | giveClassName () const =0 |
std::string | errorInfo (const char *func) const |
Error printing helper. More... | |
IML compatibility | |
FloatArray | operator* (const FloatArray &x) const |
IML compatibility, . More... | |
FloatArray | trans_mult (const FloatArray &x) const |
IML compatibility, . More... | |
Protected Attributes | |
int | nRows |
Number of rows. More... | |
int | nColumns |
Number of columns. More... | |
SparseMtrxVersionType | version |
Allows to track if receiver changes. More... | |
Base class for all matrices stored in sparse format.
Basically sparse matrix contains contribution of local element matrices. Localization of local element matrix into global (structural) matrix is determined using element code numbers. Basic methods include
Definition at line 60 of file sparsemtrx.h.
typedef long oofem::SparseMtrx::SparseMtrxVersionType |
Definition at line 63 of file sparsemtrx.h.
|
inline |
Constructor, creates (n,m) sparse matrix.
Due to sparsity character of matrix, not all coefficient are physically stored (in general, zero members are omitted).
Definition at line 87 of file sparsemtrx.h.
|
inline |
Constructor.
Definition at line 89 of file sparsemtrx.h.
|
inlinevirtual |
Destructor.
Definition at line 91 of file sparsemtrx.h.
|
inlinevirtual |
Adds x * m.
x | Value to multiply m by. |
m | Matrix to add (should of the same matrix type). |
Reimplemented in oofem::Skyline, and oofem::PetscSparseMtrx.
Definition at line 166 of file sparsemtrx.h.
References OOFEM_ERROR.
|
inlinevirtual |
Adds x * m (treats m as a diagonal matrix, stored as an array)
x | Value to multiply m by. |
m | Matrix to add. |
Reimplemented in oofem::PetscSparseMtrx.
Definition at line 172 of file sparsemtrx.h.
References oofem::FloatArray::at(), and oofem::FloatArray::giveSize().
|
pure virtual |
Assembles sparse matrix from contribution of local elements.
This method for each element adds its contribution to itself. Mapping between local element contribution and its global position is given by local code numbers of element.
loc | Location array. The values corresponding to zero loc array value are not assembled. |
mat | Contribution to be assembled using loc array. |
Implemented in oofem::CompCol, oofem::Skyline, oofem::SymCompCol, oofem::DynCompCol, oofem::DynCompRow, oofem::PetscSparseMtrx, oofem::SkylineUnsym, and oofem::SpoolesSparseMtrx.
Referenced by oofem::SurfaceTensionBoundaryCondition::assemble(), oofem::PrescribedGradientBCNeumann::assemble(), oofem::PrescribedMean::assemble(), oofem::TransportGradientNeumann::assemble(), oofem::LinearConstraintBC::assemble(), oofem::MixedGradientPressureWeakPeriodic::assemble(), oofem::PrescribedGradientBCWeak::assemble(), oofem::MixedGradientPressureNeumann::assemble(), oofem::WeakPeriodicBoundaryCondition::assemble(), oofem::EngngModel::assemble(), oofem::PrescribedGradientBCWeak::assembleExtraDisplock(), oofem::PrescribedGradientBCWeak::assembleGPContrib(), oofem::ContactDefinition::computeContactTangent(), oofem::TrabBoneNL3D::NonlocalMaterialStiffnessInterface_addIPContribution(), oofem::MisesMatNl::NonlocalMaterialStiffnessInterface_addIPContribution(), oofem::RankineMatNl::NonlocalMaterialStiffnessInterface_addIPContribution(), and oofem::IDNLMaterial::NonlocalMaterialStiffnessInterface_addIPContribution().
|
pure virtual |
Assembles sparse matrix from contribution of local elements.
This method for each element adds its contribution to itself. Mapping between local element contribution and its global position is given by row and column local code numbers.
rloc | Row location array. The values corresponding to zero loc array value are not assembled. |
cloc | Column location array. The values corresponding to zero loc array value are not assembled. |
mat | Contribution to be assembled using rloc and cloc arrays. The rloc position determines the row, the cloc position determines the corresponding column. |
Implemented in oofem::CompCol, oofem::Skyline, oofem::SymCompCol, oofem::DynCompCol, oofem::DynCompRow, oofem::PetscSparseMtrx, oofem::SkylineUnsym, and oofem::SpoolesSparseMtrx.
|
inlinevirtual |
Starts assembling the elements.
Reimplemented in oofem::PetscSparseMtrx.
Definition at line 235 of file sparsemtrx.h.
Referenced by oofem::EngngModel::assemble().
|
inlinevirtual |
Returns when assemble is completed.
Reimplemented in oofem::PetscSparseMtrx.
Definition at line 237 of file sparsemtrx.h.
Referenced by oofem::EngngModel::assemble().
|
pure virtual |
Returns coefficient at position (i,j).
Implemented in oofem::Skyline, oofem::CompCol, oofem::SymCompCol, oofem::DynCompCol, oofem::DynCompRow, oofem::PetscSparseMtrx, oofem::SkylineUnsym, and oofem::SpoolesSparseMtrx.
Referenced by oofem::NRSolver::applyConstraintsToLoadIncrement(), oofem::NRSolver::applyConstraintsToStiffness(), oofem::MicroMaterial::giveMacroStiffnessMatrix(), oofem::DiagPreconditioner::init(), oofem::InverseIteration::solve(), oofem::SubspaceIteration::solve(), and oofem::FreeWarping::updateStiffnessMatrix().
|
pure virtual |
Returns coefficient at position (i,j).
Implemented in oofem::Skyline, oofem::CompCol, oofem::SymCompCol, oofem::DynCompCol, oofem::DynCompRow, oofem::PetscSparseMtrx, oofem::SkylineUnsym, and oofem::SpoolesSparseMtrx.
|
inlinevirtual |
Computes the solution of linear system where A is receiver.
Solution vector x overwrites the right hand side vector y. Receiver must be in factorized form.
y | Right hand side on input, solution on output. |
Reimplemented in oofem::Skyline, oofem::PetscSparseMtrx, oofem::SkylineUnsym, and oofem::SpoolesSparseMtrx.
Definition at line 253 of file sparsemtrx.h.
Referenced by oofem::MicroMaterial::giveMacroStiffnessMatrix(), oofem::LDLTFactorization::solve(), and oofem::SubspaceIteration::solve().
|
inlinevirtual |
Builds internal structure of receiver based on I and J.
This call is for special purpose uses. Normal problems should use the other prealloation methods.
eModel | Pointer to corresponding engineering model. |
I | Row indices |
J | Column indices |
Reimplemented in oofem::PetscSparseMtrx.
Definition at line 185 of file sparsemtrx.h.
References OOFEM_ERROR.
Referenced by oofem::MicroMaterial::giveMacroStiffnessMatrix().
|
pure virtual |
Builds internal structure of receiver.
This method determines the internal profile of sparse matrix, allocates necessary space for storing nonzero coefficients and initializes receiver. In general, the profile of sparse matrix is determined using one (or more) loop over local code numbers of elements. This method must be called before any operation, like assembly, zeroing, or multiplication.
eModel | Pointer to corresponding engineering model. |
di | Domain index specify which domain to use. |
s | Determines unknown numbering scheme. |
Implemented in oofem::CompCol, oofem::SymCompCol, oofem::Skyline, oofem::DynCompCol, oofem::DynCompRow, oofem::PetscSparseMtrx, oofem::SkylineUnsym, and oofem::SpoolesSparseMtrx.
|
inlinevirtual |
Build internal structure of receiver.
eModel | Pointer to corresponding engineering model. |
di | Domain index specify which domain to use. |
r_s | Determines unknown numbering scheme for the rows. |
c_s | Determines unknown numbering scheme for the columns. |
Reimplemented in oofem::PetscSparseMtrx.
Definition at line 208 of file sparsemtrx.h.
References OOFEM_ERROR.
|
pure virtual |
Determines, whether receiver can be factorized.
Implemented in oofem::CompCol, oofem::Skyline, oofem::SymCompCol, oofem::DynCompCol, oofem::PetscSparseMtrx, oofem::DynCompRow, oofem::SkylineUnsym, and oofem::SpoolesSparseMtrx.
Referenced by oofem::LDLTFactorization::solve().
|
inline |
Checks size of receiver towards requested bounds.
Current implementation will call exit(1), if positions are outside bounds.
i | Required number of rows. |
j | Required number of columns. |
Definition at line 102 of file sparsemtrx.h.
References OOFEM_ERROR.
|
inlinevirtual |
Returns the norm of receiver.
Reimplemented in oofem::PetscSparseMtrx.
Definition at line 258 of file sparsemtrx.h.
References OOFEM_ERROR.
|
inline |
|
inlinevirtual |
Returns the receiver factorized.
form is used.
Reimplemented in oofem::Skyline, oofem::PetscSparseMtrx, oofem::SkylineUnsym, and oofem::SpoolesSparseMtrx.
Definition at line 245 of file sparsemtrx.h.
Referenced by oofem::MicroMaterial::giveMacroStiffnessMatrix(), oofem::LDLTFactorization::solve(), and oofem::SubspaceIteration::solve().
|
pure virtual |
Implemented in oofem::Skyline, oofem::CompCol, oofem::SymCompCol, oofem::PetscSparseMtrx, oofem::DynCompCol, oofem::SkylineUnsym, and oofem::DynCompRow.
|
inlinevirtual |
Returns a newly allocated copy of receiver.
Programmer must take care about proper deallocation of allocated space.
Reimplemented in oofem::CompCol, oofem::SymCompCol, oofem::Skyline, oofem::DynCompCol, oofem::DynCompRow, oofem::SkylineUnsym, oofem::SpoolesSparseMtrx, and oofem::PetscSparseMtrx.
Definition at line 127 of file sparsemtrx.h.
References OOFEM_ERROR.
|
inline |
Returns number of columns of receiver.
Definition at line 116 of file sparsemtrx.h.
Referenced by oofem::MicroMaterial::giveMacroStiffnessMatrix(), oofem::CompCol_ILUPreconditioner::initialize(), oofem::SLEPcSolver::solve(), oofem::SuperLUSolver::solve(), oofem::InverseIteration::solve(), oofem::SubspaceIteration::solve(), oofem::PetscSparseMtrx::times(), oofem::SpoolesSparseMtrx::times(), oofem::PetscSparseMtrx::timesT(), oofem::PetscSparseMtrx::toFloatMatrix(), and oofem::Skyline::toFloatMatrix().
|
inline |
Returns number of rows of receiver.
Definition at line 114 of file sparsemtrx.h.
Referenced by oofem::Skyline::at(), oofem::Skyline::backSubstitutionWith(), oofem::PetscSparseMtrx::createVecGlobal(), oofem::Skyline::factorized(), oofem::Skyline::GiveCopy(), oofem::DiagPreconditioner::init(), oofem::CompCol_ILUPreconditioner::initialize(), oofem::Skyline::ldl_feti_sky(), oofem::Skyline::rbmodes(), oofem::PetscSparseMtrx::scatterG2L(), oofem::SLEPcSolver::solve(), oofem::PardisoProjectOrgSolver::solve(), oofem::SuperLUSolver::solve(), oofem::SparseLinearSystemNM::solve(), oofem::SpoolesSolver::solve(), oofem::FETISolver::solve(), oofem::PetscSparseMtrx::times(), oofem::Skyline::times(), oofem::PetscSparseMtrx::timesT(), oofem::SpoolesSparseMtrx::timesT(), oofem::PetscSparseMtrx::toFloatMatrix(), and oofem::Skyline::~Skyline().
|
inlinevirtual |
Reimplemented in oofem::Skyline, and oofem::PetscSparseMtrx.
Definition at line 263 of file sparsemtrx.h.
References OOFEM_ERROR.
Referenced by oofem::StaggeredSolver::solve().
|
pure virtual |
Sparse matrix type identification.
Implemented in oofem::Skyline, oofem::CompCol, oofem::SymCompCol, oofem::PetscSparseMtrx, oofem::DynCompCol, oofem::SkylineUnsym, oofem::SpoolesSparseMtrx, and oofem::DynCompRow.
|
inline |
Return receiver version.
Definition at line 94 of file sparsemtrx.h.
Referenced by oofem::NRSolver::applyConstraintsToStiffness(), oofem::SpoolesSolver::solve(), and oofem::IMLSolver::solve().
|
inlinevirtual |
Checks whether memory is allocated at position (i,j).
Reimplemented in oofem::Skyline.
Definition at line 271 of file sparsemtrx.h.
Referenced by oofem::MicroMaterial::giveMacroStiffnessMatrix().
|
pure virtual |
Returns true if asymmetric.
Implemented in oofem::Skyline, oofem::CompCol, oofem::PetscSparseMtrx, oofem::DynCompCol, oofem::SkylineUnsym, oofem::SpoolesSparseMtrx, and oofem::DynCompRow.
|
inline |
Tests for empty matrix.
Definition at line 120 of file sparsemtrx.h.
|
inline |
Returns nonzero if receiver is square matrix.
Definition at line 118 of file sparsemtrx.h.
|
inline |
IML compatibility, .
Definition at line 291 of file sparsemtrx.h.
|
inlinevirtual |
Prints the receiver statistics (one-line) to stdout.
Reimplemented in oofem::DynCompCol, oofem::PetscSparseMtrx, oofem::DynCompRow, oofem::SkylineUnsym, and oofem::SpoolesSparseMtrx.
Definition at line 275 of file sparsemtrx.h.
References OOFEM_LOG_INFO.
Referenced by oofem::NonLinearDynamic::assemble(), and oofem::NonLinearStatic::assemble().
|
inlinevirtual |
Prints receiver to stdout. Works only for relatively small matrices.
Reimplemented in oofem::Skyline, oofem::CompCol, oofem::PetscSparseMtrx, oofem::SkylineUnsym, and oofem::SpoolesSparseMtrx.
Definition at line 277 of file sparsemtrx.h.
References OOFEM_LOG_INFO.
|
inlinevirtual |
Evaluates .
x | Array to be multiplied with receiver. |
answer | y. |
Reimplemented in oofem::CompCol, oofem::SymCompCol, oofem::Skyline, oofem::DynCompCol, oofem::DynCompRow, oofem::SkylineUnsym, oofem::SpoolesSparseMtrx, and oofem::PetscSparseMtrx.
Definition at line 137 of file sparsemtrx.h.
References OOFEM_ERROR.
Referenced by oofem::InverseIteration::solve(), oofem::SubspaceIteration::solve(), and oofem::SUPG::solveYourselfAt().
|
inlinevirtual |
Evaluates .
B | Array to be multiplied with receiver. |
answer | C. |
Reimplemented in oofem::PetscSparseMtrx.
Definition at line 149 of file sparsemtrx.h.
References OOFEM_ERROR.
|
inlinevirtual |
Multiplies receiver by scalar value.
x | Value to multiply receiver. |
Reimplemented in oofem::CompCol, oofem::SymCompCol, oofem::Skyline, oofem::DynCompCol, oofem::DynCompRow, oofem::SkylineUnsym, oofem::SpoolesSparseMtrx, and oofem::PetscSparseMtrx.
Definition at line 160 of file sparsemtrx.h.
References OOFEM_ERROR.
|
inlinevirtual |
Evaluates .
x | Array to be multiplied with transpose of the receiver. |
answer | y. |
Reimplemented in oofem::CompCol, oofem::SymCompCol, oofem::Skyline, oofem::DynCompCol, oofem::DynCompRow, oofem::SkylineUnsym, oofem::SpoolesSparseMtrx, and oofem::PetscSparseMtrx.
Definition at line 143 of file sparsemtrx.h.
References OOFEM_ERROR.
|
inlinevirtual |
Evaluates .
B | Matrix to be multiplied with receiver. |
answer | C. |
Reimplemented in oofem::PetscSparseMtrx.
Definition at line 155 of file sparsemtrx.h.
References OOFEM_ERROR.
|
inlinevirtual |
Converts receiving sparse matrix to a dense float matrix.
Reimplemented in oofem::Skyline, oofem::CompCol, oofem::PetscSparseMtrx, and oofem::SkylineUnsym.
Definition at line 273 of file sparsemtrx.h.
References OOFEM_ERROR.
|
inline |
IML compatibility, .
Definition at line 298 of file sparsemtrx.h.
|
inlinevirtual |
Helpful for debugging, writes the matrix to given file.
Reimplemented in oofem::Skyline, oofem::PetscSparseMtrx, and oofem::SkylineUnsym.
Definition at line 279 of file sparsemtrx.h.
References OOFEM_LOG_INFO.
|
pure virtual |
Zeroes the receiver.
Implemented in oofem::CompCol, oofem::Skyline, oofem::SymCompCol, oofem::DynCompCol, oofem::PetscSparseMtrx, oofem::DynCompRow, oofem::SkylineUnsym, and oofem::SpoolesSparseMtrx.
Referenced by oofem::MicroMaterial::giveMacroStiffnessMatrix().
|
protected |
Number of columns.
Definition at line 69 of file sparsemtrx.h.
Referenced by oofem::DynCompRow::at(), oofem::PetscSparseMtrx::buildInternalStructure(), oofem::SkylineUnsym::buildInternalStructure(), oofem::SpoolesSparseMtrx::buildInternalStructure(), oofem::DynCompRow::buildInternalStructure(), oofem::DynCompCol::buildInternalStructure(), oofem::Skyline::buildInternalStructure(), oofem::SymCompCol::buildInternalStructure(), oofem::CompCol::buildInternalStructure(), oofem::CompCol::CompCol(), oofem::DynCompCol::DynCompCol(), oofem::DynCompRow::DynCompRow(), oofem::PetscSparseMtrx::GiveCopy(), oofem::DynCompRow::growTo(), oofem::DynCompCol::growTo(), oofem::DynCompRow::ILUPYourself(), oofem::DynCompRow::operator()(), oofem::DynCompRow::operator=(), oofem::DynCompCol::operator=(), oofem::CompCol::operator=(), oofem::DynCompCol::printStatistics(), oofem::SkylineUnsym::setInternalStructure(), oofem::Skyline::setInternalStructure(), oofem::DynCompRow::times(), oofem::DynCompCol::times(), oofem::PetscSparseMtrx::timesT(), oofem::DynCompRow::timesT(), oofem::DynCompCol::timesT(), oofem::SkylineUnsym::writeToFile(), oofem::Skyline::writeToFile(), oofem::DynCompCol::zero(), and oofem::DynCompCol::~DynCompCol().
|
protected |
Number of rows.
Definition at line 67 of file sparsemtrx.h.
Referenced by oofem::SymCompCol::assemble(), oofem::CompCol::assemble(), oofem::DynCompRow::at(), oofem::DynCompCol::at(), oofem::SpoolesSparseMtrx::buildInternalStructure(), oofem::PetscSparseMtrx::buildInternalStructure(), oofem::SkylineUnsym::buildInternalStructure(), oofem::DynCompRow::buildInternalStructure(), oofem::DynCompCol::buildInternalStructure(), oofem::Skyline::buildInternalStructure(), oofem::SymCompCol::buildInternalStructure(), oofem::CompCol::buildInternalStructure(), oofem::CompCol::CompCol(), oofem::DynCompCol::DynCompCol(), oofem::DynCompRow::DynCompRow(), oofem::PetscSparseMtrx::GiveCopy(), oofem::DynCompRow::growTo(), oofem::DynCompCol::growTo(), oofem::DynCompRow::ILUPYourself(), oofem::DynCompRow::operator()(), oofem::DynCompCol::operator()(), oofem::DynCompRow::operator=(), oofem::DynCompCol::operator=(), oofem::CompCol::operator=(), oofem::DynCompRow::printStatistics(), oofem::SkylineUnsym::setInternalStructure(), oofem::Skyline::setInternalStructure(), oofem::PetscSparseMtrx::times(), oofem::DynCompRow::times(), oofem::DynCompCol::times(), oofem::DynCompRow::timesT(), oofem::DynCompCol::timesT(), oofem::SkylineUnsym::writeToFile(), oofem::Skyline::writeToFile(), oofem::DynCompRow::zero(), and oofem::DynCompRow::~DynCompRow().
|
protected |
Allows to track if receiver changes.
Any change on receiver should increment this variable. The factorization should not change version (nrsolver direct control relies on that). This version info is used for example by preconditioners associated to particular matrix; the preconditioner initialization can be demanding and this versioning allows to reuse initialized preconditioner for same matrix, if there is no change;
Definition at line 80 of file sparsemtrx.h.
Referenced by oofem::Skyline::add(), oofem::SpoolesSparseMtrx::assemble(), oofem::SkylineUnsym::assemble(), oofem::PetscSparseMtrx::assemble(), oofem::DynCompRow::assemble(), oofem::DynCompCol::assemble(), oofem::SymCompCol::assemble(), oofem::Skyline::assemble(), oofem::CompCol::assemble(), oofem::SkylineUnsym::at(), oofem::DynCompRow::at(), oofem::DynCompCol::at(), oofem::SymCompCol::at(), oofem::CompCol::at(), oofem::Skyline::at(), oofem::SkylineUnsym::buildInternalStructure(), oofem::DynCompRow::buildInternalStructure(), oofem::DynCompCol::buildInternalStructure(), oofem::Skyline::buildInternalStructure(), oofem::SymCompCol::buildInternalStructure(), oofem::CompCol::buildInternalStructure(), oofem::DynCompCol::DynCompCol(), oofem::DynCompRow::DynCompRow(), oofem::DynCompRow::ILUPYourself(), oofem::DynCompRow::operator()(), oofem::DynCompCol::operator()(), oofem::SymCompCol::operator()(), oofem::CompCol::operator()(), oofem::DynCompRow::operator=(), oofem::DynCompCol::operator=(), oofem::CompCol::operator=(), oofem::SkylineUnsym::setInternalStructure(), oofem::Skyline::setInternalStructure(), oofem::SkylineUnsym::times(), oofem::DynCompRow::times(), oofem::DynCompCol::times(), oofem::Skyline::times(), oofem::SymCompCol::times(), oofem::CompCol::times(), oofem::SkylineUnsym::zero(), oofem::DynCompRow::zero(), oofem::DynCompCol::zero(), oofem::SymCompCol::zero(), oofem::Skyline::zero(), and oofem::CompCol::zero().