Class xlifepp::SuTermMatrix#
-
class SuTermMatrix : public xlifepp::Term#
-
Inheritence diagram for xlifepp::SuTermMatrix:
Collaboration diagram for xlifepp::SuTermMatrix:
handles numerical representation of single unknowns bilinear forms (SuBilinearForm).
Internal class
Public Functions
-
SuTermMatrix(const SuTermMatrix&, const string_t &na = "")#
-
copy constructor
-
SuTermMatrix(const SuTermMatrix&, SpecialMatrix, const complex_t& = 1., const string_t &na = "")#
-
construct a special matrix from a SuTermMatrix (e.g _IdMatrix)
construct special matrix from a structure given by a TermMatrix A: the SuTermMatrix model sm: special matrix type (up to now, only _idMatrix is available) a: coefficient to apply na: option name
-
SuTermMatrix(const Unknown&, const GeomDomain&, const Unknown&, const GeomDomain&, const OperatorOnKernel&, const string_t &na = "")#
-
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 ux, domx is associated to matrix rows and uy, domy to matrix columns NOTE: do not use this constructor if kernel is singular on some couple of points
-
template<typename T>
SuTermMatrix(const Unknown*, const Space*, const Unknown*, const Space*, const std::vector<T>&, number_t, number_t, const string_t& = "")#
-
constructor from a row dense matrix
constructor from a row dense matrix (only scalar entry)
-
SuTermMatrix(const Unknown*, Space*, const Unknown*, Space*, MatrixEntry*, const string_t &na = "")#
-
basic constructor
-
SuTermMatrix(const Unknown*, Space*, const Unknown*, Space*, SuTermVector&, StorageType = _noStorage, AccessType = _noAccess, const string_t &na = "")#
-
constructor of a diagonal SuTermMatrix from a SuTermVector
-
SuTermMatrix(SuBilinearForm*, const Unknown*, const Unknown*, Space*, Space*, const std::vector<Space*>&, const std::vector<Space*>&, const string_t&, MatrixEntry* = nullptr)#
-
full constructor
-
SuTermMatrix(SuBilinearForm *subf, const string_t &na = "", ComputingInfo cp = ComputingInfo())#
-
basic constructor
-
SuTermMatrix(SuBilinearForm *subf = nullptr, const string_t &na = "", bool noass = false)#
-
basic constructor
-
SuTermMatrix(SuTermVector&, const string_t &na = "")#
-
constructor of a diagonal SuTermMatrix from a SuTermVector
constructor of a diagonal SuTermMatrix from a SuTermVector the diagonal Matrix is stored using SymCsStorage Note: SuTermVector cannot be passed as const because of a cast of the space pointer
-
SuTermMatrix(SuTermVector&, SuTermVector&, const SymbolicFunction &sf, const string_t &na = "")#
-
constructor of a diagonal SuTermMatrix from a SuTermVector
-
virtual ~SuTermMatrix()#
-
destructor
-
AccessType accessType() const#
-
return access type of submatrix (_row, _col, _sym, _dual)
-
MatrixEntry *actual_entries() const#
-
return the actual pointer to entries (priority to scalar entry) access to rhs matrix pointer (r/w)
-
void addDenseToStorage(const std::vector<SuBilinearForm>&, MatrixStorage*) const#
-
add dense matrix indexes to storage (extension, double intg)
-
void addDenseToStorageFE(const std::vector<SuBilinearForm>&, MatrixStorage*) const#
-
add dense matrix indexes to storage (FE with domain restriction)
-
void addDGToStorage(const std::vector<SuBilinearForm>&, MatrixStorage*&) const#
-
add DG matrix indexes to storage
-
void buildFromSuTermVectorPair(const Unknown*, Space*, const Unknown*, Space*, SuTermVector&, SuTermVector&, const SymbolicFunction&, StorageType = _noStorage, AccessType = _noAccess, const string_t &na = "")#
-
actual constructor of diagonal SuTermMatrix from a SuTermVector
-
void buildSubspaces()#
-
build the required subspaces and largest subspaces
-
inline std::vector<DofComponent> &cdofsu()#
-
access to cdofs_u vector (r/w)
-
inline std::vector<DofComponent> &cdofsv()#
-
access to cdofs_v vector (r/w)
-
inline const std::vector<DofComponent> &cdofsv() const#
-
access to cdofs_v vector (r)
-
void changeTestFunction(const Unknown &newv)#
-
change TestFuncion (row unknown) to an other (be cautious)
change test function (csay row unknown) to an other be cautious: nothing is recomputed, bilinear form is not updated, numbering remains the same.
Be sure what you are doing !!! emit warnings if something seems to be hazardous
-
void changeUnknown(const Unknown &newu)#
-
change unknown (col unknown) to an other (be cautious)
change unknown (col unknown) to an other be cautious: nothing is recomputed, bilinear form is not updated, numbering remains the same.
Be sure what you are doing !!! emit warnings if something seems to be hazardous
-
void changeUnknowns(const Unknown &newu, const Unknown &newv)#
-
change both first Unknown and TestFuncion to others (be cautious)
change unknown and test function to others be cautious: nothing is recomputed, bilinear form is not updated, numbering remains the same.
Be sure what you are doing !!! emit warnings if something seems to be hazardous
-
virtual void clear()#
-
deallocate memory used by a SuTermMatrix
-
void clearScalar()#
-
deallocate memory used by scalar part of SuTermMatrix
-
number_t colRank(const DofComponent&) const#
-
rank (>=1) in column of a DofComponent
-
virtual void compute()#
-
compute SuTermMatrix from bilinear form compute SuTermMatrix from a linear combination of SuTermMatrix’s
compute the term matrix from bilinear form and constraints data (possibly none) the process works as follows:
collect bilinear forms by type (FE, FE requiring domain extension, SP, mixed FE-SP, IE, mixed IE-SP)
if not imposed, find a suited storage with the following rule: compressed sparse for FE, else dense
construct the storage
create matrix entry with the right type
call computation routines for each collection of bilinear forms
-
void compute(const LcTerm<SuTermMatrix>&, const string_t &na = "")#
-
compute SuTermMatrix from a linear combination of SuTermMatrix’s (same unknowns) the current SuTermMatrix is assumed to be void, this function fills every things if it is not the case, every thing is clear step 1 : find structure and value type of the linear combination move SuTermMatrix’s entries to entries_p if not allocated (means move to vector representation) find largest spaces choose storage type (type, access, symmetry) step 2 : if dense storage create a new one from largest spaces if not, merge SutermMatrix storage building the colIndices (to be improved) and then create a new storage step 3 : do the linear combination step 4 : update info and go back to original representation of SuTermMatrix’s
-
template<number_t CM>
void compute(const std::vector<SuBilinearForm>&, ValueType, StrucType)#
-
utility for computing SuTermMatrix
-
template<typename T, typename K>
void computeIR(const SuLinearForm &sulf, T *mat, K &vt, const std::vector<Point> &xs, const std::vector<Vector<real_t>> *ns = nullptr)#
-
compute integral representation returns single unknown bilinear forms
FE computation of a scalar SuLinearForm with kernels on a == unique domain == Rij= intg_gamma opk(xi,y)*op(wj)(y) dy.
example: intg_gamma G(x,y)*U(y) dy
Note: U comes from a SuTermVector having same unknown as SuLinearForm
- Parameters:
-
sulf – a single unknown linear form defined on a unique domain (assumed)
vt – to pass as template the scalar type (real or complex)
xs – points where to evaluate IR
nxs – pointers to normal where to evaluate IR
mat – pointer to a matrix stored as row dense, storing the result Rij
-
void copy(const SuTermMatrix&)#
-
copy SuTermMatrix in current SuTermMatrix
-
void diagFromSuTermVector(const Unknown*, Space*, const Unknown*, Space*, SuTermVector&, StorageType = _noStorage, AccessType = _noAccess, const string_t &na = "")#
-
actual constructor of diagonal SuTermMatrix from a SuTermVector
-
inline const MatrixEntry *entries() const#
-
access to entries pointer (r)
-
FactorizationType factorization() const#
-
return factorization type if factorized
-
template<typename T>
HMatrix<T, FeDof> &getHMatrix(StrucType str = _undefStrucType) const#
-
get HMatrix representation if it exists access to cdofs_u vector (r)
get HMatrix representation if it exists, it returns a reference to the HMatrix if scalar and vector representations both exist, the scalar representation is prior except if StrucType argument is specified
-
template<typename T>
LargeMatrix<T> &getLargeMatrix(StrucType str = _undefStrucType) const#
-
get LargeMatrix representation
get internal matrix (restrictive), i.e LargeMatrix object if available if scalar and vector representations both exist, the scalar representation is prior except if StrucType argument is specified It returns a reference to the LargeMatrix
-
Value getScalarValue(number_t, number_t, dimen_t = 0, dimen_t = 0) const#
-
access to scalar matrix coefficient (i,j >= 1)
access to the matrix coefficient (i,j) >= 1 when a matrix of of m x n matrices if(k!=0) returns A(i,j)(k,l) so i,j are indexes in the vector numberings ! else i,j are indexes in the scalar numberings and function returns A(I,J)(k,l) with I=(i-1)/m + 1, k=(i-1)m +1 , J=(j-1)/n + 1, l=(j-1)n +1
-
Value getValue(number_t, number_t) const#
-
access to matrix coefficient (i,j >= 1)
access to the matrix coefficient (i,j) >= 1 in the scalar/vector row-col numbering return a scalar value if a scalar matrix, a matrix value if a matrix of matrices see getScalarValue to get scalar coefficient if (i,j) is out of storage return 0
-
bool hasDomainRestriction(const std::vector<SuBilinearForm> &subfs) const#
-
check if there are some FE subfs with domain restriction
-
void initPointers()#
-
set all pointers to 0 initialize a SuTermVector consistent with matrix rows or columns
-
inline virtual void name(const string_t &nam)#
-
update the name of a SuTermMatrix
-
real_t norm2() const#
-
Frobenius norm of a SuTermMatrix.
-
real_t norminfty() const#
-
infinite norm of a SuTermMatrix
-
template<typename T>
SuTermMatrix &operator*=(const T&)#
-
SutermMatrix *= t.
-
SuTermMatrix &operator+=(const SuTermMatrix&)#
-
SutermMatrix += SutermMatrix.
Algebraic operation += and -= the result of these operations is a SuTermMatrix with extended dof numbering, may require a new storage !!! works with entries prior to scalar entries.
-
SuTermMatrix &operator-=(const SuTermMatrix&)#
-
SutermMatrix -= SutermMatrix.
-
template<typename T>
SuTermMatrix &operator/=(const T&)#
-
SutermMatrix /= t.
-
SuTermMatrix &operator=(const SuTermMatrix&)#
-
assign operator constructor from a linear combination of SuTermMatrix’s
-
virtual void print(std::ostream&) const#
-
print SuTermMatrix
-
void print(std::ostream&, bool) const#
-
print SuTermMatrix with header option
-
void printSummary(std::ostream&) const#
-
print SuTermMatrix in brief
-
inline const MatrixEntry *rhs_matrix() const#
-
access to rhs matrix pointer (r)
-
SuTermMatrix &roundToZero(real_t aszero = 10 * theEpsilon)#
-
round to zero all coefficients close to 0 (|.
| < aszero)
-
number_t rowRank(const DofComponent&) const#
-
rank (>=1) in row of a DofComponent
-
void saveToFile(const string_t &fn, StorageType st, number_t prec = fullPrec, bool encode = false, real_t tol = theTolerance) const#
-
save TermMatrix to file (dense or Matlab format) save TermMatrix to file (dense or Matlab format)
-
inline MatrixEntry *&scalar_entries()#
-
access to entries pointer (r/w)
-
inline const MatrixEntry *scalar_entries() const#
-
access to entries pointer (r)
-
MatrixStorage *scalarStoragep() const#
-
return scalar storage of submatrix (0 if no scalar entries or storage)
-
void setCol(const Value&, number_t r1, number_t r2)#
-
set values of matrix colss c1->c2 (>= 1)
set values of a matrix cols c1->c2 (>= 1), c2=0 means nbrows, no storage change!
-
void setRow(const Value&, number_t r1, number_t r2)#
-
set values of matrix rows r1->r2 (>= 1)
set values of a matrix rows r1->r2 (>= 1), c2=0 means nbrows, value has to be consistant with matrix structure no storage change
-
void setScalarValue(number_t, number_t, const Value&, dimen_t = 0, dimen_t = 0)#
-
set matrix coefficient (i,j) >= 1 from a value when a matrix of of m x n matrices if(k!=0) set A(i,j)(k,l) so i,j are indexes in the vector numberings ! else i,j are indexes in the scalar numberings and function sets A(I,J)(k,l) with I=(i-1)/m + 1, k=(i-1)m +1 , J=(j-1)/n + 1, l=(j-1)n +1
-
void setStorage(StorageType, AccessType)#
-
set the storage of SuTermMatrix
-
void setValue(number_t, number_t, const Value&)#
-
set value to matrix coefficient (i,j >= 1) in storage set value to matrix coefficient (i,j >= 1) in storage
set value to matrix coefficient (i,j >= 1) in storage, Value must be of the same type of the matrix coefficient (i,j) must be in the matrix storage else an error occurs see getScalarValue to set a scalar coefficient of a matrix of matrices
-
inline Space *&space_up()#
-
return the pointer to the largest space in u involved in SuBilinearForm (write)
-
inline Space *space_up() const#
-
return the pointer to the largest space in u involved in SuBilinearForm
-
inline Space &space_v() const#
-
return the largest space in v involved in SuBilinearForm
-
inline Space *&space_vp()#
-
return the pointer to the largest space in v involved in SuBilinearForm (write)
-
inline Space *space_vp() const#
-
return the pointer to the largest space in v involved in SuBilinearForm
-
MatrixStorage *storagep() const#
-
return storage of submatrix (0 if no entries or storage)
-
StorageType storageType() const#
-
return storage type of submatrix (_cs,_skyline,_dense)
-
StrucType strucType() const#
-
return structure type (_scalar, _matrix) return the largest space in u involved in SuBilinearForm
-
inline SuBilinearForm *&subfp()#
-
return pointer to bilinear form
-
SuTermMatrix &toComplex()#
-
go to complex representation
-
SuTermMatrix &toConj()#
-
change to conj(U)
-
SuTermMatrix &toImag()#
-
go to real representation getting the imag part when complex
-
SuTermMatrix &toReal()#
-
go to real representation getting the real part when complex
-
void toScalar(bool keepEntries = false)#
-
go to scalar representation
create scalar representation scalar_entries_p and cdofs_u, cdofs_v numbering vectors
if one of unknowns is a vector unknown a new MatrixEntry is created using MatrixEntry::toScalar() function and linked to the MatrixEntry pointer scalar_entries_p. In that case, if keepVector is false the vector representation is destroyed
if both unknowns are scalar unknown, no new MatrixEntry: scalar_entries_p = entries_p This function builds also, in any case, the cdofs_u, cdofs_v numbering vectors (vector of ComponentDof’s)
-
void toStorage(StorageType, AccessType)#
-
change the storage of SuTermMatrix (real algorithm)
-
void toVectorUnknown()#
-
go to vector unknown representation if it has component unknown
move to vector unknown representation if SuTermMatrix has component unknown Note: the SutermMatrix is modified its unknowns too
-
virtual std::set<const Space*> unknownSpaces() const#
-
list of involved unknown spaces access to entries pointer (r/w)
-
void updateStorageType(const std::vector<SuBilinearForm>&, std::set<number_t>&, std::set<number_t>&, StorageType&) const#
-
update storage type
-
void viewStorage(std::ostream&) const#
-
print storage on stream
Friends
-
friend SuTermMatrix *mergeSuTermMatrix(const std::list<SuTermMatrix*>&)#
-
merge SuTerMatrix’s referring to the same vector unknown
-
friend SuTermMatrix operator*(const SuTermMatrix&, const SuTermMatrix&)#
-
product SuTermMatrix * SuTermMatrix
-
friend SuTermVector operator*(const SuTermMatrix&, const SuTermVector&)#
-
Product of a SuTermMatrix M and a SuTermVector V, the result is a SuTermVector R the columns dof indices of M, say (j1,j2, …jn), may be different from the dof indices of V say (k1,k2, … kp).
product SuTermMatrix * SuTermVector
The dof indices of result vector will be always the lines dof indices of the matrix, say (i1,i2, …im). It means that the product M_i,j * V_j is performed only for j=k.
when {k1,k2, … kp}={j1,j2, …jn} the product is performed using largematrix product
when {k1,k2, … kp} differs from {j1,j2, …jn}, the vector V is extended and reduced to {j1,j2, …jn} and the product is performed using largematrix product
-
friend SuTermVector operator*(const SuTermVector&, const SuTermMatrix&)#
-
Product of a SuTermVector V and a SuTermMatrix M, the result is a SuTermVector R the lines dof indices of M, say (i1,i2, …in), may be different from the dof indices of V say (k1,k2, … kp).
product SuTermVector * SuTermMatrix
The dof indices of result vector will be always the columns dof indices of the matrix, say (j1,j2, …jm). It means that the product V_i * M_i,j is performed only for i=k.
when {k1,k2, … kp}={i1,i2, …in} the product is performed using largematrix product
when {k1,k2, … kp} differs from {i1,i2, …in}, the vector V is extended and reduced to {i1,i2, …in} and the product is performed using largematrix product
exists only for particular storage type
-
SuTermMatrix(const SuTermMatrix&, const string_t &na = "")#