NeoPZ
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Static Private Member Functions | List of all members
TPZFYsmpMatrix< TVar > Class Template Reference

Implements a non symmetric sparse matrix (Yale Sparse Matrix Storage). Matrix. More...

#include <pzysmp.h>

Inheritance diagram for TPZFYsmpMatrix< TVar >:
[legend]
Collaboration diagram for TPZFYsmpMatrix< TVar >:
[legend]

Classes

struct  TPZMThread
 An auxiliary structure to hold the data of the subset
of equations used to multiply in a multi-threaded environment. More...
 

Public Member Functions

 TPZFYsmpMatrix ()
 
 TPZFYsmpMatrix (const int64_t rows, const int64_t cols)
 
 TPZFYsmpMatrix (const TPZVerySparseMatrix< TVar > &cp)
 
TPZFYsmpMatrixoperator= (const TPZFYsmpMatrix< TVar > &copy)
 
TPZFYsmpMatrixoperator= (const TPZVerySparseMatrix< TVar > &cp)
 
virtual ~TPZFYsmpMatrix ()
 
void AutoFill (int64_t nrow, int64_t ncol, int symmetric)
 Fill matrix storage with randomic values. More...
 
virtual const TVar & GetVal (const int64_t row, const int64_t col) const override
 Get the matrix entry at (row,col) without bound checking. More...
 
int64_t NumTerms ()
 
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...
 
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 MultAddMT (const TPZFMatrix< TVar > &x, const TPZFMatrix< TVar > &y, TPZFMatrix< TVar > &z, const TVar alpha=1., const TVar beta=0., const int opt=0)
 
virtual int GetSub (const int64_t sRow, const int64_t sCol, const int64_t rowSize, const int64_t colSize, TPZFMatrix< TVar > &A) const override
 Gets submatrix storing it on Target. More...
 
void GetSub (const TPZVec< int64_t > &indices, TPZFMatrix< TVar > &block) const override
 Extract the block indicated by the indices from the matrix. More...
 
virtual void SetData (int64_t *IA, int64_t *JA, TVar *A)
 Pass the data to the class. More...
 
virtual void SetData (TPZVec< int64_t > &IA, TPZVec< int64_t > &JA, TPZVec< TVar > &A)
 Pass the data to the class. More...
 
virtual void Print (const char *title, std::ostream &out=std::cout, const MatrixOutputFormat form=EFormatted) const override
 Print the matrix along with a identification title. More...
 
virtual void AddKelOld (TPZFMatrix< TVar > &elmat, TPZVec< int > &destinationindex)
 Add a contribution of a stiffness matrix putting it on destination indexes position. More...
 
virtual void AddKel (TPZFMatrix< TVar > &elmat, TPZVec< int64_t > &destinationindex) override
 Add a contribution of a stiffness matrix. More...
 
virtual void AddKel (TPZFMatrix< TVar > &elmat, TPZVec< int64_t > &sourceindex, TPZVec< int64_t > &destinationindex) override
 Add a contribution of a stiffness matrix. More...
 
void MultiplyDummy (TPZFYsmpMatrix< TVar > &B, TPZFYsmpMatrix< TVar > &Res)
 
virtual int Zero () override
 Zeroes the matrix. More...
 
int ClassId () const override
 Define the class id associated with the class. More...
 
Solvers

Linear system solvers.

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) override
 Solves the linear system using Jacobi method.
. More...
 
void SolveSOR (int64_t &numiterations, const TPZFMatrix< TVar > &rhs, TPZFMatrix< TVar > &x, TPZFMatrix< TVar > *residual, TPZFMatrix< TVar > &scratch, const REAL overrelax, REAL &tol, const int FromCurrent=0, const int direction=1) override
 Solves the linear system using Successive Over Relaxation method (Gauss Seidel).
. More...
 
Factorization

Those member functions are related to matrices factorization

virtual int Decompose_LU (std::list< int64_t > &singular) override
 Decomposes the current matrix using LU decomposition. More...
 
virtual int Decompose_LU () override
 
Substitutions

Substitutions forward and backward

virtual int Substitution (TPZFMatrix< TVar > *B) const override
 Computes Forward and Backward substitution for a "LU" decomposed matrix. More...
 
- 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 > &copy)
 
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...
 
