Class xlifepp::TermMatrix#
-
class TermMatrix : public xlifepp::Term#
-
Inheritence diagram for xlifepp::TermMatrix:
Collaboration diagram for xlifepp::TermMatrix:
end user class handling numerical representation of any bilinear form (BilinearForm).
May be multiple unknowns
Public Functions
-
inline TermMatrix(const LcTerm<TermMatrix> &lctm)#
-
constructor from a linear combination of TermMatrix’s
-
inline TermMatrix(const LcTerm<TermMatrix> &lctm, const char *na)#
-
constructor from a linear combination of TermMatrix’s
-
TermMatrix(const LcTerm<TermMatrix> &lctm, const Parameter &p1)#
-
constructor from a linear combination of TermMatrix’s
-
inline TermMatrix(const LcTerm<TermMatrix> &lctm, const string_t &na)#
-
constructor from a linear combination of TermMatrix’s
-
template<typename T>
inline explicit TermMatrix(const Matrix<T> &mat)#
-
construct TermMatrix from a standard Matrix (virtual unknown)
-
template<typename T>
inline TermMatrix(const Matrix<T> &mat, const char *na)#
-
construct TermMatrix from a standard Matrix (virtual unknown)
-
template<typename T>
TermMatrix(const Matrix<T> &mat, const Parameter &p1)#
-
construct TermMatrix from a standard Matrix (virtual unknown)
-
template<typename T>
inline TermMatrix(const Matrix<T> &mat, const string_t &na)#
-
construct TermMatrix from a standard Matrix (virtual unknown)
-
template<typename T>
TermMatrix(const std::vector<T> &mat, number_t m, number_t n = 0, const string_t &na = "")#
-
construct TermMatrix from a row dense matrix of size m x n (virtual unknown)
-
TermMatrix(const string_t &na = "?")#
-
default constructor
-
TermMatrix(const SuTermMatrix&, const string_t &na = "")#
-
constructor of a TermMatrix from a SuTermMatrix (full copy)
constructor from a SuTermMatrix, create a one block TermMatrix
-
TermMatrix(const TermMatrix &tm, const string_t &na = "")#
-
copy constructor (full copy)
-
inline TermMatrix(const TermMatrix &tm, SpecialMatrix sm)#
-
construct special matrix from a TermMatrix (e.g _IdMatrix)
-
inline TermMatrix(const TermMatrix &tm, SpecialMatrix sm, const char *na)#
-
construct special matrix from a TermMatrix (e.g _IdMatrix)
-
inline TermMatrix(const TermMatrix &tm, SpecialMatrix sm, const complex_t &a)#
-
construct special matrix from a TermMatrix (e.g _IdMatrix)
-
inline TermMatrix(const TermMatrix &tm, SpecialMatrix sm, const complex_t &a, const char *na)#
-
construct special matrix from a TermMatrix (e.g _IdMatrix)
-
TermMatrix(const TermMatrix &tm, SpecialMatrix sm, const complex_t &a, const Parameter &p1)#
-
construct special matrix from a TermMatrix (e.g _IdMatrix)
-
inline TermMatrix(const TermMatrix &tm, SpecialMatrix sm, const complex_t &a, const string_t &na)#
-
construct special matrix from a TermMatrix (e.g _IdMatrix)
-
inline TermMatrix(const TermMatrix &tm, SpecialMatrix sm, const Parameter &p1)#
-
construct special matrix from a TermMatrix (e.g _IdMatrix)
-
inline TermMatrix(const TermMatrix &tm, SpecialMatrix sm, const real_t &a)#
-
construct special matrix from a TermMatrix (e.g _IdMatrix)
-
inline TermMatrix(const TermMatrix &tm, SpecialMatrix sm, const real_t &a, const char *na)#
-
construct special matrix from a TermMatrix (e.g _IdMatrix)
-
inline TermMatrix(const TermMatrix &tm, SpecialMatrix sm, const real_t &a, const Parameter &p1)#
-
construct special matrix from a TermMatrix (e.g _IdMatrix)
-
inline TermMatrix(const TermMatrix &tm, SpecialMatrix sm, const real_t &a, const string_t &na)#
-
construct special matrix from a TermMatrix (e.g _IdMatrix)
-
inline TermMatrix(const TermMatrix &tm, SpecialMatrix sm, const string_t &na)#
-
construct special matrix from a TermMatrix (e.g _IdMatrix)
-
explicit TermMatrix(const TermVector &tv, const string_t &na = "")#
-
constructor of a diagonal TermMatrix from a TermVector
constructor of a diagonal TermMatrix from a TermVector if TermVector is a multiple unknowns vector, the matrix will be a multiple unknowns TermMatrix, each block being diagonal !
-
inline TermMatrix(const TermVector &tv1, const TermVector &tv2, const SymbolicFunction &sf)#
-
sf(tv1, tv2)
-
inline TermMatrix(const TermVector &tv1, const TermVector &tv2, const SymbolicFunction &sf, const char *na)#
-
sf(tv1, tv2)
-
TermMatrix(const TermVector &tv1, const TermVector &tv2, const SymbolicFunction &sf, const Parameter &p1)#
-
sf(tv1, tv2)
-
inline TermMatrix(const TermVector &tv1, const TermVector &tv2, const SymbolicFunction &sf, const string_t &na)#
-
sf(tv1, tv2)
-
template<typename T>
inline TermMatrix(const Unknown &u, const GeomDomain &dom, const Vector<T> &vt)#
-
special constructor of a diagonal TermMatrix from a vector
-
template<typename T>
inline TermMatrix(const Unknown &u, const GeomDomain &dom, const Vector<T> &vt, const char *na)#
-
special constructor of a diagonal TermMatrix from a vector
-
template<typename T>
TermMatrix(const Unknown &u, const GeomDomain &dom, const Vector<T> &vt, const Parameter &p1)#
-
special constructor of a diagonal TermMatrix from a vector
-
template<typename T>
inline TermMatrix(const Unknown &u, const GeomDomain &dom, const Vector<T> &vt, const string_t &na)#
-
special constructor of a diagonal TermMatrix from a vector
-
template<typename T>
TermMatrix(const Unknown &u, const TestFunction &v, const std::vector<T> &mat, number_t m, number_t n, const string_t &na = "")#
-
construct TermMatrix from a row dense matrix of size m x n (virtual unknown)
-
inline TermMatrix(const Unknown &v, const GeomDomain &domv, const Unknown &u, const GeomDomain &domu, const OperatorOnKernel &opker)#
-
constructor of matrix opker(xi,yj) (Lagrange interpolation)
-
inline TermMatrix(const Unknown &v, const GeomDomain &domv, const Unknown &u, const GeomDomain &domu, const OperatorOnKernel &opker, const char *na)#
-
constructor of matrix opker(xi,yj) (Lagrange interpolation)
-
TermMatrix(const Unknown &v, const GeomDomain &domv, const Unknown &u, const GeomDomain &domu, const OperatorOnKernel &opker, const Parameter &p1)#
-
constructor of matrix opker(xi,yj) (Lagrange interpolation)
constructor of matrix opker(xi,yj) - Lagrange interpolation - opker is evaluated on dof coordinates linked to Unknown and GeomDomains v, domv is associated to matrix rows and u, domu to matrix columns NOTE: do not use this constructor if kernel is singular on some couple of points
-
inline TermMatrix(const Unknown &v, const GeomDomain &domv, const Unknown &u, const GeomDomain &domu, const OperatorOnKernel &opker, const string_t &na)#
-
constructor of matrix opker(xi,yj) (Lagrange interpolation)
-
TermMatrix(SuTermMatrix*, const string_t &na = "")#
-
constructor of a TermMatrix from a SuTermMatrix * (pointer copy)
-
virtual ~TermMatrix()#
-
destructor
-
MatrixEntry *actual_entries() const#
-
return the actual pointer to entries (priority to scalar entry)
-
void addIndices(std::vector<std::set<number_t>> &indices, MatrixStorage *msu, const std::vector<int_t> &rowmap, const std::vector<int_t> &colmap)#
-
merging tool used by toGlobal()
-
void applyEssentialConditions()#
-
apply essential conditions to a computed TermMatrix (internal tool)
-
void applyEssentialConditions(const EssentialConditions &ecs, const ReductionMethod &rm = ReductionMethod(_pseudoReduction))#
-
apply essential conditions to a computed TermMatrix (user tool)
-
void applyEssentialConditions(const EssentialConditions &ecsu, const EssentialConditions &ecsv, const ReductionMethod &rm = ReductionMethod(_pseudoReduction))#
-
apply essential conditions to a computed TermMatrix (user tool)
-
inline it_mustm begin()#
-
iterator to the first element of the SuTermMatrix map
-
inline cit_mustm begin() const#
-
iterator to the first element of the SuTermMatrix map (const)
-
void build(const BilinearForm &bf, const EssentialConditions *ecsu, const EssentialConditions *ecsv, const std::vector<Parameter> &ps)#
-
effective constructor from bilinear form, essential conditions and options
effective constructor from BilinearForm and essential conditions
-
void build(const TermMatrix &tm, const EssentialConditions *ecsu, const EssentialConditions *ecsv, const std::vector<Parameter> &ps)#
-
effective constructor from TermMatrix, essential conditions and options
effective constructor from TermMatrix and essential conditions
-
template<typename T>
void buildFromRowDenseMatrix(const std::vector<T> &mat, number_t m, number_t n, const string_t &na)#
-
real constructor from matrix
-
void buildParam(const Parameter &p, bool &toCompute, bool &toAssemble, SymType &sym)#
-
returns the list of authorized keys
-
void buildSpecialMatrix(const TermMatrix &A, SpecialMatrix sm, const complex_t &a, const string_t &na)#
-
effective constructor of special matrix from a TermMatrix (e.g _IdMatrix)
construct special matrix from a structure given by a TermMatrix A: the matrix model sm: special matrix type (up to now, only _idMatrix is available) a: coefficient to apply na: option name
-
inline std::vector<DofComponent> &cdofsc()#
-
access to cdofs_c vector (r/w)
-
inline const std::vector<DofComponent> &cdofsc() const#
-
access to cdofs_c vector (r)
-
inline std::vector<DofComponent> &cdofsr()#
-
access to cdofs_r vector (r/w)
-
inline const std::vector<DofComponent> &cdofsr() const#
-
access to cdofs_r vector (r)
-
void changeTestFunction(const Unknown &newv)#
-
change first testFuncion (row) to an other (be cautious)
-
void changeTestFunction(const Unknown &oldv, const Unknown &newv)#
-
change oldv testfunction (row) to an other (be cautious)
-
void changeTestFunctions(const Unknowns &oldvs, const Unknowns &newvs)#
-
change set of current test functions to others (be cautious)
change collection of current test functions to others (be cautious)
-
void changeUnknown(const Unknown &oldu, const Unknown &newu)#
-
change oldu unknown (col) to an other (be cautious)
-
void changeUnknowns(const Unknown &newu, const Unknown &newv)#
-
change both first unknown and TestFuncion to others (be cautious)
-
void changeUnknowns(const Unknowns &oldus, const Unknowns &newvs)#
-
change set of current unknowns to others (be cautious)
change collection of current unknowns to others (be cautious)
-
void changeUnknowns(const Unknowns &oldus, const Unknowns &oldvs, const Unknowns &newus, const Unknowns &newvs)#
-
change set of current unknowns and test functions to others (be cautious)
change collection of current unknowns and test functions to others (be cautious)
-
virtual void clear()#
-
deallocate memory used by a TermMatrix
-
number_t colRank(const Dof &dof) const#
-
rank (>=1) in first submatrix column numbering of a Dof rank (>=1) in (u,v) submatrix row numbering of a Dof
-
number_t colRank(const Unknown &u, const Unknown &v, const Dof &dof) const#
-
rank (>=1) in (u,v) submatrix column numbering of a Dof
rank (>=1) in submatrix (u,v) row numbering of a Dof
-
TermVector column(number_t c) const#
-
return the j-th col as a TermVector (j=1,…)
-
virtual void compute()#
-
compute TermMatrix from bilinear form compute TermMatrix from a linear combination of TermMatrix’s
-
void compute(const LcTerm<TermMatrix> &lctm, const string_t &na = "")#
-
compute a LcTerm (linear combination of TermMatrix’s) to a TermMatrix all the TermMatrix’s have to be computed before and the current TermMatrix is ‘void’
NOTE: this current version works only for “standard” representation (entries pointers) but NOT for scalar representation (scalar_entries pointers) it means that TermMatrix that have been converted using toScalar or toGlobal functions cannot be used in linear combinations. Such conversions occur when essential conditions have been applied to a TermMatrix To be improved in future !
-
void copy(const TermMatrix &tm)#
-
copy a TermMatrix in current TermMatrix
-
inline it_mustm end()#
-
iterator to the last element of the SuTermMatrix map
-
inline cit_mustm end() const#
-
iterator to the last element of the SuTermMatrix map (const)
-
inline MatrixEntry *&entries()#
-
access to entries pointer (r/w)
-
inline const MatrixEntry *entries() const#
-
access to entries pointer (r)
-
FactorizationType factorization() const#
-
return factorization type if factorized
-
std::pair<StorageType, AccessType> findGlobalStorageType() const#
-
suggest storage type for global representation
regarding a multiple unknowns TermMatrix, suggest a storage according to the following rules: if dense part > 10*(cs part+skyline part) -> dense storage else -> cs storage if umfpack installed -> column access else dual row/col access if matrix has symmetry -> symmetric access
-
inline SuTermMatrix *firstSut()#
-
first SuTermMatrix as pointer, 0 if void TermMatrix
-
inline SuTermMatrix *firstSut() const#
-
first SuTermMatrix as pointer, 0 if void TermMatrix (const)
-
template<typename T>
HMatrix<T, FeDof> &getHMatrix(StrucType st = _undefStrucType) const#
-
get internal HMatrix if exists, advanced usage
get internal H-matrix (restrictive), i.e HMatrix object if available available only for single unknown TermMatrix if scalar and vector representations both exist, the scalar representation is prior except if StrucType argument is specified It returns a reference to a HMatrix Note: designed for advanced usage
-
template<typename T>
LargeMatrix<T> &getLargeMatrix(StrucType st = _undefStrucType) const#
-
get internal matrix (restrictive), advanced usage
get internal matrix (restrictive), i.e LargeMatrix object if available available only for single unknowns TermMatrix and multiple unknowns TermMatrix having a global representation if scalar and vector representations both exist, the scalar representation is prior except if StrucType argument is specified user has to specify the LargeMatrix type: TermMatrix M(intg(omega,u*v)); LargeMatrix& Alm=M.getMatrix<Real>(); It returns a reference to the LargeMatrix Note: designed for advanced usage
-
Value getScalarValue(const Unknown&, const Unknown&, number_t, number_t, dimen_t = 0, dimen_t = 0) const#
-
access to scalar (u,v) sub-matrix coefficient (i,j >= 1)
-
Value getScalarValue(number_t, number_t, dimen_t = 0, dimen_t = 0) const#
-
access to scalar first sub-matrix coefficient (i,j >= 1)
-
inline Value getValue(const Dof &du, const Dof &dv) const#
-
access to first sub-matrix coefficient related to dofs (dv,du)
-
Value getValue(const Unknown&, const Unknown&, number_t, number_t) const#
-
access to (u,v) sub-matrix coefficient (i,j >= 1)
-
inline Value getValue(const Unknown &u, const Unknown &v, const Dof &du, const Dof &dv) const#
-
access to (u,v) sub-matrix coefficient related to dofs (dv,du)
-
bool hasBlock(const Unknown &u, const Unknown &v) const#
-
true if TermMatrix has a block (u,v)
-
inline bool hasConstraints() const#
-
true if one constraints pointer is not null
-
void initFromBF(const BilinearForm &bf, const string_t &na = "", bool noass = false)#
-
initializer with bf
-
template<typename T>
void initTermVectorT(TermVector &tv, const T &t, bool col = true) const#
-
initialize a TermVector from value T, consistent with matrix rows or columns
initialize a TermVector consistent with matrix column unknown (col=true), with matrix row unknown (col=false)
-
void insert(const SuTermMatrix &sutm)#
-
insert a SuTermMatrix with full copy
-
void insert(SuTermMatrix *sutm)#
-
insert a SuTermMatrix with pointer copy initialize a TermVector consistent with matrix rows or columns
-
void insertDiagonalBlocks()#
-
insert diagonal blocks if missing (multiple unknowns)
insert 0 diagonal blocks when missing only for multiple unknowns TermMatrix and rowUnknowns()=colUnknowns()
-
bool isScalar() const#
-
true if all SuTermMatrix’s are scalar
-
inline bool isSingleUnknown() const#
-
true if a single unknown matrix
-
void ldlstarSolve(const TermVector&, TermVector&)#
-
solve linear system using LDL* factorization
-
void ldltSolve(const TermVector&, TermVector&)#
-
solve linear system using LDLt factorization
-
void luSolve(const TermVector&, TermVector&)#
-
solve linear system using LU factorization
-
void markAsComputed(bool isComputed = true)#
-
TermMatrix and SutermMatrix’s computed state = true or false.
-
MatrixEntry *matrixData()#
-
turns to scalar representation and returns pointer to data
-
void mergeBlocks()#
-
merge blocks referring to same vector unknowns
merge blocks referring to same vector unknown
when a block has a component of unknown as row or col unknown, it is convert to its vector representation
when some blocks have different components of same unknown as row or col unknown, the blocks are really merged into a block with vector unknowns see the function mergeSuTermMatrix
-
void mergeNumberingFast(std::map<const Unknown*, std::list<SuTermMatrix*>> &rcsut, std::map<SuTermMatrix*, std::vector<int_t>> &rcnum, std::vector<DofComponent> &rcdof, AccessType rc)#
-
merging tool used by toGlobal(), fast version
-
virtual void name(const string_t &nam)#
-
update the name of a TermMatrix and its SuTermMatrixs
update the name of a TermMatrix and its SuTermMatrixs SuTermMatrix name is defined as follows: name(TermMatrix)_name(row_Unknown)_name(col_Unknown)
-
real_t norm2() const#
-
Frobenius norm: sqrt(sum|a_ij|^2)
-
real_t norminfty() const#
-
infinite norm: sup|a_ij|
-
TermMatrix operator()(const Unknown &u, const Unknown &v) const#
-
access to single unknowns block matrix as TermMatrix
return single unknown term matrix as TermMatrix (copy constructor like), useful to extract single unknown term from multiple unknowns term induce a full copy !!!
-
template<typename T>
TermMatrix &operator*=(const T &t)#
-
TermMatrix *= t.
product assignment by a real or a complex (template T)
-
TermMatrix &operator+=(const TermMatrix &tM)#
-
TermMatrix += TermMatrix.
Algebraic operation += and -= the result of these operations is a TermMatrix with all unknown blocks of current TermMatrix and added Termatrix | Avu Avp | | Mvu | | Avu+Mvu Avp | | Aqu | += | Mqu Mqp | ==> | Aqu+Mqu Mqp | | Atp Ats| | Mtu Mts| | Mtu Atp Ats+Mts| in other words, the matrices are merged and summation is a real one on common unknowns blocks.
-
TermMatrix &operator-=(const TermMatrix &tM)#
-
TermMatrix -= TermMatrix.
-
template<typename T>
TermMatrix &operator/=(const T &t)#
-
TermMatrix /= t.
division assignment by a real or a complex (template T)
-
TermMatrix &operator=(const LcTerm<TermMatrix> &lctm)#
-
assign operator from a linear combination of TermMatrix’s
-
TermMatrix &operator=(const TermMatrix &tm)#
-
assign operator (full copy)
-
void penalizationReduction()#
-
reduction of essential condition using penalization method
add penalization coefficient alpha to theTermMatrix only for Dirichlet conditions like (no scalar dofs coupling)
-
virtual void print(std::ostream &out) const#
-
print TermMatrix
-
void printSummary(std::ostream &out) const#
-
print summary TermMatrix
-
void pseudoReduction()#
-
reduction of essential condition using pseudo reduction method
perform pseudo reduction of reduced essential conditions in a matrix.
Essential conditions have the following form: U_E + F*U_R = s for column unknown U V_E + G*V_R = 0 for row test function V (generally related to unknown U) where E are the eliminated unknown/test function indices and R are the reduced unknown/test function indices The pseudo reduction of matrix consists in
modifying column A.j for j in R by adding -Fjk*A.k for k in E and replacing column A.k for k in E by a 0 column
modifying row Ai. for i in R by adding -Gik*Ak. for k in E and replacing row Ak. for k in E by a 0 row
if eliminated v unknowns are dual of eliminated u unknowns, the VE rows are replaced by the equation (or a part of) U_E + F*U_R = s to delay the right hand side modification, the (A.k) columns (k in E) are stored in a new structure
At the end of the process, the eliminated system looks like U_E U_R U_S RHS
| | | | | | V_E | Id | F | 0 | | S | | | | | | |
| | | | | | V_R | 0 | ARR | ARS | | FR | | | | | | |
| | | | | | V_S | 0 | ASR | ASS | | FS | | | | | | |
In some cases (F=G=0 or local conditions u^n=g, …) the storage is not modified but in other cases (transmission condition for instance) the storage is modified
-
void realReduction()#
-
reduction of essential condition using real reduction
essential condition reduction of TermMatrix first performms pseudo-reduction and then reduce the TermMatrix by omitting rows and columns corresponding to eliminated dofs The dofs numbering vector cdofsr and cdofsc are not modified ! Expansive, in future do the real reduction
-
inline MatrixEntry *&rhs_matrix()#
-
access to rhs matrix pointer (r/w)
-
inline const MatrixEntry *rhs_matrix() const#
-
access to rhs matrix pointer (r)
-
TermMatrix &roundToZero(real_t aszero = 10 * theEpsilon)#
-
round to zero all coefficients close to 0 (|.
| < aszero)
-
TermVector row(number_t r) const#
-
return the i-th row as a TermVector (i=1,…)
-
number_t rowRank(const Dof &dof) const#
-
rank (>=1) in first submatrix row numbering of a Dof
rank (>=1) in submatrix column numbering of a Dof
-
number_t rowRank(const Unknown &u, const Unknown &v, const Dof &dof) const#
-
rank (>=1) in submatrix (u,v) row numbering of a Dof
-
inline void saveToFile(const string_t &filename, StorageType st, bool encode, real_t tol = theTolerance) const#
-
save TermMatrix to file (dense or Matlab format)
-
void saveToFile(const string_t &filename, StorageType st, number_t prec = fullPrec, bool encode = false, real_t tol = theTolerance) const#
-
save TermMatrix to file (dense or Matlab format)
-
inline MatrixEntry *&scalar_entries()#
-
access to scalar entries pointer (r/w)
-
inline const MatrixEntry *scalar_entries() const#
-
access to scalar entries pointer (r)
-
void setCol(const Unknown&, const Unknown&, const Value&, number_t c1, number_t c2)#
-
set values of sub_matrix cols c1->c2 (>= 1), no storage change
-
inline void setCol(const Unknown &u, const Unknown &v, const Value &val, const Dof &d)#
-
set values of (u,v) sub-matrix col related to dof
-
void setCol(const Value&, number_t c1, number_t c2)#
-
set values of first sub-matrix cols c1->c2 (>= 1), no storage change
-
inline void setCol(const Value &val, const Dof &d)#
-
set values of first sub-matrix col related to dof
-
void setRow(const Unknown&, const Unknown&, const Value&, number_t r1, number_t r2)#
-
set values of sub_matrix rows r1->r2 (>= 1), no storage change
-
inline void setRow(const Unknown &u, const Unknown &v, const Value &val, const Dof &d)#
-
set values of (u,v) sub-matrix row related to dof
-
void setRow(const Value&, number_t r1, number_t r2)#
-
set values of first sub-matrix rows r1->r2 (>= 1), no storage change
-
inline void setRow(const Value &val, const Dof &d)#
-
set values of first sub-matrix row related to dof
-
void setScalarValue(const Unknown&, const Unknown&, number_t, number_t, const Value&, dimen_t = 0, dimen_t = 0)#
-
set value of scalar (u,v) sub_matrix coefficient (i,j >= 1) in storage
-
void setScalarValue(number_t, number_t, const Value&, dimen_t = 0, dimen_t = 0)#
-
set value of scalar first sub_matrix coefficient (i,j >= 1) in storage
set values of first sub-matrix rows/cols rc1->rc2 (>= 1), no storage change
-
void setStorage(StorageType st, AccessType at)#
-
specify the storage of SuTermMatrix’s
-
inline void setValue(const Dof &du, const Dof &dv, const Value &val)#
-
set value of first sub-matrix coefficient related to dofs (dv,du)
-
void setValue(const Unknown&, const Unknown&, number_t, number_t, const Value&)#
-
set value of (u,v) ub-matrix coefficient (i,j >= 1) in storage
-
inline void setValue(const Unknown &u, const Unknown &v, const Dof &du, const Dof &dv, const Value &val)#
-
set value of (u,v) sub-matrix coefficient related to dofs (dv,du)
-
void setValue(number_t, number_t, const Value&)#
-
set value of first sub-matrix coefficient (i,j >= 1) in storage
-
void sorDiagonalMatrixVector(const TermVector &V, TermVector &R, const real_t w)#
-
solve linear system using SOR for block diagonal matrices
-
void sorDiagonalSolve(const TermVector &V, TermVector &R, const real_t w)#
-
solve linear system using SOR for diagonal matrices
-
void sorLowerSolve(const TermVector &V, TermVector &R, const real_t w)#
-
solve linear system using SOR for lower triangular matrices
-
void sorSolve(const TermVector &V, TermVector &R, const real_t w, SorSolverType sType)#
-
solve linear system using SOR
-
void sorUpperSolve(const TermVector &V, TermVector &R, const real_t w)#
-
solve linear system using SOR for upper triangular matrices
-
inline AccessType storageAccess()#
-
returns storage access (row, col, sym, dual)
-
inline StorageType storageType()#
-
returns storage type (cs, skyline, dense)
-
SuTermMatrix &subMatrix(const Unknown *up, const Unknown *vp)#
-
access to SuTermMatrix
-
const SuTermMatrix &subMatrix(const Unknown *up, const Unknown *vp) const#
-
access to SuTermMatrix (const)
-
SuTermMatrix *subMatrix_p(const Unknown *up, const Unknown *vp)#
-
access to SuTermMatrix pointer
-
const SuTermMatrix *subMatrix_p(const Unknown *up, const Unknown *vp) const#
-
access to SuTermMatrix pointer (const) return the i-th row or col as a TermVector (i=1,…)
-
inline std::map<uvPair, SuTermMatrix*> &suTerms()#
-
access to suterms
-
SymType symmetry() const#
-
return global symmetry property
return global symmetry property (_noSymmetry, _symmetric, _adjoint, _skewAdjoint, _skewSymmetric) that is same symmetry property of diagonal unknowns blocks it there are no non diagonal blocks if there are non diagonal blocks, their symmetry is not checked (to do later) NOTE: this function works only for allocated SCALAR entries
-
TermMatrix &toComplex()#
-
go to complex representation
-
TermMatrix &toConj()#
-
change to conj(U)
-
void toGlobal(StorageType st, AccessType at, SymType symt = _noSymmetry, bool keepSuTerms = false)#
-
create scalar global representation of TermMatrix make scalar representations of SuTermMatrix’s if not exist and merge it in a scalar global representation by allocating scalar_entries_p and setting the cdofs vectors
st, at: storage and access type of global matrix, if st = _nosStorage try to identify “best” storage symt: to force symmetry of matrix (default = _noSymmetry) keepSuterms: if false (default value), SutermMatrix entries are deleted
-
TermMatrix &toImag()#
-
go to real representation getting the imag part when complex
-
TermMatrix &toReal()#
-
go to real representation getting the real part when complex
-
void toScalar(bool keepEntries = false)#
-
create scalar representation of TermMatrix create global scalar representation of TermMatrix
create scalar representation of TermMatrix for each SuTermMatrix (if computed, entries_!=nullptr), if not exists create its scalar representation (set scalar_entries_p and cdofs vectors) if SuTermMatrix is already scalar, scalar_entries_p = entries_p NOTE: scalar_entries_p of TermMatrix is not update by this function (see toGlobal function)
-
void toSkyline()#
-
convert to skyline storage (works only for cs storage single unknown matrix) merging tool used by toGlobal()
-
void umfpackSolve(const TermVector&, TermVector&)#
-
solve linear system using umfpack
-
void viewStorage(std::ostream &out) const#
-
print storage on stream
Friends
-
friend TermMatrix imag(const TermMatrix &tm)#
-
return imag part as a real TermMatrix
-
friend TermVector &multMatrixVector(const TermMatrix&, const TermVector&, TermVector&)#
-
product TermMatrix * TermVector
-
friend TermVector &multVectorMatrix(const TermMatrix&, const TermVector&, TermVector&)#
-
product TermVector * TermMatrix
-
friend TermVector &multVectorMatrix(const TermVector&, const TermMatrix&, TermVector&)#
-
product TermVector * TermMatrix
-
friend TermMatrix operator*(const TermMatrix&, const TermMatrix&)#
-
product of TermMatrix
-
friend TermVector operator*(const TermMatrix&, const TermVector&)#
-
product TermMatrix * TermVector
product of a TermMatrix A and a TermVector X There are 2 cases:
TermMatrix has a global representation, say scalar_entries_p !=nullptr or entries_p!=nullptr in that case the product is realized as a standard matrix/vector product Note that TermVector X has to have a global representation consistent with column numbering of TermMatrix if it not the case, global representation of X is computed
TermMatrix has a local representation, say scalar_entries_p =0 and entries_p=0 Assume that A has (v1,v2, …, vm) has row unknowns and (u1,u2, …, un) has column unknowns and X has (p1,p2, …, pq) unknowns where pi may belongs to {u1,u2, …, un} but it is not mandatory the result will be always a (v1,v2, …, vm)-TermVector. Only product of matrix block (vi,uj) with vector block pk=uj will be performed For instance: the product of a (u,v)-Termatrix M with a p-TermVector X is a zero v-TermVector: |Mvu 0vp| * [0u Xp]t = [0v] the product of a [(u,p),(v,q)]-Termatrix M with a p-TermVector X is a (v,q)-TermVector: |Mvu Mvp| |0u| |Mvp*Xp| | | | | = | | |Mqu Mqp| |Xp| |Mqp*Xp| Note that the zero blocks are never created; they are not indexed in the map TermVector::suTerms_ ! The real products are made with SuTermMatrix and SuTermVector
-
friend TermVector operator*(const TermVector&, const TermMatrix&)#
-
product TermVector * TermMatrix
-
friend TermVector prepareLinearSystem(TermMatrix&, TermVector&, MatrixEntry*&, VectorEntry*&, StorageType, AccessType, bool)#
-
prepare linear system AX=B (internal tool)
-
friend TermMatrix real(const TermMatrix &tm)#
-
return real part as a real TermMatrix
-
friend TermMatrix roundToZero(const TermMatrix &tm, real_t aszero)#
-
return the 0 rounded TermMatrix)
return the 0 rounded TermMatrix
-
friend TermMatrix toComplex(const TermMatrix &tm)#
-
return the complex representation of the given TermMatrix
-
friend void updateRhs(TermMatrix&, TermVector&)#
-
prepare rhs B of linear system AX=B (internal tool)
-
inline TermMatrix(const LcTerm<TermMatrix> &lctm)#