NeoPZ
|
Full matrix class. Matrix. More...
#include <pzfmatrix.h>
Public Member Functions | |
TPZFMatrix () | |
Simple constructor. More... | |
TPZFMatrix (const int64_t rows, const int64_t columns, TVar *buf, const int64_t size) | |
Constructor with initialization parameters. More... | |
TPZFMatrix (const int64_t rows, const int64_t columns, const TVar &val) | |
Constructor with initialization parameters. More... | |
TPZFMatrix (const int64_t rows, const int64_t columns=1) | |
Constructor with initialization parameters. More... | |
TPZFMatrix (TPZVerySparseMatrix< TVar > const &A) | |
Copy constructor specialized form TPZVerySparseMatrix. More... | |
TPZFMatrix (const TPZFMatrix< TVar > &refmat) | |
Copy constructor. More... | |
TPZFMatrix (const TPZMatrix< TVar > &refmat) | |
TPZFMatrix (const std::initializer_list< TVar > &list) | |
Creates a matrix from initializer list (one column matrix) More... | |
TPZFMatrix (const std::initializer_list< std::initializer_list< TVar > > &list) | |
Creates a matrix from initializer list. More... | |
virtual | ~TPZFMatrix () |
Simple destructor. More... | |
int64_t | MemoryFootprint () const override |
Returns the approximate size of the memory footprint (amount of memory required to store this object). More... | |
TVar * | Adress () |
template<class TVar2 > | |
void | CopyFrom (TPZFMatrix< TVar2 > &orig) |
copy the values from a matrix with a different precision More... | |
virtual void | UpdateFrom (TPZAutoPointer< TPZMatrix< TVar > > mat) override |
Updates the values of the matrix based on the values of the matrix. More... | |
int | PutVal (const int64_t row, const int64_t col, const TVar &value) override |
Put values without bounds checking This method is faster than "Put" if DEBUG is defined. More... | |
const TVar & | GetVal (const int64_t row, const int64_t col) const override |
Get values without bounds checking This method is faster than "Get" if DEBUG is defined. More... | |
virtual TVar & | s (const int64_t row, const int64_t col) override |
The operators check on the bounds if the DEBUG variable is defined. More... | |
TVar & | g (const int64_t row, const int64_t col) const |
void | AddFel (TPZFMatrix< TVar > &rhs, TPZVec< int64_t > &destination) |
Performs a right hand side assemblage. More... | |
void | AddFel (TPZFMatrix< TVar > &rhs, TPZVec< int64_t > &source, TPZVec< int64_t > &destination) |
Performs a right hand side assemblage. More... | |
virtual void | MultAdd (const TPZFMatrix< TVar > &x, const TPZFMatrix< TVar > &y, TPZFMatrix< TVar > &z, const TVar alpha=1., const TVar beta=0., const int opt=0) const override |
It computes z = beta * y + alpha * opt(this)*x but z and x can not overlap in memory. More... | |
void | ZAXPY (const TVar alpha, const TPZFMatrix< TVar > &p) |
Performs an ZAXPY operation being *this += alpha * p. More... | |
void | TimesBetaPlusZ (const TVar beta, const TPZFMatrix< TVar > &z) |
Performs an operation *this = this * beta + z. More... | |
TPZFMatrix< TVar > & | operator= (const TPZMatrix< TVar > &A) |
Generic operator with matrices. More... | |
int | Resize (const int64_t newRows, const int64_t wCols) override |
Redimension a matrix, but maintain your elements. More... | |
int | SetSize (int64_t newRows, int64_t newCols) |
Redimension the matrix doing nothing with the elements. More... | |
int | Remodel (const int64_t newRows, const int64_t wCols) |
Remodel the shape of the matrix, but keeping the same dimension. More... | |
int | Redim (const int64_t newRows, const int64_t newCols) override |
Redimension a matrix and ZERO your elements. More... | |
int | Zero () override |
Makes Zero all the elements. More... | |
void | GramSchmidt (TPZFMatrix< TVar > &Orthog, TPZFMatrix< TVar > &TransfToOrthog) |
This method implements a Gram Schimidt method. this = Orthog.TransfToOrthog. More... | |
void | DeterminantInverse (TVar &determinant, TPZFMatrix< TVar > &inverse) |
void | Transpose (TPZMatrix< TVar > *const T) const override |
It makes *T the transpose of current matrix. More... | |
void | Transpose () |
int | ClassId () const override |
Routines to send and receive messages. More... | |
void | Read (TPZStream &buf, void *context) override |
read objects from the stream More... | |
void | Write (TPZStream &buf, int withclassid) const override |
Writes this object to the TPZStream buffer. Include the classid if withclassid = true. More... | |
virtual bool | Compare (TPZSavable *copy, bool override=false) override |
Compare the object for identity with the object pointed to, eventually copy the object. More... | |
virtual bool | Compare (TPZSavable *copy, bool override=false) const override |
Compare the object for identity with the object pointed to, eventually copy the object. More... | |
operator const TVar * () const | |
template<> | |
void | AddFel (TPZFMatrix< double > &rhs, TPZVec< int64_t > &source, TPZVec< int64_t > &destination) |
template<> | |
void | AddFel (TPZFMatrix< float > &rhs, TPZVec< int64_t > &source, TPZVec< int64_t > &destination) |
template<> | |
void | GramSchmidt (TPZFMatrix< TPZFlopCounter > &Orthog, TPZFMatrix< TPZFlopCounter > &TransfToOrthog) |
Generic operator with TVar type | |
TVar & | operator() (const int64_t row, const int64_t col) |
TVar & | operator() (const int64_t row) |
Operations with FULL matrices | |
virtual TPZFMatrix & | operator= (const TPZFMatrix< TVar > &A) |
Generic operator with FULL matrices. More... | |
TPZFMatrix< TVar > & | operator= (const std::initializer_list< TVar > &list) |
TPZFMatrix< TVar > & | operator= (const std::initializer_list< std::initializer_list< TVar > > &list) |
TPZFMatrix< TVar > | operator+ (const TPZFMatrix< TVar > &A) const |
TPZFMatrix< TVar > | operator- (const TPZFMatrix< TVar > &A) const |
TPZFMatrix< TVar > | operator* (TPZFMatrix< TVar > A) const |
TPZFMatrix< TVar > & | operator+= (const TPZFMatrix< TVar > &A) |
TPZFMatrix< TVar > & | operator-= (const TPZFMatrix< TVar > &A) |
Numerics | |
Numeric operations with matrices | |
TPZFMatrix< TVar > & | operator= (const TVar val) |
Numeric operator with matrices. More... | |
TPZFMatrix< TVar > | operator+ (const TVar val) const |
TPZFMatrix< TVar > | operator- (const TVar val) const |
TPZFMatrix< TVar > | operator* (const TVar val) const |
TPZFMatrix< TVar > & | operator+= (const TVar val) |
TPZFMatrix< TVar > & | operator-= (const TVar val) |
TPZFMatrix< TVar > & | operator*= (const TVar val) |
Public Member Functions inherited from TPZMatrix< TVar > | |
TPZMatrix () | |
Simple constructor. More... | |
TPZMatrix (const TPZMatrix< TVar > &cp) | |
virtual | ~TPZMatrix () |
Simple destructor. More... | |
virtual TPZMatrix< TVar > * | Clone () const =0 |
template<class TVar2 > | |
void | CopyFrom (TPZMatrix< TVar2 > ©) |
void | AutoFill (int64_t nrow, int64_t ncol, int symmetric) |
Fill matrix storage with randomic values. More... | |
virtual int | VerifySymmetry (REAL tol=1.e-13) const |
Checks if current matrix value is symmetric. More... | |
virtual int | Put (const int64_t row, const int64_t col, const TVar &value) |
Put values with bounds checking if DEBUG variable is defined. More... | |
virtual const TVar & | Get (const int64_t row, const int64_t col) const |
Get value with bound checking. More... | |
const TVar & | g (const int64_t row, const int64_t col) const |
Substitution for the () operator when const arguments are needed. More... | |
TVar & | operator() (const int64_t row, const int64_t col) |
The operators check on the bounds if the DEBUG variable is defined. More... | |
TVar & | operator() (const int64_t row) |
The operators check on the bounds if the DEBUG variable is defined. More... | |
virtual void | AddKel (TPZFMatrix< TVar > &elmat, TPZVec< int64_t > &destinationindex) |
Add a contribution of a stiffness matrix. More... | |
virtual void | AddKel (TPZFMatrix< TVar > &elmat, TPZVec< int64_t > &sourceindex, TPZVec< int64_t > &destinationindex) |
Add a contribution of a stiffness matrix. More... | |
void | SetIsDecomposed (int val) |
Sets current matrix to decomposed state. More... | |
virtual void | GetSub (const TPZVec< int64_t > &indices, TPZFMatrix< TVar > &block) const |
Extract the block indicated by the indices from the matrix. More... | |
bool | CompareValues (TPZMatrix< TVar > &M, TVar tol) |
Compare values of this to B, with a precision tolerance tol. More... | |
template<> | |
void | Print (const char *name, std::ostream &out, const MatrixOutputFormat form) const |
template<> | |
void | Print (const char *name, std::ostream &out, const MatrixOutputFormat form) const |
template<> | |
void | Print (const char *name, std::ostream &out, const MatrixOutputFormat form) const |
template<> | |
void | SolveCG (int64_t &numiterations, TPZSolver< std::complex< float > > &preconditioner, const TPZFMatrix< std::complex< float > > &F, TPZFMatrix< std::complex< float > > &result, TPZFMatrix< std::complex< float > > *residual, REAL &tol, const int FromCurrent) |
template<> | |
void | SolveCG (int64_t &numiterations, TPZSolver< std::complex< double > > &preconditioner, const TPZFMatrix< std::complex< double > > &F, TPZFMatrix< std::complex< double > > &result, TPZFMatrix< std::complex< double > > *residual, REAL &tol, const int FromCurrent) |
template<> | |
void | SolveCG (int64_t &numiterations, TPZSolver< std::complex< long double > > &preconditioner, const TPZFMatrix< std::complex< long double > > &F, TPZFMatrix< std::complex< long double > > &result, TPZFMatrix< std::complex< long double > > *residual, REAL &tol, const int FromCurrent) |
template<> | |
void | SolveGMRES (int64_t &numiterations, TPZSolver< std::complex< float > > &preconditioner, TPZFMatrix< std::complex< float > > &H, int &numvectors, const TPZFMatrix< std::complex< float > > &F, TPZFMatrix< std::complex< float > > &result, TPZFMatrix< std::complex< float > > *residual, REAL &tol, const int FromCurrent) |
template<> | |
void | SolveGMRES (int64_t &numiterations, TPZSolver< std::complex< double > > &preconditioner, TPZFMatrix< std::complex< double > > &H, int &numvectors, const TPZFMatrix< std::complex< double > > &F, TPZFMatrix< std::complex< double > > &result, TPZFMatrix< std::complex< double > > *residual, REAL &tol, const int FromCurrent) |
template<> | |
void | SolveGMRES (int64_t &numiterations, TPZSolver< std::complex< long double > > &preconditioner, TPZFMatrix< std::complex< long double > > &H, int &numvectors, const TPZFMatrix< std::complex< long double > > &F, TPZFMatrix< std::complex< long double > > &result, TPZFMatrix< std::complex< long double > > *residual, REAL &tol, const int FromCurrent) |
template<> | |
void | SolveBICG (int64_t &numiterations, TPZSolver< std::complex< float > > &preconditioner, const TPZFMatrix< std::complex< float > > &F, TPZFMatrix< std::complex< float > > &result, REAL &tol) |
template<> | |
void | SolveBICG (int64_t &numiterations, TPZSolver< std::complex< double > > &preconditioner, const TPZFMatrix< std::complex< double > > &F, TPZFMatrix< std::complex< double > > &result, REAL &tol) |
template<> | |
void | SolveBICG (int64_t &numiterations, TPZSolver< std::complex< long double > > &preconditioner, const TPZFMatrix< std::complex< long double > > &F, TPZFMatrix< std::complex< long double > > &result, REAL &tol) |
template<> | |
void | SolveBICGStab (int64_t &numiterations, TPZSolver< std::complex< float > > &preconditioner, const TPZFMatrix< std::complex< float > > &F, TPZFMatrix< std::complex< float > > &result, TPZFMatrix< std::complex< float > > *residual, REAL &tol, const int FromCurrent) |
template<> | |
void | SolveBICGStab (int64_t &numiterations, TPZSolver< std::complex< double > > &preconditioner, const TPZFMatrix< std::complex< double > > &F, TPZFMatrix< std::complex< double > > &result, TPZFMatrix< std::complex< double > > *residual, REAL &tol, const int FromCurrent) |
template<> | |
void | SolveBICGStab (int64_t &numiterations, TPZSolver< std::complex< long double > > &preconditioner, const TPZFMatrix< std::complex< long double > > &F, TPZFMatrix< std::complex< long double > > &result, TPZFMatrix< std::complex< long double > > *residual, REAL &tol, const int FromCurrent) |
template<> | |
void | SolveIR (int64_t &numiterations, TPZSolver< std::complex< float > > &preconditioner, const TPZFMatrix< std::complex< float > > &F, TPZFMatrix< std::complex< float > > &result, TPZFMatrix< std::complex< float > > *residual, REAL &tol, const int FromCurrent) |
template<> | |
void | SolveIR (int64_t &numiterations, TPZSolver< std::complex< double > > &preconditioner, const TPZFMatrix< std::complex< double > > &F, TPZFMatrix< std::complex< double > > &result, TPZFMatrix< std::complex< double > > *residual, REAL &tol, const int FromCurrent) |
template<> | |
void | SolveIR (int64_t &numiterations, TPZSolver< std::complex< long double > > &preconditioner, const TPZFMatrix< std::complex< long double > > &F, TPZFMatrix< std::complex< long double > > &result, TPZFMatrix< std::complex< long double > > *residual, REAL &tol, const int FromCurrent) |
template<> | |
bool | SolveEigenvaluesJacobi (int64_t &numiterations, REAL &tol, TPZVec< std::complex< float > > *Sort) |
template<> | |
bool | SolveEigenvaluesJacobi (int64_t &numiterations, REAL &tol, TPZVec< std::complex< double > > *Sort) |
template<> | |
bool | SolveEigenvaluesJacobi (int64_t &numiterations, REAL &tol, TPZVec< std::complex< long double > > *Sort) |
virtual void | Multiply (const TPZFMatrix< TVar > &A, TPZFMatrix< TVar > &res, int opt=0) const |
It mutiplies itself by TPZMatrix<TVar>A putting the result in res. More... | |
virtual void | Add (const TPZMatrix< TVar > &A, TPZMatrix< TVar > &res) const |
It adds itself to TPZMatrix<TVar>A putting the result in res. More... | |
virtual void | Residual (const TPZFMatrix< TVar > &x, const TPZFMatrix< TVar > &rhs, TPZFMatrix< TVar > &res) |
Computes res = rhs - this * x. More... | |
virtual void | Substract (const TPZMatrix< TVar > &A, TPZMatrix< TVar > &result) const |
It substracts A from storing the result in result. More... | |
virtual void | Identity () |
Converts the matrix in an identity matrix. More... | |
int | Inverse (TPZFMatrix< TVar > &Inv, DecomposeType dec) |
It makes Inv =[this]. IMPORTANT OBSERVATION –> The original matrix (calling object) no is more equal. It containts the some decomposition (LU or Cholesky or ...) More... | |
TVar | MatrixNorm (int p, int64_t numiter=2000000, REAL tol=1.e-10) const |
Computes the matrix norm of this. More... | |
TVar | ConditionNumber (int p, int64_t numiter=2000000, REAL tol=1.e-10) |
Computes the matrix condition number of this. More... | |
virtual int | PutSub (const int64_t sRow, const int64_t sCol, const TPZFMatrix< TVar > &Source) |
It puts submatrix Source on actual matrix structure. More... | |
virtual int | GetSub (const int64_t sRow, const int64_t sCol, const int64_t rowSize, const int64_t colSize, TPZFMatrix< TVar > &Target) const |
Gets submatrix storing it on Target. More... | |
virtual int | AddSub (const int64_t sRow, const int64_t sCol, const TPZFMatrix< TVar > &Source) |
It adds Source matrix on current matrix from position (sRow, sCol) More... | |
virtual int | InsertSub (const int64_t sRow, const int64_t sCol, const int64_t rowSize, const int64_t colSize, const int64_t pRow, const int64_t pCol, TPZMatrix< TVar > *Target) const |
Inserts a submatrix from current object on matrix *Target with no redimentioning. More... | |
virtual int | AddSub (const int64_t sRow, const int64_t sCol, const int64_t rowSize, const int64_t colSize, const int64_t pRow, const int64_t pCol, TPZMatrix< TVar > *pA) const |
Adds a submatrix from current object in *Target. More... | |
virtual int | IsSimetric () const |
Checks if the current matrix is symmetric. More... | |
int | IsSquare () const |
Checks if current matrix is square. More... | |
virtual void | Simetrize () |
Simetrizes copies upper plan to the lower plan, making its data simetric. More... | |
virtual int | IsDefPositive () const |
Checks if current matrix is definite positive. More... | |
int | IsDecomposed () const |
Checks if current matrix is already decomposed. More... | |
virtual int | Decompose_LDLt (std::list< int64_t > &singular) |
Decomposes the current matrix using LDLt. The current matrix has to be symmetric. "L" is lower triangular with 1.0 in its diagonal and "D" is a Diagonal matrix. More... | |
virtual int | Subst_Forward (TPZFMatrix< TVar > *b) const |
Computes B = Y, where A*Y = B, A is lower triangular. More... | |
virtual int | Subst_Backward (TPZFMatrix< TVar > *b) const |
Computes B = Y, where A*Y = B, A is upper triangular. More... | |
virtual int | Subst_LForward (TPZFMatrix< TVar > *b) const |
Computes B = Y, where A*Y = B, A is lower triangular with A(i,i)=1. More... | |
virtual int | Subst_LBackward (TPZFMatrix< TVar > *b) const |
Computes B = Y, where A*Y = B, A is upper triangular with A(i,i)=1. More... | |
virtual int | Subst_Diag (TPZFMatrix< TVar > *b) const |
Computes B = Y, where A*Y = B, A is diagonal matrix. More... | |
int | ClassId () const override |
Define the class id associated with the class. More... | |
void | Read (TPZStream &buf, void *context) override |
Unpacks the object structure from a stream of bytes. More... | |
void | Write (TPZStream &buf, int withclassid) const override |
Packs the object structure in a stream of bytes. More... | |
virtual void | Input (std::istream &in=std::cin) |
Input operation. More... | |
virtual void | Print (std::ostream &out) const |
virtual void | Print (const char *name, std::ostream &out=std::cout, const MatrixOutputFormat form=EFormatted) const |
It prints the matrix data in a MatrixFormat Rows X Cols. More... | |
int64_t | Rows () const |
Returns number of rows. More... | |
int64_t | Cols () const |
Returns number of cols. More... | |
virtual int64_t | Dim () const |
Returns the dimension of the matrix if the matrix is square. More... | |
virtual void | SolveJacobi (int64_t &numiterations, const TPZFMatrix< TVar > &F, TPZFMatrix< TVar > &result, TPZFMatrix< TVar > *residual, TPZFMatrix< TVar > &scratch, REAL &tol, const int FromCurrent=0) |
Solves the linear system using Jacobi method. . More... | |
virtual void | SolveSOR (int64_t &numiterations, const TPZFMatrix< TVar > &F, TPZFMatrix< TVar > &result, TPZFMatrix< TVar > *residual, TPZFMatrix< TVar > &scratch, const REAL overrelax, REAL &tol, const int FromCurrent=0, const int direction=1) |
Solves the linear system using Successive Over Relaxation method (Gauss Seidel). . More... | |
virtual void | SolveSSOR (int64_t &numiterations, const TPZFMatrix< TVar > &F, TPZFMatrix< TVar > &result, TPZFMatrix< TVar > *residual, TPZFMatrix< TVar > &scratch, const REAL overrelax, REAL &tol, const int FromCurrent=0) |
Solves the linear system using Symmetric Successive Over Relaxation method (Gauss Seidel). . More... | |
virtual void | SolveCG (int64_t &numiterations, TPZSolver< TVar > &preconditioner, const TPZFMatrix< TVar > &F, TPZFMatrix< TVar > &result, TPZFMatrix< TVar > *residual, REAL &tol, const int FromCurrent=0) |
Solves the linear system using Conjugate Gradient method. . More... | |
virtual void | SolveBICG (int64_t &numiterations, TPZSolver< TVar > &preconditioner, const TPZFMatrix< TVar > &F, TPZFMatrix< TVar > &result, REAL &tol) |
Solves the linear system using Bi-Conjugate Gradient method. . More... | |
virtual void | SolveBICGStab (int64_t &numiterations, TPZSolver< TVar > &preconditioner, const TPZFMatrix< TVar > &F, TPZFMatrix< TVar > &result, TPZFMatrix< TVar > *residual, REAL &tol, const int FromCurrent=0) |
Solves the linear system using Bi-Conjugate Gradient stabilized method. . More... | |
virtual void | SolveGMRES (int64_t &numiterations, TPZSolver< TVar > &preconditioner, TPZFMatrix< TVar > &H, int &numvectors, const TPZFMatrix< TVar > &F, TPZFMatrix< TVar > &result, TPZFMatrix< TVar > *residual, REAL &tol, const int FromCurrent) |
Solves the linear system using Generalized Minimal Residual (GMRES) method. . More... | |
virtual void | SolveIR (int64_t &numiterations, TPZSolver< TVar > &preconditioner, const TPZFMatrix< TVar > &F, TPZFMatrix< TVar > &result, TPZFMatrix< TVar > *residual, REAL &tol, const int FromCurrent=0) |
Solves the linear system using IR method. . More... | |
virtual bool | SolveEigenvaluesJacobi (int64_t &numiterations, REAL &tol, TPZVec< TVar > *Sort=0) |
Transforms this matrix in a diagonal matrix, where the diagonal values are its eigenvalues. This method is efficient only for small matrices. More... | |
virtual bool | SolveEigensystemJacobi (int64_t &numiterations, REAL &tol, TPZVec< TVar > &Eigenvalues, TPZFMatrix< TVar > &Eigenvectors) const |
Compute Eigenvalues and Eigenvectors of this matrix. This method is efficient only for small matrices. More... | |
virtual int | SolveDirect (TPZFMatrix< TVar > &F, const DecomposeType dt, std::list< int64_t > &singular) |
Solves the linear system using Direct methods. More... | |
virtual int | SolveDirect (TPZFMatrix< TVar > &F, const DecomposeType dt) |
Solves the linear system using Direct methods. More... | |
virtual int | Decompose (const DecomposeType dt, std::list< int64_t > &singular) |
decompose the system of equations acording to the decomposition scheme More... | |
int | Solve_LU (TPZFMatrix< TVar > *B, std::list< int64_t > &singular) |
Solves the linear system using LU method . More... | |
int | Solve_LU (TPZFMatrix< TVar > *B) |
Solves the linear system using LU method . More... | |
virtual int | Solve_Cholesky (TPZFMatrix< TVar > *B) |
Solves the linear system using Cholesky method . More... | |
int | Solve_Cholesky (TPZFMatrix< TVar > *B, std::list< int64_t > &singular) |
Solves the linear system using Cholesky method . More... | |
int | Solve_LDLt (TPZFMatrix< TVar > *B, std::list< int64_t > &singular) |
Solves the linear system using LDLt method . More... | |
int | Solve_LDLt (TPZFMatrix< TVar > *B) |
Solves the linear system using LDLt method . More... | |
Public Member Functions inherited from TPZSavable | |
TPZSavable () | |
virtual | ~TPZSavable () |
virtual std::list< std::map< std::string, uint64_t > > | VersionHistory () const |
virtual std::pair< std::string, uint64_t > | Version () const |
Public Member Functions inherited from TPZRegisterClassId | |
template<typename T > | |
TPZRegisterClassId (int(T::*)() const) | |
TPZRegisterClassId ()=default | |
Static Public Member Functions | |
static void | MultAdd (const TVar *ptr, int64_t rows, int64_t cols, const TPZFMatrix< TVar > &x, const TPZFMatrix< TVar > &y, TPZFMatrix< TVar > &z, const TVar alpha=1., const TVar beta=0., const int opt=0) |
static void | PrintStatic (const TVar *ptr, int64_t rows, int64_t cols, const char *name, std::ostream &out, const MatrixOutputFormat form) |
Static Public Member Functions inherited from TPZMatrix< TVar > | |
static int | Error (const char *msg, const char *msg2=0) |
Returns error messages. More... | |
static TVar | ReturnNearestValue (TVar val, TPZVec< TVar > &Vec, TVar tol) |
Retorna o valor mais proximo a "val" (exceto valores no intervalo -tol <= val <= +tol) contido no vetor Vec. More... | |
Static Public Member Functions inherited from TPZSavable | |
static std::set< TPZRestoreClassBase * > & | RestoreClassSet () |
This static function guarantees that the gMap object is available when needed. More... | |
static std::map< int, TPZRestore_t > & | ClassIdMap () |
This static function guarantees that the gMap object is available when needed. More... | |
static std::pair< std::string, uint64_t > | NeoPZVersion () |
static void | Register (TPZRestoreClassBase *restore) |
static void | RegisterClassId (int classid, TPZRestore_t fun) |
static TPZSavable * | CreateInstance (const int &classId) |
Private Member Functions | |
int | Clear () override |
It clears data structure. More... | |
Static Private Member Functions | |
static int | Error (const char *msg1, const char *msg2=0) |
Private Attributes | |
TVar * | fElem |
TVar * | fGiven |
int64_t | fSize |
Friends | |
class | TPZFMatrix< float > |
class | TPZFMatrix< double > |
virtual int | Decompose_Cholesky () override |
Cholesky Decomposition Optmized. for walks in the direction of the vector that composes the matrix. More... | |
virtual int | Decompose_Cholesky (std::list< int64_t > &singular) override |
Decomposes the current matrix using Cholesky method. More... | |
virtual int | Decompose_LU (std::list< int64_t > &singular) override |
LU Decomposition. Stores L and U matrices at the storage of the same matrix. More... | |
virtual int | Decompose_LU () override |
virtual int | Decompose_LDLt () override |
Decomposes the current matrix using LDLt. The current matrix has to be symmetric. "L" is lower triangular with 1.0 in its diagonal and "D" is a Diagonal matrix. More... | |
virtual int | Substitution (TPZFMatrix< TVar > *B) const override |
Computes Forward and Backward substitution for a "LU" decomposed matrix. More... | |
virtual int | Decompose_LU (TPZVec< int > &index) |
LU Decomposition using pivot. More... | |
virtual int | Substitution (TPZFMatrix< TVar > *B, const TPZVec< int > &index) const |
LU substitution using pivot. More... | |
static int | Substitution (const TVar *ptr, int64_t rows, TPZFMatrix< TVar > *B) |
static int | Substitution (const TVar *ptr, int64_t rows, TPZFMatrix< TVar > *B, const TPZVec< int > &index) |
LU substitution using pivot. Static version. More... | |
Additional Inherited Members | |
Protected Member Functions inherited from TPZMatrix< TVar > | |
void | PrepareZ (const TPZFMatrix< TVar > &y, TPZFMatrix< TVar > &z, const TVar beta, const int opt) const |
Is an auxiliar method used by MultiplyAdd. More... | |
TPZMatrix (const int64_t row, const int64_t col) | |
Constructor. More... | |
Static Protected Member Functions inherited from TPZMatrix< TVar > | |
static void | Swap (int64_t *a, int64_t *b) |
Swaps contents of a in b and b in a. More... | |
Protected Attributes inherited from TPZMatrix< TVar > | |
int64_t | fRow |
Number of rows in matrix. More... | |
int64_t | fCol |
Number of cols in matrix. More... | |
char | fDecomposed |
Decomposition type used to decompose the current matrix. More... | |
char | fDefPositive |
Definite Posistiveness of current matrix. More... | |
Static Protected Attributes inherited from TPZMatrix< TVar > | |
static TVar | gZero |
Initializing value to static variable. More... | |
Full matrix class. Matrix.
Definition at line 32 of file pzfmatrix.h.
|
inline |
Simple constructor.
Definition at line 74 of file pzfmatrix.h.
Referenced by TPZFMatrix< STATE >::TPZFMatrix().
|
inline |
Constructor with initialization parameters.
rows | Initial number of rows |
columns | Number of columns |
buf | Preallocated memory area which can be used by the matrix object |
size | Size of the area pointed to by buf |
Definition at line 445 of file pzfmatrix.h.
|
inline |
Constructor with initialization parameters.
rows | Initial number of rows |
columns | Number of columns |
val | Inital value fill all elements |
Definition at line 462 of file pzfmatrix.h.
|
inline |
Constructor with initialization parameters.
rows | Initial number of rows |
columns | Number of columns |
Definition at line 96 of file pzfmatrix.h.
TPZFMatrix< TVar >::TPZFMatrix | ( | TPZVerySparseMatrix< TVar > const & | A | ) |
Copy constructor specialized form TPZVerySparseMatrix.
refmat | Used as a model for current object |
Definition at line 99 of file pzfmatrix.cpp.
TPZFMatrix< TVar >::TPZFMatrix | ( | const TPZFMatrix< TVar > & | refmat | ) |
Copy constructor.
refmat | Used as a model for current object |
Definition at line 80 of file pzfmatrix.cpp.
TPZFMatrix< TVar >::TPZFMatrix | ( | const TPZMatrix< TVar > & | refmat | ) |
Definition at line 59 of file pzfmatrix.cpp.
|
inline |
Creates a matrix from initializer list (one column matrix)
list | : initializer list, usually a list of elements in curly brackets |
Definition at line 475 of file pzfmatrix.h.
|
inline |
Creates a matrix from initializer list.
list | : initializer list, usually a list of elements in curly brackets |
|
inlinevirtual |
Simple destructor.
Definition at line 556 of file pzfmatrix.h.
Referenced by TPZFMatrix< STATE >::TPZFMatrix().
void TPZFMatrix< TVar >::AddFel | ( | TPZFMatrix< TVar > & | rhs, |
TPZVec< int64_t > & | destination | ||
) |
Performs a right hand side assemblage.
rhs | Load vector |
destination | Destine index on current matrix |
Definition at line 209 of file pzfmatrix.cpp.
Referenced by TPZFrontStructMatrix< front >::AssembleElement(), TPZGradientReconstruction::AssembleGlobalMatrix(), TPZFrontStructMatrix< front >::AssembleNew(), TPZStructMatrixST::ExecuteAssemble(), TPZStructMatrixGCTP::MultiThread_Assemble(), TPZStructMatrixGCTP::Serial_Assemble(), TPZStructMatrixGC::Serial_Assemble(), TPZStructMatrixCS::Serial_Assemble(), TPZStructMatrixOT::Serial_Assemble(), TPZStructMatrixOR::Serial_Assemble(), TPZPairStructMatrix::SerialAssemble(), TPZStructMatrixOR::ThreadData::ThreadAssembly(), TPZPairStructMatrix::ThreadData::ThreadAssembly1(), TPZStructMatrixGC::ThreadData::ThreadWork(), TPZStructMatrixCS::ThreadData::ThreadWork(), TPZStructMatrixOT::ThreadData::ThreadWork(), TPZStructMatrixGC::ThreadData::ThreadWorkResidual(), TPZStructMatrixOT::ThreadData::ThreadWorkResidual(), TPZPairStructMatrix::TPZPairStructMatrix(), and TPZFMatrix< STATE >::UpdateFrom().
void TPZFMatrix< TVar >::AddFel | ( | TPZFMatrix< TVar > & | rhs, |
TPZVec< int64_t > & | source, | ||
TPZVec< int64_t > & | destination | ||
) |
Performs a right hand side assemblage.
rhs | Load vector |
source | Source index on rhs |
destination | Destine index on current matrix |
Definition at line 226 of file pzfmatrix.cpp.
void TPZFMatrix< double >::AddFel | ( | TPZFMatrix< double > & | rhs, |
TPZVec< int64_t > & | source, | ||
TPZVec< int64_t > & | destination | ||
) |
Definition at line 243 of file pzfmatrix.cpp.
void TPZFMatrix< float >::AddFel | ( | TPZFMatrix< float > & | rhs, |
TPZVec< int64_t > & | source, | ||
TPZVec< int64_t > & | destination | ||
) |
Definition at line 261 of file pzfmatrix.cpp.
|
inline |
Definition at line 140 of file pzfmatrix.h.
|
overridevirtual |
Routines to send and receive messages.
Implements TPZSavable.
Reimplemented in TPZFNMatrix< N, TVar >, TPZFNMatrix< 12, REAL >, TPZFNMatrix< 180 >, TPZFNMatrix< TGeo::NNodes *3, REAL >, TPZFNMatrix< 660, REAL >, TPZFNMatrix< 9, STATE >, TPZFNMatrix< 9, T >, TPZFNMatrix< 9, REAL >, TPZFNMatrix< 50, REAL >, TPZFNMatrix< Geo::Dimension *Geo::NNodes >, TPZFNMatrix< 6, STATE >, TPZFNMatrix< 1000, STATE >, TPZFNMatrix< 220, REAL >, TPZFNMatrix< 3, T >, TPZFNMatrix< 9 >, TPZFNMatrix< 3, STATE >, and TPZFNMatrix< 30, std::complex< double > >.
Definition at line 2371 of file pzfmatrix.cpp.
Referenced by TPZFNMatrix< 30, std::complex< double > >::ClassId(), TPZFMatrix< STATE >::operator-=(), and TPZFNMatrix< 30, std::complex< double > >::operator=().
|
overrideprivatevirtual |
It clears data structure.
Reimplemented from TPZMatrix< TVar >.
Definition at line 2175 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< TVar >::operator const TVar *().
|
overridevirtual |
Compare the object for identity with the object pointed to, eventually copy the object.
compare both objects bitwise for identity. Put an entry in the log file if different overwrite the calling object if the override flag is true
Reimplemented from TPZMatrix< TVar >.
Definition at line 2236 of file pzfmatrix.cpp.
Referenced by TPZMaterialData::Compare(), and TPZFMatrix< STATE >::operator-=().
|
overridevirtual |
Compare the object for identity with the object pointed to, eventually copy the object.
compare both objects bitwise for identity. Put an entry in the log file if different generate an interupt if the override flag is true and the objects are different
compare both objects bitwise for identity. Put an entry in the log file if different overwrite the calling object if the override flag is true
Reimplemented from TPZMatrix< TVar >.
Definition at line 2277 of file pzfmatrix.cpp.
|
inline |
copy the values from a matrix with a different precision
Definition at line 150 of file pzfmatrix.h.
|
overridevirtual |
Cholesky Decomposition Optmized. for walks in the direction of the vector that composes the matrix.
Reimplemented from TPZMatrix< TVar >.
Definition at line 1574 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::Decompose_Cholesky(), TPZEqnArray< TVar >::main(), and TPZFMatrix< STATE >::operator-=().
|
overridevirtual |
Decomposes the current matrix using Cholesky method.
singular |
escopo
Reimplemented from TPZMatrix< TVar >.
Definition at line 1628 of file pzfmatrix.cpp.
|
overridevirtual |
Decomposes the current matrix using LDLt.
The current matrix has to be symmetric. "L" is lower triangular with 1.0 in its diagonal and "D" is a Diagonal matrix.
Reimplemented from TPZMatrix< TVar >.
Definition at line 1792 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::operator-=(), and TPZFMatrix< STATE >::Substitution().
|
overridevirtual |
LU Decomposition. Stores L and U matrices at the storage of the same matrix.
Reimplemented from TPZMatrix< TVar >.
Definition at line 1246 of file pzfmatrix.cpp.
Referenced by TPZSBFemElementGroup::CalcStiff(), TPZGradientReconstruction::TPZGradientData::ComputeGradient(), and TPZSpBlockDiagPivot< TVar >::Decompose_LU().
|
overridevirtual |
Reimplemented from TPZMatrix< TVar >.
Definition at line 1307 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::Decompose_LU(), TPZFMatrix< STATE >::operator-=(), and TPZFMatrix< STATE >::Transpose().
|
virtual |
LU Decomposition using pivot.
Definition at line 1164 of file pzfmatrix.cpp.
void TPZFMatrix< TVar >::DeterminantInverse | ( | TVar & | determinant, |
TPZFMatrix< TVar > & | inverse | ||
) |
Definition at line 507 of file pzfmatrix.cpp.
Referenced by pzgeom::GPr< TFather, Topology >::Jacobian(), and TPZFMatrix< STATE >::operator-=().
|
staticprivate |
Definition at line 2161 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::Decompose_Cholesky(), TPZFMatrix< STATE >::Decompose_LDLt(), TPZFMatrix< STATE >::Decompose_LU(), TPZFMatrix< STATE >::g(), TPZFMatrix< STATE >::GetVal(), TPZFMatrix< STATE >::MultAdd(), TPZFMatrix< TVar >::operator const TVar *(), TPZFMatrix< STATE >::operator()(), TPZFMatrix< STATE >::operator*(), TPZFMatrix< STATE >::operator+(), TPZFMatrix< STATE >::operator+=(), TPZFMatrix< STATE >::operator-(), TPZFMatrix< STATE >::operator-=(), TPZFMatrix< STATE >::operator=(), TPZFMatrix< STATE >::Redim(), TPZFMatrix< STATE >::Resize(), TPZFMatrix< STATE >::SetSize(), TPZFMatrix< STATE >::Substitution(), and TPZFMatrix< STATE >::TPZFMatrix().
|
inline |
Definition at line 594 of file pzfmatrix.h.
Referenced by TPZSBFemVolume::AdjustAxes3D(), TPZMatElasticity2D::ComputeSigma(), Dot(), TPZSYsmpMatrix< TVar >::MultAdd(), TPZSkylNSymMatrix< TVar >::MultAdd(), TPZFYsmpMatrix< TVar >::MultAdd(), TPZSpMatrix< TVar >::MultAdd(), TPZFMatrix< STATE >::MultAdd(), TPZSkylMatrix< TVar >::MultAdd(), TPZFYsmpMatrix< TVar >::MultAddMT(), TPZTransfer< TVar >::MultAddScalar(), TPZMatrix< STATE >::PrepareZ(), TPZFYsmpMatrix< TVar >::SolveSOR(), and TPZFMatrix< STATE >::UpdateFrom().
|
inlineoverridevirtual |
Get values without bounds checking
This method is faster than "Get" if DEBUG is defined.
Reimplemented from TPZMatrix< TVar >.
Definition at line 566 of file pzfmatrix.h.
Referenced by TPZFYsmpMatrix< TVar >::AddKel(), TPZMatrix< STATE >::AddSub(), TPZProjectEllipse::AlmostZeroToZero(), TPZDohrAssembly< STATE >::Assemble(), TPZDohrAssembly< STATE >::AssembleCoarse(), TPZAxesTools< TVar >::Axes2XYZ(), TPZNonLinMultGridAnalysis::CalcResidual(), TPZCondensedCompEl::CalcResidual(), TPZCondensedCompEl::CalcStiff(), TPZElasticityAxiMaterial::Contribute(), TPZDohrSubstruct< TVar >::Contribute_v3(), TPZElasticityAxiMaterial::ContributeInterface(), TPZInterpolationSpace::Convert2Axes(), TPZTransform< T >::CopyFrom(), TPZFMatrix< STATE >::Decompose_LDLt(), TPZFMatrix< STATE >::Decompose_LU(), TPZDohrAssembly< STATE >::Extract(), TPZDohrAssembly< STATE >::ExtractCoarse(), TPZSparseBlockDiagonal< TVar >::Gather(), TPZEquationFilter::Gather(), pzgeom::TPZGeoCube::GradX(), pzgeom::TPZGeoQuad::GradX(), pzgeom::TPZQuadraticPrism::GradX(), pzgeom::TPZGeoLinear::GradX(), pzgeom::TPZGeoPrism::GradX(), pzgeom::TPZCylinderMap< TGeo >::GradX(), pzgeom::TPZQuadraticPyramid::GradX(), pzgeom::TPZGeoTetrahedra::GradX(), pzgeom::TPZGeoPyramid::GradX(), pzgeom::TPZQuadraticCube::GradX(), pzgeom::TPZQuadraticTetra::GradX(), pzgeom::TPZGeoTriangle::GradX(), pzgeom::TPZWavyLine::GradX(), pzgeom::TPZEllipse3D::GradX(), pzgeom::TPZQuadraticTrig::GradX(), pzgeom::TPZQuadraticLine::GradX(), pzgeom::TPZQuadraticQuad::GradX(), TPZFMatrix< STATE >::GramSchmidt(), TEulerDiffusivity::InverseJacob(), pzgeom::TPZGeoTriangle::Jacobian(), pzgeom::GPr< TFather, Topology >::Jacobian(), TPZGeoEl::Jacobian(), TPZMathTools::JacobianConv(), ConvTest::JacobianConv(), TPZGeoElMapped< TBase >::JacobianConv(), TPZGeoEl::JacobianXYZ(), TPZGeoElMapped< TBase >::KsiBar(), TPZSubMeshAnalysis::LoadSolution(), TPZSubMeshFrontalAnalysis::LoadSolution(), TPZCompMesh::LoadSolution(), TPZDohrSubstruct< TVar >::LoadWeightedResidual(), TPZStencilMatrix< TVar >::MultAdd(), TPZBlockDiagonal< STATE >::MultAdd(), TPZSBMatrix< TVar >::MultAdd(), TPZTransfer< TVar >::MultAdd(), TPZFBMatrix< TVar >::MultAdd(), TPZVerySparseMatrix< TVar >::MultAdd(), TPZMatrix< STATE >::MultAdd(), TPZRefPatternTools::NodesHunter(), TPZDohrSubstructCondense< TTVar >::PermuteGather(), TPZDohrSubstructCondense< TTVar >::PermuteScatter(), TPZSBFemVolume::Print(), TPZProjectEllipse::PrintAxes(), TPZMatrix< STATE >::PutSub(), PYBIND11_MODULE(), TPZSparseBlockDiagonal< TVar >::Scatter(), TPZEquationFilter::Scatter(), TPZMatRed< TTVar, TPZFMatrix< TTVar > >::SetF(), TPZElasticityAxiMaterial::Solution(), TPZMatrix< STATE >::SolveGMRES(), TPZFYsmpMatrix< TVar >::SolveJacobi(), TPZStencilMatrix< TVar >::SolveSOR(), TPZMatrix< STATE >::SolveSOR(), TPZSFMatrix< TVar >::Subst_Backward(), TPZMatrix< STATE >::Subst_Backward(), TPZSSpMatrix< TVar >::Subst_Diag(), TPZSBMatrix< TVar >::Subst_Diag(), TPZSFMatrix< TVar >::Subst_Diag(), TPZMatrix< STATE >::Subst_Diag(), TPZSSpMatrix< TVar >::Subst_Forward(), TPZSFMatrix< TVar >::Subst_Forward(), TPZMatrix< STATE >::Subst_Forward(), TPZSBMatrix< TVar >::Subst_LBackward(), TPZSFMatrix< TVar >::Subst_LBackward(), TPZMatrix< STATE >::Subst_LBackward(), TPZSSpMatrix< TVar >::Subst_LForward(), TPZSBMatrix< TVar >::Subst_LForward(), TPZSFMatrix< TVar >::Subst_LForward(), TPZMatrix< STATE >::Subst_LForward(), TPZMatrix< STATE >::Substitution(), TPZGeoElMapped< TBase >::TKsiBar(), TPZBlockDiagonal< STATE >::TPZBlockDiagonal(), TPZTensor< STATE >::TPZTensor(), TPZVerySparseMatrix< TVar >::TPZVerySparseMatrix(), TPZTransfer< TVar >::TransferResidual(), TPZTransfer< TVar >::TransferSolution(), TPZMatRed< TTVar, TPZFMatrix< TTVar > >::UGlobal(), TPZFMatrix< STATE >::UpdateFrom(), TPZAxesTools< TVar >::VerifyAxes(), pzgeom::TPZGeoQuad::X(), pzgeom::TPZGeoCube::X(), pzgeom::TPZQuadraticPrism::X(), pzgeom::TPZGeoLinear::X(), pzgeom::TPZGeoPrism::X(), pzgeom::TPZGeoTriangle::X(), pzgeom::TPZGeoTetrahedra::X(), pzgeom::TPZQuadraticPyramid::X(), pzgeom::TPZGeoPyramid::X(), pzgeom::TPZQuadraticCube::X(), pzgeom::TPZQuadraticTetra::X(), pzgeom::TPZGeoPoint::X(), pzgeom::TPZQuadraticTrig::X(), pzgeom::GPr< TFather, Topology >::X(), pzgeom::TPZQuadraticLine::X(), pzgeom::TPZQuadraticQuad::X(), pzgeom::TPZCylinderMap< TGeo >::X(), pzgeom::TPZEllipse3D::X(), and pzgeom::TPZArc3D::X().
void TPZFMatrix< TVar >::GramSchmidt | ( | TPZFMatrix< TVar > & | Orthog, |
TPZFMatrix< TVar > & | TransfToOrthog | ||
) |
This method implements a Gram Schimidt method.
this = Orthog.TransfToOrthog.
Orthog | [out] each column represents a vector orthogonalized with respect to the first vector (first column of *this). Vectors are normalized |
TransfToOrthog | [out] is the basis change from *this to Orthog |
Making a copy of *this (Ortog = *this)
Definition at line 341 of file pzfmatrix.cpp.
Referenced by pzgeom::TPZArc3D::ComputeAtributes(), pzgeom::TPZGeoTriangle::Jacobian(), pzgeom::GPr< TFather, Topology >::Jacobian(), pzgeom::TPZGeoBlend< TGeo >::Jacobian(), NormalVector(), TPZFMatrix< STATE >::operator-(), and TPZFMatrix< STATE >::operator-=().
void TPZFMatrix< TPZFlopCounter >::GramSchmidt | ( | TPZFMatrix< TPZFlopCounter > & | Orthog, |
TPZFMatrix< TPZFlopCounter > & | TransfToOrthog | ||
) |
Definition at line 500 of file pzfmatrix.cpp.
|
inlineoverridevirtual |
Returns the approximate size of the memory footprint (amount of memory required to store this object).
Reimplemented from TPZMatrix< TVar >.
Definition at line 135 of file pzfmatrix.h.
|
overridevirtual |
It computes z = beta * y + alpha * opt(this)*x but z and x can not overlap in memory.
x | Is x on the above operation |
y | Is y on the above operation |
z | Is z on the above operation |
alpha | Is alpha on the above operation |
beta | Is beta on the above operation |
opt | Indicates if is Transpose or not |
Reimplemented from TPZMatrix< TVar >.
Definition at line 757 of file pzfmatrix.cpp.
Referenced by AdjustSolutionDerivatives(), TPZCompElPostProc< TCOMPEL >::CalcResidual(), TPZSBFemElementGroup::ComputeMassMatrix(), pzgeom::TPZEllipse3D::GradX(), TPZFMatrix< STATE >::MultAdd(), TPZTransfer< TVar >::MultAddScalar(), TPZFMatrix< STATE >::operator*(), TPZDohrSubstruct< TVar >::PrepareSystems(), Shape(), TPZFMatrix< STATE >::UpdateFrom(), and pzgeom::TPZEllipse3D::X().
|
static |
Definition at line 533 of file pzfmatrix.cpp.
|
inline |
Definition at line 422 of file pzfmatrix.h.
References TPZFMatrix< TVar >::Clear(), TPZFMatrix< TVar >::Error(), TPZFMatrix< TVar >::fElem, and TPZFMatrix< TVar >::PrintStatic().
|
inline |
Definition at line 577 of file pzfmatrix.h.
Referenced by TPZFMatrix< STATE >::AddFel(), TPZFMatrix< STATE >::Decompose_Cholesky(), TPZFMatrix< STATE >::Decompose_LU(), TPZFMatrix< STATE >::s(), and TPZFMatrix< STATE >::UpdateFrom().
|
inline |
Definition at line 605 of file pzfmatrix.h.
|
inline |
Definition at line 535 of file pzfmatrix.h.
Referenced by TPZFMatrix< STATE >::UpdateFrom().
TPZFMatrix< TVar > TPZFMatrix< TVar >::operator* | ( | const TVar | val | ) | const |
Definition at line 996 of file pzfmatrix.cpp.
TPZFMatrix< TVar > & TPZFMatrix< TVar >::operator*= | ( | const TVar | val | ) |
Definition at line 1006 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::GramSchmidt(), and TPZFMatrix< STATE >::operator-=().
TPZFMatrix< TVar > TPZFMatrix< TVar >::operator+ | ( | const TPZFMatrix< TVar > & | A | ) | const |
Definition at line 284 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::operator-(), and TPZFMatrix< STATE >::UpdateFrom().
TPZFMatrix< TVar > TPZFMatrix< TVar >::operator+ | ( | const TVar | val | ) | const |
Definition at line 975 of file pzfmatrix.cpp.
TPZFMatrix< TVar > & TPZFMatrix< TVar >::operator+= | ( | const TPZFMatrix< TVar > & | A | ) |
Definition at line 842 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::operator-=(), and TPZFMatrix< STATE >::UpdateFrom().
TPZFMatrix< TVar > & TPZFMatrix< TVar >::operator+= | ( | const TVar | val | ) |
Definition at line 962 of file pzfmatrix.cpp.
TPZFMatrix< TVar > TPZFMatrix< TVar >::operator- | ( | const TPZFMatrix< TVar > & | A | ) | const |
Definition at line 303 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::UpdateFrom().
TPZFMatrix< TVar > TPZFMatrix< TVar >::operator- | ( | const TVar | val | ) | const |
Definition at line 987 of file pzfmatrix.cpp.
TPZFMatrix< TVar > & TPZFMatrix< TVar >::operator-= | ( | const TPZFMatrix< TVar > & | A | ) |
Definition at line 856 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::UpdateFrom().
|
inline |
Definition at line 265 of file pzfmatrix.h.
|
virtual |
Generic operator with FULL matrices.
Reimplemented in TPZFNMatrix< N, TVar >, TPZFNMatrix< 12, REAL >, TPZFNMatrix< 180 >, TPZFNMatrix< TGeo::NNodes *3, REAL >, TPZFNMatrix< 660, REAL >, TPZFNMatrix< 9, STATE >, TPZFNMatrix< 9, T >, TPZFNMatrix< 9, REAL >, TPZFNMatrix< 50, REAL >, TPZFNMatrix< Geo::Dimension *Geo::NNodes >, TPZFNMatrix< 6, STATE >, TPZFNMatrix< 1000, STATE >, TPZFNMatrix< 220, REAL >, TPZFNMatrix< 3, T >, TPZFNMatrix< 9 >, TPZFNMatrix< 3, STATE >, and TPZFNMatrix< 30, std::complex< double > >.
Definition at line 131 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::Compare(), TPZFNMatrix< 30, std::complex< double > >::operator=(), TPZFNMatrix< 30, std::complex< double > >::TPZFNMatrix(), and TPZFMatrix< STATE >::UpdateFrom().
TPZFMatrix< TVar > & TPZFMatrix< TVar >::operator= | ( | const std::initializer_list< TVar > & | list | ) |
Definition at line 159 of file pzfmatrix.cpp.
TPZFMatrix< TVar > & TPZFMatrix< TVar >::operator= | ( | const std::initializer_list< std::initializer_list< TVar > > & | list | ) |
Definition at line 172 of file pzfmatrix.cpp.
TPZFMatrix< TVar > & TPZFMatrix< TVar >::operator= | ( | const TPZMatrix< TVar > & | A | ) |
Generic operator with matrices.
Definition at line 920 of file pzfmatrix.cpp.
TPZFMatrix< TVar > & TPZFMatrix< TVar >::operator= | ( | const TVar | val | ) |
Numeric operator with matrices.
Definition at line 948 of file pzfmatrix.cpp.
|
static |
Definition at line 2321 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::Compare(), and TPZFMatrix< TVar >::operator const TVar *().
|
inlineoverridevirtual |
Put values without bounds checking
This method is faster than "Put" if DEBUG is defined.
Reimplemented from TPZMatrix< TVar >.
Definition at line 548 of file pzfmatrix.h.
Referenced by TPZProjectEllipse::AlmostZeroToZero(), TPZElasticCriterion::ApplyStrainComputeSigma(), pztopology::TPZCube::ComputeDirections(), pztopology::TPZPrism::ComputeDirections(), pztopology::TPZPyramid::ComputeDirections(), pztopology::TPZTetrahedron::ComputeDirections(), pztopology::TPZQuadrilateral::ComputeDirections(), TPZElasticityAxiMaterial::Contribute(), TPZElasticityAxiMaterial::ContributeInterface(), TPZPorousElasticResponse::De_G_constant(), TPZPorousElasticResponse::De_Poisson_constant(), TPZFMatrix< STATE >::Decompose_LDLt(), TPZFMatrix< STATE >::Decompose_LU(), Hdiv2dPaper201504::ErrorH1(), hdivCurvedJCompAppMath::ErrorH1(), Hdiv3dPaper201504::ErrorH1(), Hdiv2dPaper201504::ErrorPrimalDual(), hdivCurvedJCompAppMath::ErrorPrimalDual(), Hdiv3dPaper201504::ErrorPrimalDual(), TPZMatrix< STATE >::GetSub(), TPZProjectEllipse::LeastSquaresToGetEllipse(), TPZProjectEllipse::LeastSquaresToGetSimpleEllipse(), TPZProjectEllipse::LeastSquaresToGetVerySimpleEllipse(), TPZSBMatrix< TVar >::MultAdd(), TPZFBMatrix< TVar >::MultAdd(), TPZVerySparseMatrix< TVar >::MultAdd(), TPZMatrix< STATE >::MultAdd(), TPZRefPatternTools::NodesHunter(), TPZMatrixMarket::Read(), pzshape::TPZShapeDisc::Shape0D(), TPZElasticityAxiMaterial::Solution(), TPZMatrix< STATE >::SolveEigensystemJacobi(), TPZSFMatrix< TVar >::Subst_Backward(), TPZMatrix< STATE >::Subst_Backward(), TPZSSpMatrix< TVar >::Subst_Diag(), TPZSBMatrix< TVar >::Subst_Diag(), TPZSFMatrix< TVar >::Subst_Diag(), TPZMatrix< STATE >::Subst_Diag(), TPZSSpMatrix< TVar >::Subst_Forward(), TPZSFMatrix< TVar >::Subst_Forward(), TPZMatrix< STATE >::Subst_Forward(), TPZSFMatrix< TVar >::Subst_LBackward(), TPZMatrix< STATE >::Subst_LBackward(), TPZSSpMatrix< TVar >::Subst_LForward(), TPZSBMatrix< TVar >::Subst_LForward(), TPZSFMatrix< TVar >::Subst_LForward(), TPZMatrix< STATE >::Subst_LForward(), TPZMatrix< STATE >::Substitution(), TPZFMatrix< STATE >::TPZFMatrix(), TPZMatRed< TTVar, TPZFMatrix< TTVar > >::UGlobal(), and TPZFMatrix< STATE >::UpdateFrom().
|
overridevirtual |
read objects from the stream
Reimplemented from TPZSavable.
Reimplemented in TPZFNMatrix< N, TVar >, TPZFNMatrix< 12, REAL >, TPZFNMatrix< 180 >, TPZFNMatrix< TGeo::NNodes *3, REAL >, TPZFNMatrix< 660, REAL >, TPZFNMatrix< 9, STATE >, TPZFNMatrix< 9, T >, TPZFNMatrix< 9, REAL >, TPZFNMatrix< 50, REAL >, TPZFNMatrix< Geo::Dimension *Geo::NNodes >, TPZFNMatrix< 6, STATE >, TPZFNMatrix< 1000, STATE >, TPZFNMatrix< 220, REAL >, TPZFNMatrix< 3, T >, TPZFNMatrix< 9 >, TPZFNMatrix< 3, STATE >, and TPZFNMatrix< 30, std::complex< double > >.
Definition at line 2205 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::Clear(), TPZFMatrix< STATE >::operator-=(), TPZFNMatrix< 30, std::complex< double > >::operator=(), TPZMaterialTest3D::Read(), TPZAnalysis::Read(), TPZGeoEl::Read(), TPZCompMesh::Read(), and TPZFNMatrix< 30, std::complex< double > >::Read().
|
inlineoverridevirtual |
Redimension a matrix and ZERO your elements.
Reimplemented from TPZMatrix< TVar >.
Definition at line 616 of file pzfmatrix.h.
Referenced by TPZInterpolatedElement::AdjustPreferredSideOrder(), TPZMulticamadaOrthotropic::AnalyticTensor(), TPZCompElHDivPressure< TSHAPE >::Append(), TPZCompElHDiv< TSHAPE >::Append(), Append(), TPZElementMatrix::ApplyConstraints(), TPZTransientAnalysis< TRANSIENTCLASS >::Assemble(), TPZSubMeshAnalysis::Assemble(), TPZDohrStructMatrix::Assemble(), TPZAnalysis::Assemble(), TPZAnalysis::AssembleResidual(), TPZAxesTools< TVar >::Axes2XYZ(), TPZEulerAnalysis::BufferLastStateAssemble(), TPZInterpolatedElement::CalcIntegral(), TPZNonLinMultGridAnalysis::CalcResidual(), TPZSBFemElementGroup::CalcStiff(), TPZAgglomerateElement::CalcStiff(), TPZCondensedCompEl::CalcStiff(), TPZSubCompMesh::CalcStiff(), pzshape::TPZShapeDisc::ChebyshevWithoutScale(), CheckConvergence(), TPZAnalysis::CleanUp(), TPZCompMesh::CleanUp(), TPZMatElastoPlastic< T, TMEM >::ComputeDeltaStrainVector(), TPZCompMeshTools::ComputeDifferenceNorm(), pztopology::TPZCube::ComputeDirections(), pztopology::TPZPrism::ComputeDirections(), pztopology::TPZPyramid::ComputeDirections(), pztopology::TPZLine::ComputeDirections(), pztopology::TPZTetrahedron::ComputeDirections(), pztopology::TPZPoint::ComputeDirections(), pztopology::TPZQuadrilateral::ComputeDirections(), TPZCompMesh::ComputeFillIn(), TPZTransientAnalysis< TRANSIENTCLASS >::ComputeFluxOnly(), TPZMaterialData::ComputeFluxValues(), TPZMaterialData::ComputeFunctionDivergence(), TPZGradientReconstruction::TPZGradientData::ComputeGradient(), TPZAxesTools< TVar >::ComputeGradX(), TPZSBFemVolume::ComputeKMatrices(), TPZTransientAnalysis< TRANSIENTCLASS >::ComputeLinearTangentMatrix(), TPZSBFemElementGroup::ComputeMassMatrix(), TPZTransientAnalysis< TRANSIENTCLASS >::ComputeMassMatrix(), TPZElasticity3D::ComputeStrainTensor(), TPZElasticity3D::ComputeStrainVector(), TPZElasticity3D::ComputeStressTensor(), TPZMatElastoPlastic< T, TMEM >::ComputeStressVector(), TPZElasticity3D::ComputeStressVector(), TPZNonLinearAnalysis::ComputeTangent(), TPZYCTrescaRegularized::ComputeTangent(), TPZElastoPlasticAnalysis::ComputeTangent(), TPZPlasticTest::ComputeTangent(), TPZYCWillamWarnke::ComputeTangent(), TPZYCVonMises::ComputeTangent(), TPZYCMohrCoulomb::ComputeTangent(), TPZLadeKimThermoForceA::ComputeTangent(), TPZYCTresca::ComputeTangent(), TPZYCDruckerPrager::ComputeTangent(), TPZLadeNelsonElasticResponse::ComputeTangent(), TPZYCLadeKim::ComputeTangent(), TPZLadeKim::ComputeTangent(), TPZPlasticStep< YC_t, TF_t, ER_t >::ComputeTangent(), TPZYCSandlerDimaggio::ComputeTangent(), TPZSubMeshFrontalAnalysis::CondensedSolution(), TPZViscoelastic::Contribute(), TPZDohrSubstruct< TVar >::Contribute_v3_local(), TPZSwelling::ContributeBC(), TPZSBandStructMatrix::CreateAssemble(), TPBSpStructMatrix::CreateAssemble(), TPZSpStructMatrix::CreateAssemble(), TPZBandStructMatrix::CreateAssemble(), TPZStructMatrixBase::CreateAssemble(), TPZSymetricSpStructMatrix::CreateAssemble(), TPZParSkylineStructMatrix::CreateAssemble(), TPZBlockDiagonalStructMatrix::CreateAssemble(), TPZStructMatrixST::CreateAssemble(), TPZStructMatrixGCTP::CreateAssemble(), TPZStructMatrixGC::CreateAssemble(), TPZStructMatrixCS::CreateAssemble(), TPZStructMatrixTBBFlow::CreateAssemble(), TPZStructMatrixOT::CreateAssemble(), TPZParFrontStructMatrix< front >::CreateAssemble(), TPZFMatrix< STATE >::DeterminantInverse(), TPZDiffusionConsLaw::Divergence(), TPZMixedElasticityMaterial::ElasticityModulusTensor(), TPZMatConvectionProblem::Errors(), TPZAnalysisError::EvaluateError(), TPZEulerAnalysis::EvaluateFluxEpsilon(), TPZCompMesh::ExpandSolution(), TPZFrontSym< TVar >::ExtractFrontMatrix(), TPZFrontNonSym< TVar >::ExtractFrontMatrix(), TPZBlockDiagonal< STATE >::GetBlock(), TPZSandlerExtended::GradF1SigmaTrial(), TPZSandlerExtended::GradF2SigmaTrial(), pzgeom::TPZGeoBlend< TGeo >::GradX(), TPZInterfaceElement::InitializeElementMatrix(), TPZReducedSpace::InitializeElementMatrix(), TPZMultiphysicsInterfaceElement::InitializeElementMatrix(), TPZInterpolationSpace::InitializeElementMatrix(), TPZCompElLagrange::InitializeElementMatrix(), TPZElementGroup::InitializeElementMatrix(), TPZMultiphysicsCompEl< TGeometry >::InitializeElementMatrix(), TPZInterfaceElement::InitMaterialData(), TPZSBFemVolume::InitMaterialData(), TPZInterpolationSpace::InitMaterialData(), TPZMultiphysicsInterfaceElement::InitMaterialData(), pzshape::TPZShapeDisc::IntegratedLegendre(), TPZMatrix< STATE >::Inverse(), TPZNonLinearAnalysis::IterativeProcess(), TPZElastoPlasticAnalysis::IterativeProcess(), TPZDiffusionConsLaw::JacobFlux(), TPZGeoElMapped< TBase >::KsiBar(), TPZProjectEllipse::LeastSquaresToGetEllipse(), TPZProjectEllipse::LeastSquaresToGetSimpleEllipse(), TPZProjectEllipse::LeastSquaresToGetVerySimpleEllipse(), pzshape::TPZShapeDisc::Legendre(), pzshape::TPZShapeDisc::LegendreWithoutScale(), TPZSkylNSymMatrix< TVar >::MultAdd(), TPZFYsmpMatrix< TVar >::MultAdd(), TPZDohrMatrix< TVar, TSubStruct >::MultAdd(), TPZFMatrix< STATE >::MultAdd(), TPZSkylMatrix< TVar >::MultAdd(), TPZMatrix< STATE >::Multiply(), TPZTransfer< TVar >::MultiplyScalar(), operator*(), TPZFMatrix< STATE >::operator*(), operator+(), TPZFMatrix< STATE >::operator+(), operator-(), TPZFMatrix< STATE >::operator-(), TPZFMatrix< STATE >::operator-=(), TPZDiffusionConsLaw::PointOperator(), pzshape::TPZShapeDisc::Polynomial(), pzshape::TPZShapeDisc::PolynomialWithoutScale(), TPZGradientReconstruction::ProjectionL2GradientReconstructed(), TPZCompMesh::ProjectSolution(), TPZGradientReconstruction::TPZGradientData::QRFactorization(), TPZMatrixMarket::Read(), TPZNonLinearAnalysis::Residual(), TPZYCTrescaRegularized::Residual(), TPZElastoPlasticAnalysis::Residual(), TPZPlasticTest::Residual(), TPZYCWillamWarnke::Residual(), TPZLadeKimThermoForceA::Residual(), TPZYCVonMises::Residual(), TPZYCMohrCoulomb::Residual(), TPZYCDruckerPrager::Residual(), TPZLadeNelsonElasticResponse::Residual(), TPZYCTresca::Residual(), TPZYCLadeKim::Residual(), TPZLadeKim::Residual(), TPZPlasticStep< YC_t, TF_t, ER_t >::Residual(), TPZYCSandlerDimaggio::Residual(), TPZTransientAnalysis< TRANSIENTCLASS >::RunTransient(), TPZSubMeshAnalysis::SetCompMesh(), TPZMatLaplacian::SetDimension(), TPZTransientAnalysis< TRANSIENTCLASS >::SetInitialSolutionAsZero(), TPZMat2dLin::SetMaterial(), TPZCompElHDivBound2< TSHAPE >::Shape(), TPZCompElHDiv< TSHAPE >::Shape(), TPZSBFemVolume::Shape(), pzshape::TPZShapeDisc::Shape0D(), pzshape::TPZShapeDisc::Shape2D(), pzshape::TPZShapeDisc::Shape2DFull(), pzshape::TPZShapeDisc::Shape3D(), TPZMixedPoisson::Solution(), TPZSequenceSolver< TVar >::Solve(), TPZMGAnalysis::Solve(), TPZMGSolver< TVar >::Solve(), TPZStepSolver< TVar >::Solve(), TPZEulerAnalysis::Solve(), TPZAnalysis::Solve(), TPZCheckConvergence(), TPZDiffusionConsLaw::TPZDiffusionConsLaw(), TPZElastoPlasticAnalysis::TPZElastoPlasticAnalysis(), TPZMatMixedPoisson3D::TPZMatMixedPoisson3D(), TPZMixedPoisson::TPZMixedPoisson(), TPZMultPlaca::TPZMultPlaca(), TPZPoroElastoPlasticAnalysis::TPZPoroElastoPlasticAnalysis(), TPZSubMeshAnalysis::TPZSubMeshAnalysis(), TPZSwelling::TPZSwelling(), TPZTransfer< TVar >::TransferResidual(), TPZPostProcAnalysis::TransferSolution(), TPZTransfer< TVar >::TransferSolution(), TPZDohrSubstructCondense< TTVar >::UGlobal(), TPZMatRed< TTVar, TPZFMatrix< TTVar > >::UGlobal(), TPZMatRed< TTVar, TPZFMatrix< TTVar > >::UGlobal2(), TPZCompMesh::UpdatePreviousState(), and TPZViscoelastic::UpdateQsi().
int TPZFMatrix< TVar >::Remodel | ( | const int64_t | newRows, |
const int64_t | wCols | ||
) |
Remodel the shape of the matrix, but keeping the same dimension.
Definition at line 1063 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::operator-=().
|
overridevirtual |
Redimension a matrix, but maintain your elements.
Reimplemented from TPZMatrix< TVar >.
Definition at line 1016 of file pzfmatrix.cpp.
Referenced by AdjustSolutionDerivatives(), TPZCompElDisc::AppendExternalShapeFunctions(), TPZCompMesh::AssembleError(), TPZAxesTools< TVar >::Axes2XYZ(), TPZSparseBlockDiagonal< TVar >::BuildFromMatrix(), TPZSBFemElementGroup::CalcStiff(), TPZDohrSubstruct< TVar >::ComputeCoarseStiffness(), pztopology::TPZCube::ComputeDirections(), pztopology::TPZPrism::ComputeDirections(), pztopology::TPZPyramid::ComputeDirections(), pztopology::TPZTetrahedron::ComputeDirections(), TPZMixedElasticityMaterial::ComputeDivergenceOnDeformed(), TPZEulerEquation::ComputeEulerFlux(), TPZReducedSpace::ComputeRequiredData(), TPZCompElDisc::ComputeShape(), TPZCompElHDivPressure< TSHAPE >::ComputeSolution(), TPZCompElDisc::ComputeSolution(), TPZGeoEl::ComputeXInverse(), TPZDohrSubstructCondense< TTVar >::Contribute_v1_local(), TPZDohrSubstruct< TVar >::Contribute_v1_local(), TPZDohrSubstructCondense< TTVar >::Contribute_v2_local(), TPZDohrSubstruct< TVar >::Contribute_v2_local(), TPZDohrSubstructCondense< TTVar >::ContributeDiagonalLocal(), TPZDohrSubstruct< TVar >::ContributeDiagonalLocal(), TPZBiharmonicEstimator::ContributeErrorsSimple(), TPZDohrSubstruct< TVar >::ContributeTestV1(), TPZInterpolationSpace::Convert2Axes(), TPZTransform< T >::CopyFrom(), TPZFMatrix< STATE >::CopyFrom(), TPZPlasticStepPV< YC_t, ER_t >::CopyFromFMatrixToTensor(), TPZTensor< STATE >::CopyToTensor(), TPZYCCamClayPV::D2DistanceToSurface(), TPZYCDruckerPragerPV::D2DistanceToSurfaceF1(), TPZSandlerExtended::D2DistFunc1(), TPZSandlerExtended::DDistFunc1(), TPZFMatrix< STATE >::DeterminantInverse(), TPZSandlerExtended::DF1Cart(), TPZSandlerExtended::DF2Cart(), TPZDohrAssembly< STATE >::Extract(), TPZDohrAssembly< STATE >::ExtractCoarse(), TPZMatRed< TTVar, TPZFMatrix< TTVar > >::F1Red(), TPZMixedElasticityMaterial::FillVecShapeIndex(), Hdiv2dPaper201504::ForcingH1(), TPZMatrix< STATE >::GetSub(), pzgeom::TPZGeoQuad::GradX(), pzgeom::TPZGeoCube::GradX(), pzgeom::TPZQuadraticPrism::GradX(), pzgeom::TPZGeoLinear::GradX(), pzgeom::TPZGeoPrism::GradX(), pzgeom::TPZGeoTetrahedra::GradX(), pzgeom::TPZQuadraticPyramid::GradX(), pzgeom::TPZGeoPyramid::GradX(), pzgeom::TPZQuadraticCube::GradX(), pzgeom::TPZQuadraticTetra::GradX(), pzgeom::TPZGeoPoint::GradX(), pzgeom::TPZGeoTriangle::GradX(), pzgeom::TPZWavyLine::GradX(), pzgeom::TPZQuadraticTrig::GradX(), pzgeom::TPZQuadraticLine::GradX(), pzgeom::TPZQuadraticQuad::GradX(), TPZGeoElRefLess< TGeo >::GradX(), TPZFMatrix< STATE >::GramSchmidt(), Hdiv3dPaper201504::Hdiv3dPaper201504(), TPZCompElHDiv< TSHAPE >::InitMaterialData(), pzgeom::TPZGeoTriangle::Jacobian(), pzgeom::GPr< TFather, Topology >::Jacobian(), pzgeom::TPZGeoBlend< TGeo >::Jacobian(), TPZGeoElSide::Jacobian(), TPZGeoEl::Jacobian(), TPZSandlerExtended::JacobianCoVertex(), TPZSandlerExtended::Jacobianf1(), TPZSandlerExtended::Jacobianf2(), TPZSandlerExtended::Jacobianf2CoVertex(), TPZGeoEl::JacobianXYZ(), TPZMultiphase::K(), TPZMultiphase::LoadKMap(), TPZCompMesh::LoadSolution(), TPZMatPoisson3d::LocalNeumanContribute(), pzgeom::TPZGeoBlend< TGeo >::MapToNeighSide(), pztopology::TPZPyramid::MapToSide(), pztopology::TPZCube::MapToSide(), pztopology::TPZPrism::MapToSide(), pztopology::TPZTetrahedron::MapToSide(), pztopology::TPZPoint::MapToSide(), pztopology::TPZTriangle::MapToSide(), pztopology::TPZQuadrilateral::MapToSide(), TPZSparseBlockDiagonal< TVar >::MultAdd(), TPZYCSandlerDimaggio::NewtonF1(), TPZGeoEl::NodesCoordinates(), TPZFMatrix< STATE >::operator-=(), TPZFMatrix< STATE >::operator=(), TPZBiharmonicEstimator::OrderSolution(), TPZBiharmonicEstimator::OrderSolutionLeft(), TPZBiharmonicEstimator::OrderSolutionRight(), TPZMohrCoulombNeto::PieceWise(), TPZMohrCoulombPV::PieceWise(), TPZAnalysis::PostProcessErrorSerial(), TPZInterpolationSpace::ProjectFlux(), TPZFMatrix< STATE >::Read(), TPZElementMatrix::Reset(), TPZTensor< T >::TPZDecomposed::ResidualCheckConv(), RotationMatrix(), Hdiv2dPaper201504::Run(), hdivCurvedJCompAppMath::Run(), Hdiv3dPaper201504::Run(), TPZAnalysis::SetCompMesh(), TPZCompMesh::SetElementSolution(), TPZMatMixedPoisson3D::SetPermeability(), TPZSBFemVolume::SetPhiEigVal(), pzshape::SPr< TFather >::Shape(), TPZCompElHDivBound2< TSHAPE >::Shape(), pzshape::TPZShapeLinear::Shape(), pzshape::SPr< TFather >::ShapeCorner(), pzshape::TPZShapeTriang::ShapeInternal(), pzshape::TPZShapeQuad::ShapeInternal(), pzshape::TPZShapeCube::ShapeInternal(), TPZReducedSpace::ShapeX(), SolExata(), Hdiv2dPaper201504::SolExata(), hdivCurvedJCompAppMath::SolExata(), Hdiv3dPaper201504::SolExata(), Hdiv2dPaper201504::SolExataH1(), hdivCurvedJCompAppMath::SolExataH1(), Hdiv3dPaper201504::SolExataH1(), TPZMatrix< STATE >::SolveEigensystemJacobi(), TPZDohrSubstruct< TVar >::SolveSystemZi(), TPZFMatrix< STATE >::Substitution(), TPZTensor< T >::TPZDecomposed::TangentCheckConv(), TPZMatMixedPoisson3D::TPZMatMixedPoisson3D(), TPZMixedDarcyFlow::TPZMixedDarcyFlow(), TPZMixedPoisson::TPZMixedPoisson(), TPZMultPlaca::TPZMultPlaca(), TPZProjectEllipse::TPZProjectEllipse(), pzgeom::TPZGeoBlend< TGeo >::X(), TPZAxesTools< TVar >::XYZ2Axes(), TPZAnalyticSolution::~TPZAnalyticSolution(), and TPZTransform< T >::~TPZTransform().
|
inlineoverridevirtual |
The operators check on the bounds if the DEBUG variable is defined.
row | Row number. |
col | Column number. |
Reimplemented from TPZMatrix< TVar >.
Definition at line 588 of file pzfmatrix.h.
Referenced by TPZFBMatrix< TVar >::Transpose(), and TPZFMatrix< STATE >::UpdateFrom().
int TPZFMatrix< TVar >::SetSize | ( | int64_t | newRows, |
int64_t | newCols | ||
) |
Redimension the matrix doing nothing with the elements.
Definition at line 2376 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::operator-=().
|
static |
Definition at line 1323 of file pzfmatrix.cpp.
Referenced by TPZSBFemElementGroup::CalcStiff(), TPZGradientReconstruction::TPZGradientData::ComputeGradient(), TPZFMatrix< STATE >::operator-=(), TPZBlockDiagonal< STATE >::Substitution(), TPZFMatrix< STATE >::Substitution(), and TPZSpBlockDiagPivot< TVar >::Substitution2().
|
overridevirtual |
Computes Forward and Backward substitution for a "LU" decomposed matrix.
B | right hand side and result after all |
Reimplemented from TPZMatrix< TVar >.
Definition at line 1364 of file pzfmatrix.cpp.
|
virtual |
LU substitution using pivot.
Definition at line 1496 of file pzfmatrix.cpp.
|
static |
LU substitution using pivot. Static version.
Definition at line 1676 of file pzfmatrix.cpp.
void TPZFMatrix< TVar >::TimesBetaPlusZ | ( | const TVar | beta, |
const TPZFMatrix< TVar > & | z | ||
) |
Performs an operation *this = this * beta + z.
beta | Being beta on above opereation |
z | Being z on above operation |
Definition at line 898 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::UpdateFrom(), and TPZFMatrix< STATE >::ZAXPY().
|
overridevirtual |
It makes *T the transpose of current matrix.
Reimplemented from TPZMatrix< TVar >.
Definition at line 1073 of file pzfmatrix.cpp.
Referenced by TPZSBFemVolume::AdjustAxes3D(), AdjustSolutionDerivatives(), TPZAxesTools< TVar >::Axes2XYZ(), TPZMixedElasticityMaterial::ComputeDivergenceOnDeformed(), TPZGradientReconstruction::TPZGradientData::ComputeGradient(), pztopology::TPZPrism::ComputeHCurlDirections(), pztopology::TPZCube::ComputeHCurlDirections(), pztopology::TPZTetrahedron::ComputeHCurlDirections(), TPZSBFemVolume::ComputeKMatrices(), TPZGeoEl::ComputeXInverse(), TPZMatElastoPlastic2D< T, TMEM >::Contribute(), TPZMatElastoPlastic< T, TMEM >::Contribute(), TPZMatPorous< T, TMEM >::Contribute(), TPZPlaca::Errors(), TPZGeoElSide::GradX(), TPZFMatrix< STATE >::GramSchmidt(), pzgeom::TPZGeoTriangle::Jacobian(), pzgeom::GPr< TFather, Topology >::Jacobian(), pzgeom::TPZGeoBlend< TGeo >::Jacobian(), TPZGeoElSide::Jacobian(), TPZProjectEllipse::LeastSquaresToGetEllipse(), TPZProjectEllipse::LeastSquaresToGetSimpleEllipse(), TPZProjectEllipse::LeastSquaresToGetVerySimpleEllipse(), TPZDiffusionConsLaw::LS(), TPZDiffusionConsLaw::PointOperator(), TPZDohrSubstruct< TVar >::PrepareSystems(), TPZGradientReconstruction::TPZGradientData::QRFactorization(), TPZMatPlaca2::SetNAxes(), TPZMatOrthotropic::Solution(), TPZMatPlaca2::TPZMatPlaca2(), TPZMultPlaca::TPZMultPlaca(), TPZPlaca::TPZPlaca(), and TPZAxesTools< TVar >::XYZ2Axes().
void TPZFMatrix< TVar >::Transpose | ( | ) |
Definition at line 1086 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::operator-=(), and TPZFMatrix< STATE >::Transpose().
|
inlineoverridevirtual |
Updates the values of the matrix based on the values of the matrix.
Reimplemented from TPZMatrix< TVar >.
Definition at line 161 of file pzfmatrix.h.
|
overridevirtual |
Writes this object to the TPZStream buffer. Include the classid if withclassid = true.
Writes this object to the TPZStream buffer. Include the classid if withclassid = true
Reimplemented from TPZSavable.
Reimplemented in TPZFNMatrix< N, TVar >, TPZFNMatrix< 12, REAL >, TPZFNMatrix< 180 >, TPZFNMatrix< TGeo::NNodes *3, REAL >, TPZFNMatrix< 660, REAL >, TPZFNMatrix< 9, STATE >, TPZFNMatrix< 9, T >, TPZFNMatrix< 9, REAL >, TPZFNMatrix< 50, REAL >, TPZFNMatrix< Geo::Dimension *Geo::NNodes >, TPZFNMatrix< 6, STATE >, TPZFNMatrix< 1000, STATE >, TPZFNMatrix< 220, REAL >, TPZFNMatrix< 3, T >, TPZFNMatrix< 9 >, TPZFNMatrix< 3, STATE >, and TPZFNMatrix< 30, std::complex< double > >.
Definition at line 2225 of file pzfmatrix.cpp.
Referenced by TPZFMatrix< STATE >::Clear(), main(), TPZFMatrix< STATE >::operator-=(), TPZFNMatrix< 30, std::complex< double > >::operator=(), TPZMaterialTest3D::Write(), TPZAnalysis::Write(), TPZGeoEl::Write(), TPZCompMesh::Write(), and TPZFNMatrix< 30, std::complex< double > >::Write().
void TPZFMatrix< TVar >::ZAXPY | ( | const TVar | alpha, |
const TPZFMatrix< TVar > & | p | ||
) |
Performs an ZAXPY operation being *this += alpha * p.
alpha | Being alpha on above opereation |
p | Being p on above operation |
Definition at line 879 of file pzfmatrix.cpp.
Referenced by TPZSandlerDimaggio< SANDLERDIMAGGIOPARENT >::ComputeDep(), TPZElastoPlasticAnalysis::ManageIterativeProcess(), TPZFMatrix< STATE >::operator-=(), and TPZFMatrix< STATE >::UpdateFrom().
|
inlineoverridevirtual |
Makes Zero all the elements.
Reimplemented from TPZMatrix< TVar >.
Definition at line 651 of file pzfmatrix.h.
Referenced by TPZPoroElastoPlasticAnalysis::AcceptSolution(), TPZElastoPlasticAnalysis::AcceptSolution(), TPZCompElDisc::AppendExternalShapeFunctions(), TPZCompMesh::AssembleError(), TPZAxesTools< TVar >::Axes2XYZ(), TPZEulerAnalysis::BufferLastStateAssemble(), TPZInterpolationSpace::BuildTransferMatrix(), TPZInterpolatedElement::BuildTransferMatrix(), TPZCompElDisc::BuildTransferMatrix(), TPZSBFemElementGroup::CalcStiff(), pzshape::TPZShapeLinear::Chebyshev(), TPZPlasticDiagnostic::CheckGlobal(), TPZCompElDisc::ComputeShape(), TPZMulticamadaOrthotropic::ComputeSolution(), TPZInterpolatedElement::ComputeSolution(), TPZInterfaceElement::ComputeSolution(), TPZCompElDisc::ComputeSolution(), TPZCompElHDiv< TSHAPE >::ComputeSolutionHDiv(), TPZDohrPrecond< TVar, TSubStruct >::ComputeV1(), TPZBCTension::Contribute(), TPZMixedPoissonParabolic::Contribute(), TPZMixedPoisson::Contribute(), TPZMatMixedPoisson3D::Contribute(), TPZDohrSubstruct< TVar >::Contribute_v2_local(), TPZDohrSubstruct< TVar >::Contribute_v3(), TPZDohrSubstruct< TVar >::ContributeKU(), TPZDohrSubstruct< TVar >::ContributeKULocal(), TPZMatMixedPoisson3D::ContributeWithoutSecondIntegration(), TPZElasticResponse::De(), TPZDiffusionConsLaw::Divergence(), TPZMGAnalysis::ElementError(), TPZEulerAnalysis::EvaluateFluxEpsilon(), TPZDummyFunction< TVar >::Execute(), TPZEquationFilter::Gather(), TPZMixedPoisson::GetPermeability(), pzgeom::TPZGeoQuad::GradX(), pzgeom::TPZGeoCube::GradX(), pzgeom::TPZQuadraticPrism::GradX(), pzgeom::TPZGeoLinear::GradX(), pzgeom::TPZGeoPrism::GradX(), pzgeom::TPZGeoTetrahedra::GradX(), pzgeom::TPZQuadraticPyramid::GradX(), pzgeom::TPZGeoPyramid::GradX(), pzgeom::TPZQuadraticCube::GradX(), pzgeom::TPZQuadraticTetra::GradX(), pzgeom::TPZGeoTriangle::GradX(), pzgeom::TPZWavyLine::GradX(), pzgeom::TPZQuadraticTrig::GradX(), pzgeom::TPZQuadraticLine::GradX(), pzgeom::TPZQuadraticQuad::GradX(), TPZFMatrix< STATE >::GramSchmidt(), Input::InsertElasticityCubo(), InsertElasticityCubo(), InsertViscoElasticityCubo(), TPZMatrix< STATE >::Inverse(), TPZElastoPlasticAnalysis::IterativeProcess(), TPZElastoPlasticAnalysis::IterativeProcessPrecomputedMatrix(), pzgeom::TPZGeoTriangle::Jacobian(), pzgeom::TPZGeoBlend< TGeo >::Jacobian(), TPZGeoElSide::Jacobian(), TPZGeoEl::Jacobian(), TPZGeoEl::JacobianXYZ(), TPZMultiphase::K(), TPZProjectEllipse::LeastSquaresToGetEllipse(), TPZProjectEllipse::LeastSquaresToGetSimpleEllipse(), TPZProjectEllipse::LeastSquaresToGetVerySimpleEllipse(), TPZSBFemElementGroup::LoadEigenVector(), TCedricTest::LoadInterpolation(), TPZMultiphase::LoadKMap(), TPZAnalysis::LoadShape(), TPZSubCompMesh::main(), TEulerDiffusivity::MatrixDiff(), TPZSparseBlockDiagonal< TVar >::MultAdd(), TPZSYsmpMatrix< TVar >::MultAdd(), TPZFYsmpMatrix< TVar >::MultAdd(), TPZFMatrix< STATE >::MultAdd(), TPZGeoEl::NodesCoordinates(), NullForce(), TPZNonLinMultGridAnalysis::OneGridAlgorithm(), TPZFMatrix< STATE >::operator-=(), TPZDiffusionConsLaw::PointOperator(), TPZDohrSubstruct< TVar >::PrepareSystems(), TPZCompMesh::ProjectSolution(), TPZFMatrix< STATE >::Redim(), TPZInterpolatedElement::RestrainSide(), TPZEulerAnalysis::Run(), TPZEquationFilter::Scatter(), TPZTransientAnalysis< TRANSIENTCLASS >::SetInitialSolutionAsZero(), TPZMatLaplacian::SetParameters(), TPZMixedPoisson::SetPermeability(), TPZMixedDarcyFlow::SetPermeability(), TPZMatLaplacian::SetPermeability(), pzshape::TPZShapePrism::Shape(), pzshape::TPZShapePiram::Shape(), pzshape::TPZShapeTetra::Shape(), pzshape::TPZShapeCube::Shape(), TPZSBFemVolume::Shape(), pzshape::TPZShapeDisc::Shape2D(), pzshape::TPZShapeDisc::Shape2DFull(), pzshape::TPZShapeDisc::Shape3D(), TPZBuildMultiphysicsMesh::ShowShape(), TPZAnalysis::ShowShape(), TPZNonLinMultGridAnalysis::SmoothingSolution(), TPZNonLinMultGridAnalysis::SmoothingSolution2(), TPZSequenceSolver< TVar >::Solve(), TPZMatrix< STATE >::SolveEigensystemJacobi(), TPZMatrix< STATE >::SolveJacobi(), TPZStencilMatrix< TVar >::SolveSOR(), TPZFYsmpMatrix< TVar >::SolveSOR(), TPZMatrix< STATE >::SolveSOR(), TPZSkylMatrix< TVar >::SolveSOR(), TPZDohrSubstruct< TVar >::SolveSystemPhi(), TPZDohrSubstruct< TVar >::SolveSystemZi(), TPZSpBlockDiagPivot< TVar >::Substitution(), TPZSparseBlockDiagonal< TVar >::Substitution(), TPZElastoPlasticAnalysis::TPZElastoPlasticAnalysis(), TPZNonLinearAnalysis::TPZNonLinearAnalysis(), TPZTransform< T >::TPZTransform(), pztopology::TPZPyramid::TransformElementToSide(), pztopology::TPZCube::TransformElementToSide(), pztopology::TPZPrism::TransformElementToSide(), pztopology::TPZLine::TransformElementToSide(), pztopology::TPZTetrahedron::TransformElementToSide(), pztopology::TPZTriangle::TransformElementToSide(), pztopology::TPZQuadrilateral::TransformElementToSide(), pztopology::TPZPyramid::TransformSideToElement(), pztopology::TPZCube::TransformSideToElement(), pztopology::TPZLine::TransformSideToElement(), pztopology::TPZPrism::TransformSideToElement(), pztopology::TPZTetrahedron::TransformSideToElement(), pztopology::TPZTriangle::TransformSideToElement(), pztopology::TPZQuadrilateral::TransformSideToElement(), TPZNonLinMultGridAnalysis::TwoGridAlgorithm(), TPZEulerAnalysis::UpdateHistory(), TEulerDiffusivity::ValJacobFlux(), and TPZAxesTools< TVar >::XYZ2Axes().
|
friend |
Definition at line 146 of file pzfmatrix.h.
|
friend |
Definition at line 145 of file pzfmatrix.h.
|
private |
Definition at line 431 of file pzfmatrix.h.
Referenced by TPZFMatrix< STATE >::Adress(), TPZFMatrix< STATE >::Clear(), TPZFMatrix< STATE >::Compare(), TPZFMatrix< STATE >::CopyFrom(), TPZFMatrix< STATE >::Decompose_Cholesky(), TPZFMatrix< STATE >::Decompose_LDLt(), TPZFMatrix< STATE >::g(), TPZFMatrix< STATE >::GetVal(), TPZFMatrix< STATE >::MultAdd(), TPZFMatrix< TVar >::operator const TVar *(), TPZFMatrix< STATE >::operator()(), TPZFMatrix< STATE >::operator*=(), TPZFMatrix< STATE >::operator+(), TPZFMatrix< STATE >::operator+=(), TPZFMatrix< STATE >::operator-(), TPZFMatrix< STATE >::operator-=(), TPZFMatrix< STATE >::operator=(), TPZFMatrix< STATE >::PutVal(), TPZFMatrix< STATE >::Read(), TPZFMatrix< STATE >::Redim(), TPZFMatrix< STATE >::Resize(), TPZFMatrix< STATE >::SetSize(), TPZFMatrix< STATE >::Substitution(), TPZFMatrix< STATE >::TimesBetaPlusZ(), TPZFMatrix< STATE >::TPZFMatrix(), TPZFMatrix< STATE >::Transpose(), TPZFMatrix< STATE >::Write(), TPZFMatrix< STATE >::ZAXPY(), TPZFMatrix< STATE >::Zero(), and TPZFMatrix< STATE >::~TPZFMatrix().
|
private |
Definition at line 432 of file pzfmatrix.h.
Referenced by TPZFMatrix< STATE >::Clear(), TPZFMatrix< STATE >::operator=(), TPZFMatrix< STATE >::Redim(), TPZFMatrix< STATE >::Resize(), TPZFMatrix< STATE >::SetSize(), and TPZFMatrix< STATE >::~TPZFMatrix().
|
private |
Definition at line 433 of file pzfmatrix.h.
Referenced by TPZFMatrix< STATE >::operator=(), TPZFMatrix< STATE >::Redim(), TPZFMatrix< STATE >::Resize(), TPZFMatrix< STATE >::SetSize(), and TPZFMatrix< STATE >::~TPZFMatrix().