virtual TVar & s (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...
 
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...
 
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...
 
virtual void Transpose (TPZMatrix< TVar > *const T) const
 It makes *T the transpose of current 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 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 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_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 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...
 
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
 
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 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 Member Functions

void InitializeData ()
 Implements a initialization method for the sparse structure. It sets the initial value for the fIA and fJA. More...
 
- 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...
 

Protected Attributes

TPZVec< int64_t > fIA
 
TPZVec< int64_t > fJA
 
TPZVec< TVar > fA
 
TPZVec< TVar > fDiag
 
int fSymmetric
 
- 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 ComputeDiagonal ()
 
void RowLUUpdate (int64_t sourcerow, int64_t destrow)
 

Static Private Member Functions

static void * ExecuteMT (void *entrydata)
 

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 TPZSavableCreateInstance (const int &classId)
 
- 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...
 

Detailed Description

template<class TVar>
class TPZFYsmpMatrix< TVar >

Implements a non symmetric sparse matrix (Yale Sparse Matrix Storage). Matrix.

Defines operations on general sparse matrices stored in the (old) Yale Sparse Matrix Package format.

Definition at line 45 of file pzysmp.h.

Constructor & Destructor Documentation

◆ TPZFYsmpMatrix() [1/3]

template<class TVar >
TPZFYsmpMatrix< TVar >::TPZFYsmpMatrix ( )

Definition at line 70 of file pzysmp.cpp.

◆ TPZFYsmpMatrix() [2/3]

template<class TVar >
TPZFYsmpMatrix< TVar >::TPZFYsmpMatrix ( const int64_t  rows,
const int64_t  cols 
)

◆ TPZFYsmpMatrix() [3/3]

template<class TVar >
TPZFYsmpMatrix< TVar >::TPZFYsmpMatrix ( const TPZVerySparseMatrix< TVar > &  cp)

Definition at line 59 of file pzysmp.cpp.

◆ ~TPZFYsmpMatrix()

template<class TVar >
TPZFYsmpMatrix< TVar >::~TPZFYsmpMatrix ( )
virtual

Definition at line 349 of file pzysmp.cpp.

Member Function Documentation

◆ AddKel() [1/2]

template<class TVar >
void TPZFYsmpMatrix< TVar >::AddKel ( TPZFMatrix< TVar > &  elmat,
TPZVec< int64_t > &  destinationindex 
)
overridevirtual

Add a contribution of a stiffness matrix.

Parameters
elmatElement matrix to be contributed
destinationindexContains destine indexes on current matrix

Reimplemented from TPZMatrix< TVar >.

Definition at line 162 of file pzysmp.cpp.

References TPZFYsmpMatrix< TVar >::fA, TPZFYsmpMatrix< TVar >::fIA, TPZFYsmpMatrix< TVar >::fJA, TPZFMatrix< TVar >::GetVal(), and TPZMatrix< TVar >::Rows().

Referenced by TPZFYsmpMatrix< TVar >::NumTerms().

◆ AddKel() [2/2]

template<class TVar >
void TPZFYsmpMatrix< TVar >::AddKel ( TPZFMatrix< TVar > &  elmat,
TPZVec< int64_t > &  sourceindex,
TPZVec< int64_t > &  destinationindex 
)
overridevirtual

Add a contribution of a stiffness matrix.

Parameters
elmatElement matrix to be contributed
sourceindexContains source indexes on current matrix
destinationindexContains destine indexes on current matrix

Reimplemented from TPZMatrix< TVar >.

Definition at line 205 of file pzysmp.cpp.

References TPZFYsmpMatrix< TVar >::fA, TPZFYsmpMatrix< TVar >::fIA, TPZFYsmpMatrix< TVar >::fJA, TPZFMatrix< TVar >::GetVal(), and TPZVec< T >::NElements().

◆ AddKelOld()

template<class TVar >
void TPZFYsmpMatrix< TVar >::AddKelOld ( TPZFMatrix< TVar > &  elmat,
TPZVec< int > &  destinationindex 
)
virtual

Add a contribution of a stiffness matrix putting it on destination indexes position.

Parameters
elmatMember stiffness matrix beeing added
destinationindexPositioning of such members on global stiffness matrix

Definition at line 248 of file pzysmp.cpp.

References TPZVec< T >::begin(), TPZFYsmpMatrix< TVar >::fA, TPZFYsmpMatrix< TVar >::fIA, TPZFYsmpMatrix< TVar >::fJA, and TPZVec< T >::NElements().

Referenced by TPZFYsmpMatrix< TVar >::NumTerms().

◆ AutoFill()

template<class TVar >
void TPZFYsmpMatrix< TVar >::AutoFill ( int64_t  nrow,
int64_t  ncol,
int  symmetric 
)

Fill matrix storage with randomic values.

This method use GetVal and PutVal which are implemented by each type matrices

Definition at line 899 of file pzysmp.cpp.

References TPZMatrix< TVar >::AutoFill(), DebugStop, TPZVec< T >::end(), TPZStack< T, NumExtAlloc >::Push(), TPZMatrix< TVar >::Resize(), TPZFYsmpMatrix< TVar >::SetData(), and TPZVec< T >::size().

◆ ClassId()

template<class TVar >
int TPZFYsmpMatrix< TVar >::ClassId ( ) const
overridevirtual

Define the class id associated with the class.

This id has to be unique for all classes A non unique id is flagged at the startup of the program

Reimplemented from TPZMatrix< TVar >.

Definition at line 1042 of file pzysmp.cpp.

References TPZMatrix< TVar >::ClassId(), and Hash().

Referenced by TPZFYsmpMatrix< TVar >::NumTerms().

◆ ComputeDiagonal()

template<class TVar >
void TPZFYsmpMatrix< TVar >::ComputeDiagonal ( )
private

◆ Decompose_LU() [1/2]

template<class TVar >
int TPZFYsmpMatrix< TVar >::Decompose_LU ( std::list< int64_t > &  singular)
overridevirtual

Decomposes the current matrix using LU decomposition.

Decomposes the current matrix using LU decomposition.

Reimplemented from TPZMatrix< TVar >.

Definition at line 941 of file pzysmp.cpp.

References TPZFYsmpMatrix< TVar >::Decompose_LU().

◆ Decompose_LU() [2/2]

template<class TVar >
int TPZFYsmpMatrix< TVar >::Decompose_LU ( )
overridevirtual

◆ ExecuteMT()

template<class TVar >
void * TPZFYsmpMatrix< TVar >::ExecuteMT ( void *  entrydata)
staticprivate

◆ GetSub() [1/2]

template<class TVar >
int TPZFYsmpMatrix< TVar >::GetSub ( const int64_t  sRow,
const int64_t  sCol,
const int64_t  rowSize,
const int64_t  colSize,
TPZFMatrix< TVar > &  Target 
) const
overridevirtual

Gets submatrix storing it on Target.

Parameters
sRowSpecifies starting row on current object
sColSpecifies starting column on current object.
rowSizeSpecifies the amount of rows from sRow
colSizeSpecifies the amount of columns from sCol
TargetThe matrix to be aquired.

Reimplemented from TPZMatrix< TVar >.

Definition at line 803 of file pzysmp.cpp.

References TPZFYsmpMatrix< TVar >::fA, TPZFYsmpMatrix< TVar >::fIA, FindCol(), and TPZFYsmpMatrix< TVar >::fJA.

Referenced by TPZFYsmpMatrix< TVar >::NumTerms().

◆ GetSub() [2/2]

template<class TVar >
void TPZFYsmpMatrix< TVar >::GetSub ( const TPZVec< int64_t > &  indices,
TPZFMatrix< TVar > &  block 
) const
overridevirtual

Extract the block indicated by the indices from the matrix.

Reimplemented from TPZMatrix< TVar >.

Definition at line 820 of file pzysmp.cpp.

References dist(), TPZFYsmpMatrix< TVar >::fA, TPZFYsmpMatrix< TVar >::fIA, FindCol(), TPZFYsmpMatrix< TVar >::fJA, and TPZVec< T >::NElements().

◆ GetVal()

template<class TVar >
const TVar & TPZFYsmpMatrix< TVar >::GetVal ( const int64_t  row,
const int64_t  col 
) const
overridevirtual

Get the matrix entry at (row,col) without bound checking.

Reimplemented from TPZMatrix< TVar >.

Definition at line 363 of file pzysmp.cpp.

References TPZFYsmpMatrix< TVar >::fA, TPZFYsmpMatrix< TVar >::fIA, TPZFYsmpMatrix< TVar >::fJA, and TPZMatrix< TVar >::gZero.

Referenced by TPZFYsmpMatrix< TVar >::ComputeDiagonal(), and TPZFYsmpMatrix< TVar >::MultiplyDummy().

◆ InitializeData()

template<class TVar >
void TPZFYsmpMatrix< TVar >::InitializeData ( )
protected

Implements a initialization method for the sparse structure. It sets the initial value for the fIA and fJA.

-fIA will contain the initial positions for all the equations -fJA will contain (-1) on all its positions -fA will contain 0 on all its value

Definition at line 32 of file pzysmp.cpp.

◆ MultAdd()

template<class TVar >
void TPZFYsmpMatrix< TVar >::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
overridevirtual

It computes z = beta * y + alpha * opt(this)*x but z and x can not overlap in memory.

Parameters
xIs x on the above operation
yIs y on the above operation
zIs z on the above operation
alphaIs alpha on the above operation
betaIs beta on the above operation
optIndicates if is Transpose or not

Reimplemented from TPZMatrix< TVar >.

Definition at line 486 of file pzysmp.cpp.

References TPZMatrix< TVar >::Cols(), TPZFYsmpMatrix< TVar >::ExecuteMT(), TPZFMatrix< TVar >::g(), substruct_tst14.test::numthreads, PZ_PTHREAD_CREATE, PZ_PTHREAD_JOIN, TPZFMatrix< TVar >::Redim(), test::res, TPZMatrix< TVar >::Rows(), and TPZFMatrix< TVar >::Zero().

Referenced by TPZFYsmpMatrix< TVar >::NumTerms().

◆ MultAddMT()

template<class TVar >
void TPZFYsmpMatrix< TVar >::MultAddMT ( const TPZFMatrix< TVar > &  x,
const TPZFMatrix< TVar > &  y,
TPZFMatrix< TVar > &  z,
const TVar  alpha = 1.,
const TVar  beta = 0.,
const int  opt = 0 
)
virtual

◆ MultiplyDummy()

template<class TVar >
void TPZFYsmpMatrix< TVar >::MultiplyDummy ( TPZFYsmpMatrix< TVar > &  B,
TPZFYsmpMatrix< TVar > &  Res 
)

◆ NumTerms()

template<class TVar>
int64_t TPZFYsmpMatrix< TVar >::NumTerms ( )
inline

◆ operator=() [1/2]

template<class TVar >
TPZFYsmpMatrix< TVar > & TPZFYsmpMatrix< TVar >::operator= ( const TPZFYsmpMatrix< TVar > &  copy)

◆ operator=() [2/2]

template<class TVar >
TPZFYsmpMatrix< TVar > & TPZFYsmpMatrix< TVar >::operator= ( const TPZVerySparseMatrix< TVar > &  cp)

◆ Print()

template<class TVar>
void TPZFYsmpMatrix< TVar >::Print ( const char *  title,
std::ostream &  out = std::cout,
const MatrixOutputFormat  form = EFormatted 
) const
overridevirtual

Print the matrix along with a identification title.

Reimplemented from TPZMatrix< TVar >.

Definition at line 588 of file pzysmp.cpp.

References TPZMatrix< TVar >::Cols(), EInputFormat, TPZFYsmpMatrix< TVar >::fA, TPZFYsmpMatrix< TVar >::fIA, TPZFYsmpMatrix< TVar >::fJA, and TPZMatrix< TVar >::Rows().

Referenced by TPZFYsmpMatrix< TVar >::NumTerms().

◆ PutVal()

template<class TVar >
int TPZFYsmpMatrix< TVar >::PutVal ( const int64_t  ,
const int64_t  ,
const TVar &  val 
)
overridevirtual

Put values without bounds checking
This method is faster than "Put" if DEBUG is defined.

Reimplemented from TPZMatrix< TVar >.

Definition at line 140 of file pzysmp.cpp.

References DebugStop, TPZFYsmpMatrix< TVar >::fA, TPZFYsmpMatrix< TVar >::fIA, and TPZFYsmpMatrix< TVar >::fJA.

Referenced by TPZFYsmpMatrix< TVar >::MultiplyDummy(), and TPZFYsmpMatrix< TVar >::NumTerms().

◆ RowLUUpdate()

template<class TVar >
void TPZFYsmpMatrix< TVar >::RowLUUpdate ( int64_t  sourcerow,
int64_t  destrow 
)
private

◆ SetData() [1/2]

template<class TVar >
void TPZFYsmpMatrix< TVar >::SetData ( int64_t *  IA,
int64_t *  JA,
TVar *  A 
)
inlinevirtual

◆ SetData() [2/2]

template<class TVar >
void TPZFYsmpMatrix< TVar >::SetData ( TPZVec< int64_t > &  IA,
TPZVec< int64_t > &  JA,
TPZVec< TVar > &  A 
)
inlinevirtual

◆ SolveJacobi()

template<class TVar >
void TPZFYsmpMatrix< TVar >::SolveJacobi ( int64_t &  numiterations,
const TPZFMatrix< TVar > &  F,
TPZFMatrix< TVar > &  result,
TPZFMatrix< TVar > *  residual,
TPZFMatrix< TVar > &  scratch,
REAL &  tol,
const int  FromCurrent = 0 
)
overridevirtual

Solves the linear system using Jacobi method.
.

For symmetric decompositions lower triangular matrix is used.
Solves a system A*X = B returning X in B

Parameters
numiterationsThe number of interations for the process.
FThe right hand side of the system.
resultThe solution.
residualReturns F - A*U which is the solution residual.
scratchAvailable manipulation area on memory.
tolThe tolerance value.
FromCurrentIt starts the solution based on FromCurrent. Obtaining solution FromCurrent + 1.

Solves the linear system using Jacobi method.

Parameters
numiterationsThe number of interations for the process.
FThe right hand side of the system.
resultThe solution.
residualReturns F - A*U which is the solution residual.
scratchAvailable manipulation area on memory.
tolThe tolerance value.
FromCurrentIt starts the solution based on FromCurrent. Obtaining solution FromCurrent + 1.

Reimplemented from TPZMatrix< TVar >.

Definition at line 683 of file pzysmp.cpp.

References TPZMatrix< TVar >::Cols(), TPZFYsmpMatrix< TVar >::fDiag, TPZFMatrix< TVar >::GetVal(), Norm(), test::res, TPZMatrix< TVar >::Residual(), TPZMatrix< TVar >::Rows(), TPZVec< T >::size(), sqrt, and pzgeom::tol.

Referenced by TPZFYsmpMatrix< TVar >::NumTerms().

◆ SolveSOR()

template<class TVar >
void TPZFYsmpMatrix< TVar >::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 
)
overridevirtual

