|
OOFEM 3.0
|
#include <skylineu.h>
Public Member Functions | |
| SkylineUnsym (int n=0) | |
| SkylineUnsym (const SkylineUnsym &s) | |
| virtual | ~SkylineUnsym () |
| Destructor. | |
| std::unique_ptr< SparseMtrx > | clone () const override |
| void | times (const FloatArray &x, FloatArray &answer) const override |
| void | timesT (const FloatArray &x, FloatArray &answer) const override |
| void | times (double x) override |
| int | buildInternalStructure (EngngModel *, int, const UnknownNumberingScheme &s) override |
| int | setInternalStructure (IntArray &a) |
| int | assemble (const IntArray &loc, const FloatMatrix &mat) override |
| int | assemble (const IntArray &rloc, const IntArray &cloc, const FloatMatrix &mat) override |
| bool | canBeFactorized () const override |
| Determines, whether receiver can be factorized. | |
| SparseMtrx * | factorized () override |
| FloatArray * | backSubstitutionWith (FloatArray &) const override |
| void | zero () override |
| Zeroes the receiver. | |
| double & | at (int i, int j) override |
| Returns coefficient at position (i,j). | |
| double | at (int i, int j) const override |
| Returns coefficient at position (i,j). | |
| void | toFloatMatrix (FloatMatrix &answer) const override |
| Converts receiving sparse matrix to a dense float matrix. | |
| void | printYourself () const override |
| Prints receiver to stdout. Works only for relatively small matrices. | |
| void | printStatistics () const override |
| Prints the receiver statistics (one-line) to stdout. | |
| void | writeToFile (const char *fname) const override |
| Helpful for debugging, writes the matrix to given file. | |
| SparseMtrxType | giveType () const override |
| Sparse matrix type identification. | |
| bool | isAsymmetric () const override |
| Returns true if asymmetric. | |
| const char * | giveClassName () const override |
| Public Member Functions inherited from oofem::SparseMtrx | |
| SparseMtrx (int n=0, int m=0) | |
| virtual | ~SparseMtrx () |
| Destructor. | |
| SparseMtrxVersionType | giveVersion () |
| Return receiver version. | |
| void | checkBounds (int i, int j) const |
| int | giveNumberOfRows () const |
| Returns number of rows of receiver. | |
| int | giveNumberOfColumns () const |
| Returns number of columns of receiver. | |
| bool | isSquare () const |
| Returns nonzero if receiver is square matrix. | |
| bool | isNotEmpty () const |
| Tests for empty matrix. | |
| virtual void | times (const FloatMatrix &B, FloatMatrix &answer) const |
| virtual void | timesT (const FloatMatrix &B, FloatMatrix &answer) const |
| virtual void | add (double x, SparseMtrx &m) |
| virtual void | addDiagonal (double x, FloatArray &m) |
| virtual int | buildInternalStructure (EngngModel *eModel, int n, int m, const IntArray &I, const IntArray &J) |
| virtual int | buildInternalStructure (EngngModel *eModel, int di, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s) |
| virtual int | assembleBegin () |
| Starts assembling the elements. | |
| virtual int | assembleEnd () |
| Returns when assemble is completed. | |
| virtual double | computeNorm () const |
| Returns the norm of receiver. | |
| virtual std::unique_ptr< SparseMtrx > | giveSubMatrix (const IntArray &rows, const IntArray &cols) |
| virtual bool | isAllocatedAt (int i, int j) const |
| Checks whether memory is allocated at position (i,j). | |
| std::string | errorInfo (const char *func) const |
| Error printing helper. | |
| FloatArray | operator* (const FloatArray &x) const |
| IML compatibility, \( A \cdot x\). | |
| FloatArray | trans_mult (const FloatArray &x) const |
| IML compatibility, \( A^{\mathrm{T}} \cdot x\). | |
Protected Member Functions | |
| void | checkSizeTowards (const IntArray &) |
| void | checkSizeTowards (const IntArray &rloc, const IntArray &cloc) |
| void | growTo (int) |
| SkylineUnsym (int n, std::vector< RowColumn > data, bool isFact) | |
Protected Attributes | |
| std::vector< RowColumn > | rowColumns |
| Row column segments. | |
| int | isFactorized |
| Factorization flag. | |
| Protected Attributes inherited from oofem::SparseMtrx | |
| int | nRows |
| Number of rows. | |
| int | nColumns |
| Number of columns. | |
| SparseMtrxVersionType | version |
Additional Inherited Members | |
| Public Types inherited from oofem::SparseMtrx | |
| typedef long | SparseMtrxVersionType |
This class implements a nonsymmetric matrix stored in a compacted (skyline) form. Its shape is symmetric, but not its coefficients. The Skyline contains 'size' row-column segments, each of them of any size; these are stored in the attribute 'rowColumns'.
Definition at line 54 of file skylineu.h.
| oofem::SkylineUnsym::SkylineUnsym | ( | int | n = 0 | ) |
Constructor. Before any operation an internal profile must be built.
| n | Size of matrix |
Definition at line 60 of file skylineu.C.
References isFactorized, and oofem::SparseMtrx::SparseMtrx().
Referenced by SkylineUnsym().
| oofem::SkylineUnsym::SkylineUnsym | ( | const SkylineUnsym & | s | ) |
Definition at line 66 of file skylineu.C.
References isFactorized, oofem::SparseMtrx::nColumns, oofem::SparseMtrx::nRows, rowColumns, SkylineUnsym(), and oofem::SparseMtrx::SparseMtrx().
|
inlinevirtual |
Destructor.
Definition at line 71 of file skylineu.h.
|
protected |
|
overridevirtual |
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 124 of file skylineu.C.
References oofem::FloatMatrix::at(), oofem::IntArray::at(), at(), checkSizeTowards(), oofem::FloatMatrix::giveNumberOfRows(), oofem::IntArray::giveSize(), OOFEM_ERROR, and oofem::SparseMtrx::version.
|
overridevirtual |
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 158 of file skylineu.C.
References oofem::FloatMatrix::at(), oofem::IntArray::at(), at(), checkSizeTowards(), oofem::IntArray::giveSize(), and oofem::SparseMtrx::version.
|
overridevirtual |
Returns coefficient at position (i,j).
Implements oofem::SparseMtrx.
Definition at line 111 of file skylineu.C.
References rowColumns.
|
overridevirtual |
Returns coefficient at position (i,j).
Implements oofem::SparseMtrx.
Definition at line 96 of file skylineu.C.
References rowColumns, and oofem::SparseMtrx::version.
Referenced by assemble(), assemble(), and toFloatMatrix().
|
overridevirtual |
Computes the solution of linear system \( A\cdot x = y \) 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 462 of file skylineu.C.
References oofem::FloatArray::at(), oofem::diag(), oofem::SparseMtrx::giveNumberOfColumns(), oofem::FloatArray::giveSize(), OOFEM_ERROR, rowColumns, and SkylineUnsym_TINY_PIVOT.
|
overridevirtual |
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 214 of file skylineu.C.
References oofem::IntArray::at(), oofem::Domain::giveBc(), oofem::EngngModel::giveDomain(), oofem::Domain::giveElement(), oofem::ActiveBoundaryCondition::giveLocationArrays(), oofem::Domain::giveNumberOfBoundaryConditions(), oofem::EngngModel::giveNumberOfDomainEquations(), oofem::Domain::giveNumberOfElements(), oofem::IntArray::giveSize(), growTo(), oofem::min(), oofem::SparseMtrx::nColumns, oofem::SparseMtrx::nRows, printStatistics(), rowColumns, oofem::EngngModel::useNonlocalStiffnessOption(), and oofem::SparseMtrx::version.
|
inlineoverridevirtual |
Determines, whether receiver can be factorized.
Implements oofem::SparseMtrx.
Definition at line 82 of file skylineu.h.
|
protected |
Definition at line 368 of file skylineu.C.
References growTo(), oofem::max(), and rowColumns.
Referenced by assemble(), and assemble().
|
protected |
Definition at line 183 of file skylineu.C.
References oofem::SparseMtrx::giveNumberOfColumns(), growTo(), oofem::max(), and rowColumns.
|
overridevirtual |
Returns a newly allocated copy of receiver. Programmer must take care about proper deallocation of allocated space.
Reimplemented from oofem::SparseMtrx.
Definition at line 74 of file skylineu.C.
|
overridevirtual |
Returns the receiver factorized. \( L^{\mathrm{T}} \cdot D \cdot L \) form is used.
Reimplemented from oofem::SparseMtrx.
Definition at line 388 of file skylineu.C.
References oofem::FloatArray::at(), oofem::diag(), oofem::Timer::getUtime(), oofem::SparseMtrx::giveNumberOfColumns(), isFactorized, oofem::max(), OOFEM_LOG_DEBUG, OOFEM_WARNING, oofem::FloatArray::resize(), rowColumns, SkylineUnsym_TINY_PIVOT, oofem::Timer::startTimer(), oofem::Timer::stopTimer(), and oofem::FloatArray::zero().
|
inlineoverridevirtual |
Implements oofem::SparseMtrx.
Definition at line 94 of file skylineu.h.
|
inlineoverridevirtual |
Sparse matrix type identification.
Implements oofem::SparseMtrx.
Definition at line 92 of file skylineu.h.
References oofem::SMT_SkylineU.
|
protected |
Definition at line 541 of file skylineu.C.
References oofem::SparseMtrx::nColumns, oofem::SparseMtrx::nRows, and rowColumns.
Referenced by buildInternalStructure(), checkSizeTowards(), checkSizeTowards(), and setInternalStructure().
|
inlineoverridevirtual |
|
overridevirtual |
Prints the receiver statistics (one-line) to stdout.
Reimplemented from oofem::SparseMtrx.
Definition at line 552 of file skylineu.C.
References oofem::SparseMtrx::giveNumberOfRows(), OOFEM_LOG_INFO, and rowColumns.
Referenced by buildInternalStructure(), and setInternalStructure().
|
overridevirtual |
Prints receiver to stdout. Works only for relatively small matrices.
Reimplemented from oofem::SparseMtrx.
Definition at line 580 of file skylineu.C.
References oofem::FloatMatrix::printYourself(), rowColumns, and toFloatMatrix().
| int oofem::SkylineUnsym::setInternalStructure | ( | IntArray & | a | ) |
Definition at line 332 of file skylineu.C.
References oofem::IntArray::at(), oofem::IntArray::giveSize(), growTo(), oofem::SparseMtrx::nColumns, oofem::SparseMtrx::nRows, printStatistics(), oofem::IntArray::resize(), rowColumns, oofem::SparseMtrx::version, and oofem::IntArray::zero().
|
overridevirtual |
Evaluates \( y = A \cdot x \)
| x | Array to be multiplied with receiver. |
| answer | y. |
Reimplemented from oofem::SparseMtrx.
Definition at line 501 of file skylineu.C.
References oofem::FloatArray::at(), oofem::SparseMtrx::giveNumberOfColumns(), oofem::SparseMtrx::giveNumberOfRows(), oofem::FloatArray::giveSize(), OOFEM_ERROR, oofem::FloatArray::resize(), rowColumns, and oofem::FloatArray::zero().
|
overridevirtual |
Multiplies receiver by scalar value.
| x | Value to multiply receiver. |
Reimplemented from oofem::SparseMtrx.
Definition at line 523 of file skylineu.C.
References oofem::SparseMtrx::giveNumberOfColumns(), rowColumns, and oofem::SparseMtrx::version.
|
overridevirtual |
Evaluates \( y = A^{\mathrm{T}} \cdot x \)
| x | Array to be multiplied with transpose of the receiver. |
| answer | y. |
Reimplemented from oofem::SparseMtrx.
Definition at line 603 of file skylineu.C.
References oofem::FloatArray::at(), oofem::SparseMtrx::giveNumberOfColumns(), oofem::SparseMtrx::giveNumberOfRows(), oofem::FloatArray::giveSize(), OOFEM_ERROR, oofem::FloatArray::resize(), rowColumns, and oofem::FloatArray::zero().
|
overridevirtual |
Converts receiving sparse matrix to a dense float matrix.
Reimplemented from oofem::SparseMtrx.
Definition at line 81 of file skylineu.C.
References oofem::FloatMatrix::at(), at(), oofem::SparseMtrx::giveNumberOfColumns(), oofem::SparseMtrx::giveNumberOfRows(), oofem::FloatMatrix::resize(), rowColumns, and oofem::FloatMatrix::zero().
Referenced by printYourself(), and writeToFile().
|
overridevirtual |
Helpful for debugging, writes the matrix to given file.
Reimplemented from oofem::SparseMtrx.
Definition at line 564 of file skylineu.C.
References oofem::FloatMatrix::at(), oofem::SparseMtrx::nColumns, oofem::SparseMtrx::nRows, and toFloatMatrix().
|
overridevirtual |
Zeroes the receiver.
Implements oofem::SparseMtrx.
Definition at line 591 of file skylineu.C.
References isFactorized, rowColumns, and oofem::SparseMtrx::version.
|
protected |
Factorization flag.
Definition at line 60 of file skylineu.h.
Referenced by factorized(), SkylineUnsym(), SkylineUnsym(), and zero().
|
protected |
Row column segments.
Definition at line 58 of file skylineu.h.
Referenced by at(), at(), backSubstitutionWith(), buildInternalStructure(), checkSizeTowards(), checkSizeTowards(), factorized(), growTo(), printStatistics(), printYourself(), setInternalStructure(), SkylineUnsym(), times(), times(), timesT(), toFloatMatrix(), and zero().