NeoPZ
|
Implements a matrix whose nonzero elements are stored in binary tree. Matrix. More...
#include <pzfmatrix.h>
Public Member Functions | |
TPZVerySparseMatrix () | |
TPZVerySparseMatrix (int64_t rows, int64_t cols) | |
TPZVerySparseMatrix (int64_t rows, int64_t cols, TVar val) | |
virtual | ~TPZVerySparseMatrix () |
TPZVerySparseMatrix (const TPZVerySparseMatrix< TVar > ©) | |
TPZVerySparseMatrix (const TPZFMatrix< TVar > &cp) | |
void | Simetrize () override |
Simetrizes copies the data of the matrix to make its data simetric. More... | |
int | PutVal (const int64_t row, const int64_t col, const TVar &val) override |
Put values checking bounds. More... | |
virtual const TVar & | GetVal (const int64_t row, const int64_t col) const override |
Get values checking bounds. 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... | |
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... | |
virtual void | Transpose (TPZVerySparseMatrix< TVar > *T) const |
It makes *T the transpose of current matrix. More... | |
virtual void | Transpose (TPZMatrix< TVar > *const T) const override |
It makes *T the transpose of current matrix. More... | |
int | ClassId () const override |
Saveable methods. More... | |
void | Write (TPZStream &buf, int withclassid) const override |
Packs the object structure in a stream of bytes. More... | |
void | Read (TPZStream &buf, void *context) override |
Unpacks the object structure from a stream of bytes. More... | |
std::map< std::pair< int64_t, int64_t >, TVar >::const_iterator | MapBegin () const |
std::map< std::pair< int64_t, int64_t >, TVar >::const_iterator | MapEnd () const |
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 |
virtual int64_t | MemoryFootprint () const |
Returns the approximate size of the memory footprint (amount of memory required to store this object). More... | |
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 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... | |
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 void | UpdateFrom (TPZAutoPointer< TPZMatrix< TVar > >) |
Updates the values of the matrix based on the values of the matrix. More... | |
virtual int | IsSimetric () const |
Checks if the current matrix is symmetric. More... | |
int | IsSquare () const |
Checks if current matrix is square. 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_LU (std::list< int64_t > &singular) |
Decomposes the current matrix using LU decomposition. More... | |
virtual int | Decompose_LU () |
virtual int | Decompose_Cholesky () |
Decomposes the current matrix using Cholesky method. The current matrix has to be symmetric. More... | |
virtual int | Decompose_Cholesky (std::list< int64_t > &singular) |
Decomposes the current matrix using Cholesky method. 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 | Decompose_LDLt () |
Decomposes the current matrix using LDLt. More... | |
virtual int | Substitution (TPZFMatrix< TVar > *B) const |
Computes Forward and Backward substitution for a "LU" decomposed 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... | |
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 int | Resize (const int64_t newRows, const int64_t newCols) |
Redimensions a matriz keeping the previous values. More... | |
virtual int | Redim (const int64_t newRows, const int64_t newCols) |
Redimensions the matrix reinitializing it with zero. More... | |
virtual int | Zero () |
Zeroes the matrix. 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 | |
Protected Attributes | |
std::map< std::pair< int64_t, int64_t >, TVar > | fExtraSparseData |
Save elements different from zero, of Sparse matrix. 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... | |
Private Member Functions | |
void | WriteMap (TPZStream &buf, int withclassid, const std::map< std::pair< int64_t, int64_t >, TVar > &TheMap) const |
Auxiliary functions only reading and writing a map as the third paremeter. More... | |
void | ReadMap (TPZStream &buf, void *context, std::map< std::pair< int64_t, int64_t >, TVar > &TheMap) |
Friends | |
class | TPZFYsmpMatrix< TVar > |
Additional Inherited Members | |
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) |
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... | |
virtual int | Clear () |
It clears data structure. 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... | |
Static Protected Attributes inherited from TPZMatrix< TVar > | |
static TVar | gZero |
Initializing value to static variable. More... | |
Implements a matrix whose nonzero elements are stored in binary tree. Matrix.
Definition at line 33 of file pzfmatrix.h.
TPZVerySparseMatrix< TVar >::TPZVerySparseMatrix | ( | ) |
Definition at line 11 of file tpzverysparsematrix.cpp.
Referenced by TPZVerySparseMatrix< TVar >::TPZVerySparseMatrix().
|
inline |
Definition at line 36 of file tpzverysparsematrix.h.
|
inline |
Definition at line 41 of file tpzverysparsematrix.h.
References TPZVerySparseMatrix< TVar >::~TPZVerySparseMatrix().
|
virtual |
Definition at line 18 of file tpzverysparsematrix.cpp.
Referenced by TPZVerySparseMatrix< TVar >::TPZVerySparseMatrix().
|
inline |
Definition at line 49 of file tpzverysparsematrix.h.
References TPZVerySparseMatrix< TVar >::GetVal(), TPZVerySparseMatrix< TVar >::PutVal(), TPZVerySparseMatrix< TVar >::Simetrize(), TPZVerySparseMatrix< TVar >::TPZVerySparseMatrix(), and val().
TPZVerySparseMatrix< TVar >::TPZVerySparseMatrix | ( | const TPZFMatrix< TVar > & | cp | ) |
Definition at line 49 of file tpzverysparsematrix.cpp.
References TPZMatrix< TVar >::fCol, TPZMatrix< TVar >::fRow, TPZFMatrix< TVar >::GetVal(), IsZero(), and TPZVerySparseMatrix< TVar >::PutVal().
|
overridevirtual |
Saveable methods.
Reimplemented from TPZMatrix< TVar >.
Definition at line 127 of file tpzverysparsematrix.h.
References TPZMatrix< TVar >::ClassId(), and Hash().
Referenced by TPZVerySparseMatrix< TVar >::Transpose().
|
overridevirtual |
Get values checking bounds.
Reimplemented from TPZMatrix< TVar >.
Definition at line 117 of file tpzverysparsematrix.cpp.
References TPZMatrix< TVar >::fCol, TPZVerySparseMatrix< TVar >::fExtraSparseData, TPZMatrix< TVar >::fRow, and TPZMatrix< TVar >::gZero.
Referenced by TPZVerySparseMatrix< TVar >::TPZVerySparseMatrix().
|
inline |
Definition at line 111 of file tpzverysparsematrix.h.
References TPZVerySparseMatrix< TVar >::fExtraSparseData.
Referenced by TPZFMatrix< STATE >::TPZFMatrix().
|
inline |
Definition at line 112 of file tpzverysparsematrix.h.
References TPZVerySparseMatrix< TVar >::fExtraSparseData, TPZVerySparseMatrix< TVar >::ReadMap(), and TPZVerySparseMatrix< TVar >::WriteMap().
Referenced by TPZFMatrix< STATE >::TPZFMatrix().
|
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 138 of file tpzverysparsematrix.cpp.
References TPZMatrix< TVar >::Cols(), TPZVerySparseMatrix< TVar >::fExtraSparseData, TPZFMatrix< TVar >::GetVal(), TPZMatrix< TVar >::PrepareZ(), TPZFMatrix< TVar >::PutVal(), TPZMatrix< TVar >::Rows(), and val().
Referenced by TPZVerySparseMatrix< TVar >::s().
|
overridevirtual |
Put values checking bounds.
Reimplemented from TPZMatrix< TVar >.
Definition at line 23 of file tpzverysparsematrix.cpp.
References TPZMatrix< TVar >::fCol, TPZVerySparseMatrix< TVar >::fExtraSparseData, TPZMatrix< TVar >::fRow, and val().
Referenced by TPZVerySparseMatrix< TVar >::TPZVerySparseMatrix().
|
overridevirtual |
Unpacks the object structure from a stream of bytes.
buf | The buffer containing the object in a packed form |
context |
Reimplemented from TPZMatrix< TVar >.
Definition at line 226 of file tpzverysparsematrix.cpp.
References TPZMatrix< TVar >::fCol, TPZMatrix< TVar >::fDecomposed, TPZMatrix< TVar >::fDefPositive, TPZVerySparseMatrix< TVar >::fExtraSparseData, TPZMatrix< TVar >::fRow, TPZMatrix< TVar >::gZero, TPZStream::Read(), and TPZVerySparseMatrix< TVar >::ReadMap().
Referenced by TPZVerySparseMatrix< TVar >::Transpose().
|
private |
Definition at line 237 of file tpzverysparsematrix.cpp.
References TPZStream::Read().
Referenced by TPZVerySparseMatrix< TVar >::MapEnd(), and TPZVerySparseMatrix< TVar >::Read().
|
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 71 of file tpzverysparsematrix.h.
References CLONEDEF, TPZMatrix< TVar >::Cols(), DebugStop, TPZMatrix< TVar >::Error(), TPZVerySparseMatrix< TVar >::fExtraSparseData, TPZVerySparseMatrix< TVar >::MultAdd(), TPZMatrix< TVar >::Rows(), and TPZVerySparseMatrix< TVar >::Transpose().
|
overridevirtual |
Simetrizes copies the data of the matrix to make its data simetric.
Reimplemented from TPZMatrix< TVar >.
Definition at line 62 of file tpzverysparsematrix.cpp.
References TPZMatrix< TVar >::Cols(), TPZVerySparseMatrix< TVar >::fExtraSparseData, and TPZMatrix< TVar >::Rows().
Referenced by TPZVerySparseMatrix< TVar >::TPZVerySparseMatrix().
|
virtual |
It makes *T the transpose of current matrix.
Definition at line 102 of file tpzverysparsematrix.cpp.
References TPZMatrix< TVar >::Cols(), TPZVerySparseMatrix< TVar >::fExtraSparseData, TPZMatrix< TVar >::Resize(), and TPZMatrix< TVar >::Rows().
Referenced by TPZVerySparseMatrix< TVar >::s().
|
inlineoverridevirtual |
It makes *T the transpose of current matrix.
Reimplemented from TPZMatrix< TVar >.
Definition at line 100 of file tpzverysparsematrix.h.
References TPZVerySparseMatrix< TVar >::ClassId(), TPZVerySparseMatrix< TVar >::Read(), TPZMatrix< TVar >::Transpose(), and TPZVerySparseMatrix< TVar >::Write().
|
overridevirtual |
Packs the object structure in a stream of bytes.
buf | Buffer which will receive the bytes |
withclassid |
Reimplemented from TPZMatrix< TVar >.
Definition at line 198 of file tpzverysparsematrix.cpp.
References TPZMatrix< TVar >::fCol, TPZMatrix< TVar >::fDecomposed, TPZMatrix< TVar >::fDefPositive, TPZVerySparseMatrix< TVar >::fExtraSparseData, TPZMatrix< TVar >::fRow, TPZMatrix< TVar >::gZero, TPZStream::Write(), and TPZVerySparseMatrix< TVar >::WriteMap().
Referenced by TPZVerySparseMatrix< TVar >::Transpose().
|
private |
Auxiliary functions only reading and writing a map as the third paremeter.
Definition at line 209 of file tpzverysparsematrix.cpp.
References TPZStream::Write().
Referenced by TPZVerySparseMatrix< TVar >::MapEnd(), and TPZVerySparseMatrix< TVar >::Write().
|
friend |
Definition at line 32 of file tpzverysparsematrix.h.
|
protected |
Save elements different from zero, of Sparse matrix.
Definition at line 122 of file tpzverysparsematrix.h.
Referenced by TPZVerySparseMatrix< TVar >::GetVal(), TPZVerySparseMatrix< TVar >::MapBegin(), TPZVerySparseMatrix< TVar >::MapEnd(), TPZVerySparseMatrix< TVar >::MultAdd(), TPZFYsmpMatrix< TVar >::operator=(), TPZVerySparseMatrix< TVar >::PutVal(), TPZVerySparseMatrix< TVar >::Read(), TPZVerySparseMatrix< TVar >::s(), TPZVerySparseMatrix< TVar >::Simetrize(), TPZVerySparseMatrix< TVar >::Transpose(), and TPZVerySparseMatrix< TVar >::Write().