Solves the linear system using Successive Over Relaxation method (Gauss Seidel).
.

Parameters
numiterationsThe number of interations for the process.
FThe right hand side of the system.
resultThe solution.
residualReturns F - A*U which is the solution residual.
scratchAvailable manipulation area on memory.
overrelaxThe over relaxation parameter
tolThe tolerance value..
FromCurrentIt starts the solution based on FromCurrent. Obtaining solution FromCurrent + 1.
directionIndicates interaction direction, from first to last (default 1) or from last to first (-1)

Reimplemented from TPZMatrix< TVar >.

Definition at line 629 of file pzysmp.cpp.

References TPZFYsmpMatrix< TVar >::fA, TPZFYsmpMatrix< TVar >::fDiag, TPZFYsmpMatrix< TVar >::fIA, TPZFYsmpMatrix< TVar >::fJA, TPZFMatrix< TVar >::g(), TPZMatrix< TVar >::Residual(), TPZMatrix< TVar >::Rows(), sqrt, pzgeom::tol, and TPZFMatrix< TVar >::Zero().

Referenced by TPZFYsmpMatrix< TVar >::NumTerms().

◆ Substitution()

template<class TVar >
int TPZFYsmpMatrix< TVar >::Substitution ( TPZFMatrix< TVar > *  B) const
overridevirtual

