OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
Class implementing sparse matrix stored in skyline form. More...
#include <skyline.h>
Public Member Functions | |
Skyline (int n) | |
Constructor. More... | |
Skyline () | |
Constructor. More... | |
virtual | ~Skyline () |
Destructor. 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 (double x) |
Multiplies receiver by scalar value. More... | |
virtual void | add (double x, SparseMtrx &m) |
Adds x * m. More... | |
virtual int | buildInternalStructure (EngngModel *, int, const UnknownNumberingScheme &) |
Builds internal structure of receiver. More... | |
virtual SparseMtrx * | giveSubMatrix (const IntArray &rows, const IntArray &cols) |
int | setInternalStructure (IntArray &a) |
Allocates and builds internal structure according to given array holding addresses of diagonal members values (adr). More... | |
virtual int | assemble (const IntArray &loc, const FloatMatrix &mat) |
Assembles sparse matrix from contribution of local elements. More... | |
virtual int | assemble (const IntArray &rloc, const IntArray &cloc, const FloatMatrix &mat) |
Assembles sparse matrix from contribution of local elements. More... | |
virtual bool | canBeFactorized () const |
Determines, whether receiver can be factorized. More... | |
virtual SparseMtrx * | factorized () |
Returns the receiver factorized. More... | |
virtual FloatArray * | backSubstitutionWith (FloatArray &) const |
Computes the solution of linear system where A is receiver. More... | |
virtual void | zero () |
Zeroes the receiver. More... | |
void | rbmodes (FloatMatrix &r, int &nse, IntArray &se, double limit, int tc) |
Splits the receiver to LDLT form, and computes the rigid body motions. More... | |
void | ldl_feti_sky (FloatArray &x, FloatArray &y, int nse, double limit, IntArray &se) |
Solves the singular system of equations, the receiver should be factorized using rbmodes service. More... | |
virtual double & | at (int, int) |
Returns coefficient at position (i,j). More... | |
virtual double | at (int i, int j) const |
Returns coefficient at position (i,j). More... | |
virtual bool | isAllocatedAt (int i, int j) const |
Returns 0 if the memory is not allocated at position (i,j). More... | |
int | giveNumberOfNonZeros () const |
virtual void | toFloatMatrix (FloatMatrix &answer) const |
Converts receiving sparse matrix to a dense float matrix. More... | |
virtual void | printYourself () const |
Prints receiver to stdout. More... | |
virtual void | writeToFile (const char *fname) const |
Helpful for debugging, writes the matrix to given file. More... | |
int | giveAllocatedSize () |
virtual SparseMtrxType | giveType () const |
Sparse matrix type identification. More... | |
virtual bool | isAsymmetric () const |
Returns true if asymmetric. More... | |
virtual const char * | giveClassName () const |
Public Member Functions inherited from oofem::SparseMtrx | |
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 void | times (const FloatMatrix &B, FloatMatrix &answer) const |
Evaluates . More... | |
virtual void | timesT (const FloatMatrix &B, FloatMatrix &answer) const |
Evaluates . 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 &r_s, const UnknownNumberingScheme &c_s) |
Build internal structure of receiver. More... | |
virtual int | assembleBegin () |
Starts assembling the elements. More... | |
virtual int | assembleEnd () |
Returns when assemble is completed. More... | |
virtual double | computeNorm () const |
Returns the norm of receiver. More... | |
virtual void | printStatistics () const |
Prints the receiver statistics (one-line) to stdout. More... | |
std::string | errorInfo (const char *func) const |
Error printing helper. More... | |
FloatArray | operator* (const FloatArray &x) const |
IML compatibility, . More... | |
FloatArray | trans_mult (const FloatArray &x) const |
IML compatibility, . More... | |
Protected Member Functions | |
Skyline (int, int, double *, const IntArray &) | |
Protected Attributes | |
int | nwk |
Total number of nonzero coefficients stored. More... | |
IntArray | adr |
Integer array holding addresses of diagonal members. More... | |
double * | mtrx |
Vector of stored coefficients. More... | |
int | isFactorized |
Flag indicating whether factorized. More... | |
Protected Attributes inherited from oofem::SparseMtrx | |
int | nRows |
Number of rows. More... | |
int | nColumns |
Number of columns. More... | |
SparseMtrxVersionType | version |
Allows to track if receiver changes. More... | |
Additional Inherited Members | |
Public Types inherited from oofem::SparseMtrx | |
typedef long | SparseMtrxVersionType |
Class implementing sparse matrix stored in skyline form.
This class assumes symmetric form of matrix to be stored, i.e., only upper half part is therefore stored. Coefficients are stored in one dimensional float array, containing for each column the values from the diagonal to the last nonzero entry. This requires to remember the addresses of diagonal members in such array (adr attribute)
Attribute 'isFactorized' is True if the skyline is already in U(T).D.U factorized form, else it is False. Attribute adr is array of diagonal members, it's size is size+1 (adr[0]=neni, adr[1]=1) Attribute mtrx is double pointer to skyline stored in a array form (but we start from index 1)
Tasks:
oofem::Skyline::Skyline | ( | int | n | ) |
Constructor.
Before any operation an internal profile must be built.
Definition at line 62 of file skyline.C.
References isFactorized, mtrx, and nwk.
oofem::Skyline::Skyline | ( | ) |
Constructor.
Before any operation an internal profile must be built.
Definition at line 73 of file skyline.C.
References isFactorized, mtrx, and nwk.
Referenced by GiveCopy(), and giveSubMatrix().
|
virtual |
Destructor.
Definition at line 83 of file skyline.C.
References oofem::SparseMtrx::giveNumberOfRows(), and mtrx.
|
protected |
|
virtual |
Adds x * m.
x | Value to multiply m by. |
m | Matrix to add (should of the same matrix type). |
Reimplemented from oofem::SparseMtrx.
Definition at line 633 of file skyline.C.
References mtrx, nwk, and oofem::SparseMtrx::version.
|
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. |
Implements oofem::SparseMtrx.
Definition at line 224 of file skyline.C.
References adr, oofem::IntArray::at(), oofem::FloatMatrix::at(), oofem::FloatMatrix::giveNumberOfRows(), oofem::IntArray::giveSize(), mtrx, OOFEM_ERROR, and oofem::SparseMtrx::version.
|
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. |
Implements oofem::SparseMtrx.
Definition at line 266 of file skyline.C.
References oofem::IntArray::at(), at(), oofem::FloatMatrix::at(), oofem::FloatMatrix::giveNumberOfColumns(), oofem::FloatMatrix::giveNumberOfRows(), and oofem::SparseMtrx::version.
|
virtual |
Returns coefficient at position (i,j).
Implements oofem::SparseMtrx.
Definition at line 93 of file skyline.C.
References adr, oofem::IntArray::at(), oofem::SparseMtrx::giveNumberOfRows(), mtrx, OOFEM_ERROR, and oofem::SparseMtrx::version.
Referenced by assemble(), and giveSubMatrix().
|
virtual |
Returns coefficient at position (i,j).
Implements oofem::SparseMtrx.
Definition at line 134 of file skyline.C.
References adr, oofem::IntArray::at(), oofem::SparseMtrx::giveNumberOfRows(), mtrx, and OOFEM_ERROR.
|
virtual |
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 from oofem::SparseMtrx.
Definition at line 289 of file skyline.C.
References adr, oofem::IntArray::at(), oofem::SparseMtrx::giveNumberOfRows(), and mtrx.
|
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. |
Implements oofem::SparseMtrx.
Definition at line 362 of file skyline.C.
References adr, oofem::IntArray::at(), oofem::IntArray::clear(), oofem::Domain::giveBcs(), oofem::ContactManager::giveContactDefinition(), oofem::ContactDefinition::giveContactElement(), oofem::Domain::giveContactManager(), oofem::EngngModel::giveDomain(), oofem::Domain::giveElements(), oofem::ContactElement::giveLocationArray(), oofem::ActiveBoundaryCondition::giveLocationArrays(), oofem::ContactManager::giveNumberOfContactDefinitions(), oofem::EngngModel::giveNumberOfDomainEquations(), oofem::ContactDefinition::giveNumbertOfContactElements(), oofem::UnknownNumberingScheme::giveRequiredNumberOfDomainEquation(), oofem::Domain::hasContactManager(), oofem::UnknownNumberingScheme::isDefault(), oofem::min(), mtrx, oofem::SparseMtrx::nColumns, oofem::SparseMtrx::nRows, nwk, OOFEM_ERROR, oofem::IntArray::resize(), and oofem::SparseMtrx::version.
|
inlinevirtual |
Determines, whether receiver can be factorized.
Implements oofem::SparseMtrx.
|
virtual |
Returns the receiver factorized.
form is used.
Reimplemented from oofem::SparseMtrx.
Definition at line 499 of file skyline.C.
References adr, oofem::IntArray::at(), oofem::Timer::getUtime(), oofem::SparseMtrx::giveNumberOfRows(), isFactorized, mtrx, nwk, OOFEM_LOG_DEBUG, oofem::Timer::startTimer(), and oofem::Timer::stopTimer().
|
inlinevirtual |
Implements oofem::SparseMtrx.
|
virtual |
Returns a newly allocated copy of receiver.
Programmer must take care about proper deallocation of allocated space.
Reimplemented from oofem::SparseMtrx.
Definition at line 682 of file skyline.C.
References adr, oofem::SparseMtrx::giveNumberOfRows(), mtrx, nwk, OOFEM_ERROR, and Skyline().
|
inline |
Definition at line 148 of file skyline.h.
Referenced by giveSubMatrix().
|
virtual |
Reimplemented from oofem::SparseMtrx.
Definition at line 720 of file skyline.C.
References oofem::IntArray::at(), at(), giveNumberOfNonZeros(), oofem::IntArray::giveSize(), isAllocatedAt(), OOFEM_ERROR, and Skyline().
|
inlinevirtual |
Sparse matrix type identification.
Implements oofem::SparseMtrx.
Definition at line 155 of file skyline.h.
References oofem::SMT_Skyline.
|
virtual |
Returns 0 if the memory is not allocated at position (i,j).
Reimplemented from oofem::SparseMtrx.
Definition at line 174 of file skyline.C.
References adr, and oofem::IntArray::at().
Referenced by giveSubMatrix().
|
inlinevirtual |
void oofem::Skyline::ldl_feti_sky | ( | FloatArray & | x, |
FloatArray & | y, | ||
int | nse, | ||
double | limit, | ||
IntArray & | se | ||
) |
Solves the singular system of equations, the receiver should be factorized using rbmodes service.
x | solution vector |
y | right hand side |
nse | number of rigid body motions |
limit | determines the linear dependence or independence |
se | indexes of singular equations |
Definition at line 932 of file skyline.C.
References adr, oofem::IntArray::at(), oofem::FloatArray::at(), oofem::SparseMtrx::giveNumberOfRows(), and mtrx.
Referenced by oofem::FETISolver::solve().
|
virtual |
Prints receiver to stdout.
Reimplemented from oofem::SparseMtrx.
Definition at line 644 of file skyline.C.
References oofem::FloatMatrix::printYourself(), and toFloatMatrix().
void oofem::Skyline::rbmodes | ( | FloatMatrix & | r, |
int & | nse, | ||
IntArray & | se, | ||
double | limit, | ||
int | tc | ||
) |
Splits the receiver to LDLT form, and computes the rigid body motions.
r | matrix containing the rigid body motions base vectors. |
nse | number of rigid body motions |
se | array containing indexes of singular eqs. |
limit | determines linear dependence or independence |
tc | type of calculation = 1 = Decomposition to LDL form = 2 = ? konstruovani baze prostoru ? Ker A = 3 = Decomposition to LDL form and ? konstruovani baze prostoru ? Ker A |
Definition at line 785 of file skyline.C.
References adr, oofem::IntArray::at(), oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatMatrix::clear(), oofem::SparseMtrx::giveNumberOfRows(), mtrx, nwk, OOFEM_LOG_INFO, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
Referenced by oofem::FETISolver::solve().
int oofem::Skyline::setInternalStructure | ( | IntArray & | a | ) |
Allocates and builds internal structure according to given array holding addresses of diagonal members values (adr).
Definition at line 338 of file skyline.C.
References adr, oofem::IntArray::at(), oofem::IntArray::giveSize(), mtrx, oofem::SparseMtrx::nColumns, oofem::SparseMtrx::nRows, nwk, OOFEM_ERROR, and oofem::SparseMtrx::version.
|
virtual |
Evaluates .
x | Array to be multiplied with receiver. |
answer | y. |
Reimplemented from oofem::SparseMtrx.
Definition at line 578 of file skyline.C.
References adr, oofem::IntArray::at(), oofem::FloatArray::at(), oofem::SparseMtrx::giveNumberOfRows(), oofem::FloatArray::giveSize(), mtrx, OOFEM_ERROR, oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
virtual |
Multiplies receiver by scalar value.
x | Value to multiply receiver. |
Reimplemented from oofem::SparseMtrx.
Definition at line 621 of file skyline.C.
References mtrx, nwk, and oofem::SparseMtrx::version.
|
inlinevirtual |
Evaluates .
x | Array to be multiplied with transpose of the receiver. |
answer | y. |
Reimplemented from oofem::SparseMtrx.
|
virtual |
Converts receiving sparse matrix to a dense float matrix.
Reimplemented from oofem::SparseMtrx.
Definition at line 193 of file skyline.C.
References adr, oofem::IntArray::at(), oofem::FloatMatrix::at(), oofem::SparseMtrx::giveNumberOfColumns(), oofem::IntArray::giveSize(), mtrx, OOFEM_ERROR, oofem::FloatMatrix::resize(), oofem::FloatMatrix::symmetrized(), and oofem::FloatMatrix::zero().
Referenced by printYourself(), and writeToFile().
|
virtual |
Helpful for debugging, writes the matrix to given file.
Reimplemented from oofem::SparseMtrx.
Definition at line 654 of file skyline.C.
References oofem::FloatMatrix::at(), oofem::SparseMtrx::nColumns, oofem::SparseMtrx::nRows, and toFloatMatrix().
|
virtual |
Zeroes the receiver.
Implements oofem::SparseMtrx.
Definition at line 669 of file skyline.C.
References isFactorized, mtrx, nwk, and oofem::SparseMtrx::version.
|
protected |
Integer array holding addresses of diagonal members.
Definition at line 74 of file skyline.h.
Referenced by assemble(), at(), backSubstitutionWith(), buildInternalStructure(), factorized(), GiveCopy(), isAllocatedAt(), ldl_feti_sky(), rbmodes(), setInternalStructure(), Skyline(), times(), and toFloatMatrix().
|
protected |
Flag indicating whether factorized.
Definition at line 78 of file skyline.h.
Referenced by factorized(), Skyline(), and zero().
|
protected |
Vector of stored coefficients.
Definition at line 76 of file skyline.h.
Referenced by add(), assemble(), at(), backSubstitutionWith(), buildInternalStructure(), factorized(), GiveCopy(), ldl_feti_sky(), rbmodes(), setInternalStructure(), Skyline(), times(), toFloatMatrix(), zero(), and ~Skyline().
|
protected |
Total number of nonzero coefficients stored.
Definition at line 72 of file skyline.h.
Referenced by add(), buildInternalStructure(), factorized(), GiveCopy(), rbmodes(), setInternalStructure(), Skyline(), times(), and zero().