Computes Forward and Backward substitution for a "LU" decomposed matrix.

Parameters
Bright hand side and result after all

Reimplemented from TPZMatrix< TVar >.

Definition at line 979 of file pzysmp.cpp.

References TPZMatrix< TVar >::Cols(), TPZFYsmpMatrix< TVar >::fA, TPZFYsmpMatrix< TVar >::fIA, FindCol(), TPZFYsmpMatrix< TVar >::fJA, and TPZMatrix< TVar >::Rows().

Referenced by TPZFYsmpMatrix< TVar >::NumTerms().

◆ Zero()

template<class TVar >
int TPZFYsmpMatrix< TVar >::Zero ( )
overridevirtual

Zeroes the matrix.

Reimplemented from TPZMatrix< TVar >.

Definition at line 664 of file pzysmp.cpp.

References TPZFYsmpMatrix< TVar >::fA, TPZFYsmpMatrix< TVar >::fDiag, and TPZVec< T >::Fill().

Referenced by TPZFYsmpMatrix< TVar >::NumTerms().

Member Data Documentation

◆ fA

template<class TVar>
TPZVec<TVar> TPZFYsmpMatrix< TVar >::fA
protected

◆ fDiag

template<class TVar>
TPZVec<TVar> TPZFYsmpMatrix< TVar >::fDiag
protected

◆ fIA

template<class TVar>
TPZVec<int64_t> TPZFYsmpMatrix< TVar >::fIA
protected

◆ fJA

template<class TVar>
TPZVec<int64_t> TPZFYsmpMatrix< TVar >::fJA
protected

◆ fSymmetric

template<class TVar>
int TPZFYsmpMatrix< TVar >::fSymmetric
protected

Definition at line 212 of file pzysmp.h.

Referenced by TPZFYsmpMatrix< TVar >::TPZFYsmpMatrix().


The documentation for this class was generated from the following files: