NeoPZ
|
Defines a non symmetric banded matrix. Matrix. More...
#include <pzbndmat.h>
Public Member Functions | |
virtual int | Substitution (TPZFMatrix< TVar > *B) const override |
Computes Forward and Backward substitution for a "LU" decomposed matrix. More... | |
TPZFBMatrix () | |
Simple constructor. More... | |
TPZFBMatrix (const int64_t dim, const int64_t band_width=0) | |
Constructor. More... | |
TPZFBMatrix (const TPZFBMatrix< TVar > &) | |
Copy constructor. More... | |
~TPZFBMatrix () | |
Simple destructor. More... | |
template<class TVar2 > | |
void | CopyFrom (TPZFBMatrix< TVar2 > &orig) |
copy the values from a matrix with a different precision More... | |
void | AutoFill (int64_t nrow, int64_t ncol, int symmetric) |
int | Put (const int64_t row, const int64_t col, const TVar &value) override |
Put values with bounds checking if DEBUG variable is defined. More... | |
const TVar & | Get (const int64_t row, const int64_t col) const override |
Get value with bound checking. More... | |
TVar & | operator() (const int64_t row, const int64_t col) |
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... | |
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... | |
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... | |
TPZFBMatrix & | operator= (const TPZFBMatrix< TVar > &A) |
TPZFBMatrix | operator+ (const TPZFBMatrix< TVar > &A) const |
TPZFBMatrix | operator- (const TPZFBMatrix< TVar > &A) const |
TPZFBMatrix & | operator+= (const TPZFBMatrix< TVar > &A) |
TPZFBMatrix & | operator-= (const TPZFBMatrix< TVar > &A) |
TPZFBMatrix | operator* (const TVar val) const |
TPZFBMatrix & | operator*= (const TVar val) |
TPZFBMatrix | operator- () const |
int64_t | Dim () const override |
Returns the dimension of the matrix if the matrix is square. More... | |
int64_t | GetBandLower () const |
Returns band size. More... | |
int64_t | GetBandUpper () const |
int64_t | GetBand () const |
int | SetBand (const int64_t newBand) |
Sets band size. More... | |
int | Resize (const int64_t newRows, const int64_t newCols) override |
Redimension the matrix preserving its elements. More... | |
int | Redim (const int64_t newRows, const int64_t newCols) override |
Redimension the matrix and make zero its elements. More... | |
int | Zero () override |
Zeroes the matrix. More... | |
void | Transpose (TPZMatrix< TVar > *const T) const override |
It makes *T the transpose of current matrix. More... | |
int | ClassId () const override |
Define the class id associated with the class. More... | |
template<> | |
TPZFBMatrix< std::complex< float > > & | operator*= (const std::complex< float > value) |
template<> | |
TPZFBMatrix< std::complex< double > > & | operator*= (const std::complex< double > value) |
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... | |
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 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_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 | 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 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 | |
Private Member Functions | |
int64_t | Index (int64_t i, int64_t j) const |
int | Clear () override |
It clears data structure. More... | |
Private Attributes | |
TPZVec< TVar > | fElem |
int64_t | fBandLower |
int64_t | fBandUpper |
Friends | |
class | TPZFBMatrix< float > |
class | TPZFBMatrix< double > |
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... | |
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... | |
Defines a non symmetric banded matrix. Matrix.
Definition at line 25 of file pzbndmat.h.
TPZFBMatrix< TVar >::TPZFBMatrix | ( | ) |
Simple constructor.
Definition at line 27 of file pzbndmat.cpp.
References TPZFBMatrix< TVar >::fBandLower, and TPZFBMatrix< TVar >::fBandUpper.
TPZFBMatrix< TVar >::TPZFBMatrix | ( | const int64_t | dim, |
const int64_t | band_width = 0 |
||
) |
Constructor.
dim | Initial dimension of banded matrix |
band_width | Initial band width of banded matrix |
Definition at line 39 of file pzbndmat.cpp.
TPZFBMatrix< TVar >::TPZFBMatrix | ( | const TPZFBMatrix< TVar > & | A | ) |
Copy constructor.
Definition at line 51 of file pzbndmat.cpp.
TPZFBMatrix< TVar >::~TPZFBMatrix | ( | ) |
Simple destructor.
Definition at line 62 of file pzbndmat.cpp.
void TPZFBMatrix< TVar >::AutoFill | ( | int64_t | nrow, |
int64_t | ncol, | ||
int | symmetric | ||
) |
Definition at line 351 of file pzbndmat.cpp.
References DebugStop, TPZFBMatrix< TVar >::GetVal(), TPZFBMatrix< TVar >::PutVal(), and val().
Referenced by TPZFBMatrix< TVar >::CopyFrom().
|
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
Implements TPZSavable.
Definition at line 634 of file pzbndmat.cpp.
References TPZMatrix< TVar >::ClassId(), and Hash().
Referenced by TPZFBMatrix< TVar >::GetBand().
|
overrideprivatevirtual |
It clears data structure.
Reimplemented from TPZMatrix< TVar >.
Definition at line 645 of file pzbndmat.cpp.
References TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, TPZMatrix< TVar >::fCol, TPZFBMatrix< TVar >::fElem, TPZMatrix< TVar >::fRow, and TPZVec< T >::resize().
Referenced by TPZFBMatrix< TVar >::Index().
|
inline |
copy the values from a matrix with a different precision
Definition at line 50 of file pzbndmat.h.
References TPZFBMatrix< TVar >::AutoFill(), TPZMatrix< TVar >::CopyFrom(), TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, TPZFBMatrix< TVar >::fElem, TPZFBMatrix< TVar >::Get(), TPZFBMatrix< TVar >::GetVal(), TPZFBMatrix< TVar >::MultAdd(), TPZFBMatrix< TVar >::operator()(), TPZFBMatrix< TVar >::operator*(), TPZFBMatrix< TVar >::operator*=(), TPZFBMatrix< TVar >::operator+(), TPZFBMatrix< TVar >::operator+=(), TPZFBMatrix< TVar >::operator-(), TPZFBMatrix< TVar >::operator-=(), TPZFBMatrix< TVar >::operator=(), TPZFBMatrix< TVar >::Put(), TPZFBMatrix< TVar >::PutVal(), TPZVec< T >::resize(), TPZFBMatrix< TVar >::s(), TPZVec< T >::size(), and val().
|
inlineoverridevirtual |
Returns the dimension of the matrix if the matrix is square.
If the matrix is not square, returns an error
Reimplemented from TPZMatrix< TVar >.
Definition at line 103 of file pzbndmat.h.
References TPZMatrix< TVar >::Rows().
Referenced by TPZFBMatrix< TVar >::Get(), TPZFBMatrix< TVar >::MultAdd(), TPZFBMatrix< TVar >::operator+(), TPZFBMatrix< TVar >::operator+=(), TPZFBMatrix< TVar >::operator-(), TPZFBMatrix< TVar >::operator-=(), TPZFBMatrix< TVar >::Put(), TPZFBMatrix< TVar >::SetBand(), and TPZFBMatrix< TVar >::Transpose().
|
overridevirtual |
Get value with bound checking.
row | Row number. |
col | Column number. |
Reimplemented from TPZMatrix< TVar >.
Definition at line 93 of file pzbndmat.cpp.
References TPZFBMatrix< TVar >::Dim(), TPZMatrix< TVar >::Error(), and TPZFBMatrix< TVar >::GetVal().
Referenced by TPZFBMatrix< TVar >::CopyFrom().
|
inline |
Definition at line 113 of file pzbndmat.h.
References TPZFBMatrix< TVar >::ClassId(), DebugStop, TPZMatrix< TVar >::Decompose_LU(), TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, TPZFBMatrix< TVar >::Redim(), TPZFBMatrix< TVar >::Resize(), TPZFBMatrix< TVar >::SetBand(), TPZFBMatrix< TVar >::Transpose(), and TPZFBMatrix< TVar >::Zero().
|
inline |
Returns band size.
Definition at line 105 of file pzbndmat.h.
References TPZFBMatrix< TVar >::fBandLower.
|
inline |
Definition at line 109 of file pzbndmat.h.
References TPZFBMatrix< TVar >::fBandUpper.
|
inlineoverridevirtual |
Get values without bounds checking
This method is faster than "Get" if DEBUG is defined.
Reimplemented from TPZMatrix< TVar >.
Definition at line 190 of file pzbndmat.h.
References DebugStop, TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, TPZMatrix< TVar >::fCol, TPZFBMatrix< TVar >::fElem, TPZMatrix< TVar >::fRow, TPZFBMatrix< TVar >::Index(), and TPZFBMatrix< TVar >::Zero().
Referenced by TPZFBMatrix< TVar >::AutoFill(), TPZFBMatrix< TVar >::CopyFrom(), TPZFBMatrix< TVar >::Get(), TPZFBMatrix< TVar >::MultAdd(), and TPZFBMatrix< TVar >::Transpose().
|
inlineprivate |
Definition at line 148 of file pzbndmat.h.
References TPZFBMatrix< TVar >::Clear(), TPZFBMatrix< TVar >::fBandLower, and TPZFBMatrix< TVar >::fBandUpper.
Referenced by TPZFBMatrix< TVar >::GetVal(), TPZFBMatrix< TVar >::operator()(), and TPZFBMatrix< TVar >::PutVal().
|
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 221 of file pzbndmat.cpp.
References TPZMatrix< TVar >::Cols(), TPZFBMatrix< TVar >::Dim(), TPZMatrix< TVar >::Error(), TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, TPZFBMatrix< TVar >::GetVal(), TPZFMatrix< TVar >::GetVal(), MAX, MIN, TPZMatrix< TVar >::PrepareZ(), TPZFMatrix< TVar >::PutVal(), TPZMatrix< TVar >::Rows(), and val().
Referenced by TPZFBMatrix< TVar >::CopyFrom().
|
inline |
Definition at line 206 of file pzbndmat.h.
References DebugStop, TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, TPZFBMatrix< TVar >::fElem, TPZFBMatrix< TVar >::Index(), and TPZFBMatrix< TVar >::Zero().
Referenced by TPZFBMatrix< TVar >::CopyFrom(), and TPZFBMatrix< TVar >::s().
TPZFBMatrix< TVar > TPZFBMatrix< TVar >::operator* | ( | const TVar | val | ) | const |
Definition at line 285 of file pzbndmat.cpp.
References TPZFBMatrix< TVar >::fElem, test::res, and TPZVec< T >::size().
Referenced by TPZFBMatrix< TVar >::CopyFrom().
TPZFBMatrix< TVar > & TPZFBMatrix< TVar >::operator*= | ( | const TVar | val | ) |
Definition at line 337 of file pzbndmat.cpp.
References TPZFBMatrix< TVar >::fElem, and TPZVec< T >::size().
Referenced by TPZFBMatrix< TVar >::CopyFrom().
TPZFBMatrix< std::complex< float > > & TPZFBMatrix< std::complex< float > >::operator*= | ( | const std::complex< float > | value | ) |
Definition at line 305 of file pzbndmat.cpp.
References TPZFBMatrix< TVar >::fElem, and TPZVec< T >::size().
TPZFBMatrix< std::complex< double > > & TPZFBMatrix< std::complex< double > >::operator*= | ( | const std::complex< double > | value | ) |
Definition at line 320 of file pzbndmat.cpp.
References TPZFBMatrix< TVar >::fElem, and TPZVec< T >::size().
TPZFBMatrix< TVar > TPZFBMatrix< TVar >::operator+ | ( | const TPZFBMatrix< TVar > & | A | ) | const |
Definition at line 129 of file pzbndmat.cpp.
References TPZFBMatrix< TVar >::Dim(), TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, TPZFBMatrix< TVar >::fElem, test::res, and TPZVec< T >::size().
Referenced by TPZFBMatrix< TVar >::CopyFrom().
TPZFBMatrix< TVar > & TPZFBMatrix< TVar >::operator+= | ( | const TPZFBMatrix< TVar > & | A | ) |
Definition at line 177 of file pzbndmat.cpp.
References TPZFBMatrix< TVar >::Dim(), TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, TPZFBMatrix< TVar >::fElem, and TPZVec< T >::size().
Referenced by TPZFBMatrix< TVar >::CopyFrom().
TPZFBMatrix< TVar > TPZFBMatrix< TVar >::operator- | ( | const TPZFBMatrix< TVar > & | A | ) | const |
Definition at line 150 of file pzbndmat.cpp.
References TPZFBMatrix< TVar >::Dim(), TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, TPZFBMatrix< TVar >::fElem, test::res, and TPZVec< T >::size().
TPZFBMatrix< TVar > TPZFBMatrix< TVar >::operator- | ( | ) | const |
Definition at line 269 of file pzbndmat.cpp.
Referenced by TPZFBMatrix< TVar >::CopyFrom().
TPZFBMatrix< TVar > & TPZFBMatrix< TVar >::operator-= | ( | const TPZFBMatrix< TVar > & | A | ) |
Definition at line 197 of file pzbndmat.cpp.
References TPZFBMatrix< TVar >::Dim(), TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, TPZFBMatrix< TVar >::fElem, and TPZVec< T >::size().
Referenced by TPZFBMatrix< TVar >::CopyFrom().
TPZFBMatrix< TVar > & TPZFBMatrix< TVar >::operator= | ( | const TPZFBMatrix< TVar > & | A | ) |
Definition at line 111 of file pzbndmat.cpp.
References TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, and TPZFBMatrix< TVar >::fElem.
Referenced by TPZFBMatrix< TVar >::CopyFrom().
|
overridevirtual |
Put values with bounds checking if DEBUG variable is defined.
row | Row number. |
col | Column number. |
value | Value being put. |
Reimplemented from TPZMatrix< TVar >.
Definition at line 73 of file pzbndmat.cpp.
References DebugStop, TPZFBMatrix< TVar >::Dim(), TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, and TPZFBMatrix< TVar >::PutVal().
Referenced by TPZFBMatrix< TVar >::CopyFrom().
|
inlineoverridevirtual |
Put values without bounds checking
This method is faster than "Put" if DEBUG is defined.
Reimplemented from TPZMatrix< TVar >.
Definition at line 170 of file pzbndmat.h.
References DebugStop, TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, TPZFBMatrix< TVar >::fElem, TPZFBMatrix< TVar >::Index(), and IsZero().
Referenced by TPZFBMatrix< TVar >::AutoFill(), TPZFBMatrix< TVar >::CopyFrom(), TPZFBMatrix< TVar >::Put(), and TPZMatrixMarket::Read().
|
overridevirtual |
Redimension the matrix and make zero its elements.
Reimplemented from TPZMatrix< TVar >.
Definition at line 412 of file pzbndmat.cpp.
References TPZMatrix< TVar >::Error(), TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, TPZFBMatrix< TVar >::fElem, TPZMatrix< TVar >::Redim(), TPZVec< T >::Resize(), and TPZFBMatrix< TVar >::Zero().
Referenced by TPZFBMatrix< TVar >::GetBand(), TPZMatrixMarket::Read(), and TPZFBMatrix< TVar >::Resize().
|
overridevirtual |
Redimension the matrix preserving its elements.
Reimplemented from TPZMatrix< TVar >.
Definition at line 396 of file pzbndmat.cpp.
References TPZMatrix< TVar >::Error(), and TPZFBMatrix< TVar >::Redim().
Referenced by TPZFBMatrix< TVar >::GetBand().
|
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 216 of file pzbndmat.h.
References TPZFBMatrix< TVar >::operator()().
Referenced by TPZFBMatrix< TVar >::CopyFrom().
int TPZFBMatrix< TVar >::SetBand | ( | const int64_t | newBand | ) |
Sets band size.
newBand | New band size |
Definition at line 458 of file pzbndmat.cpp.
References TPZFBMatrix< TVar >::Dim(), TPZMatrix< TVar >::Error(), TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, TPZFBMatrix< TVar >::fElem, TPZVec< T >::resize(), and TPZFBMatrix< TVar >::Zero().
Referenced by TPZFBMatrix< TVar >::GetBand(), and TPZMatrixMarket::Read().
|
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 624 of file pzbndmat.cpp.
References DebugStop, ELU, TPZMatrix< TVar >::fDecomposed, and TPZMatrix< TVar >::Substitution().
Referenced by TPZFBMatrix< TVar >::Transpose().
|
overridevirtual |
It makes *T the transpose of current matrix.
Reimplemented from TPZMatrix< TVar >.
Definition at line 477 of file pzbndmat.cpp.
References TPZMatrix< TVar >::Cols(), DebugStop, TPZMatrix< TVar >::Decompose_LU(), TPZFBMatrix< TVar >::Dim(), ELU, TPZMatrix< TVar >::Error(), TPZFBMatrix< TVar >::fBandLower, TPZFBMatrix< TVar >::fBandUpper, TPZMatrix< TVar >::fDecomposed, TPZFBMatrix< TVar >::fElem, TPZFBMatrix< TVar >::GetVal(), MAX, MIN, TPZMatrix< TVar >::PutVal(), TPZMatrix< TVar >::Resize(), TPZMatrix< TVar >::Rows(), TPZFMatrix< TVar >::s(), and TPZFBMatrix< TVar >::Substitution().
Referenced by TPZFBMatrix< TVar >::GetBand().
|
overridevirtual |
Zeroes the matrix.
Reimplemented from TPZMatrix< TVar >.
Definition at line 439 of file pzbndmat.cpp.
References TPZMatrix< TVar >::fDecomposed, TPZFBMatrix< TVar >::fElem, and TPZVec< T >::size().
Referenced by TPZFBMatrix< TVar >::GetBand(), TPZFBMatrix< TVar >::GetVal(), TPZFBMatrix< TVar >::operator()(), TPZFBMatrix< TVar >::Redim(), and TPZFBMatrix< TVar >::SetBand().
|
friend |
Definition at line 46 of file pzbndmat.h.
|
friend |
Definition at line 45 of file pzbndmat.h.
|
private |
Definition at line 155 of file pzbndmat.h.
Referenced by TPZFBMatrix< TVar >::Clear(), TPZFBMatrix< TVar >::CopyFrom(), TPZFBMatrix< TVar >::GetBand(), TPZFBMatrix< TVar >::GetBandLower(), TPZFBMatrix< TVar >::GetVal(), TPZFBMatrix< TVar >::Index(), TPZFBMatrix< TVar >::MultAdd(), TPZFBMatrix< TVar >::operator()(), TPZFBMatrix< TVar >::operator+(), TPZFBMatrix< TVar >::operator+=(), TPZFBMatrix< TVar >::operator-(), TPZFBMatrix< TVar >::operator-=(), TPZFBMatrix< TVar >::operator=(), TPZFBMatrix< TVar >::Put(), TPZFBMatrix< TVar >::PutVal(), TPZFBMatrix< TVar >::Redim(), TPZFBMatrix< TVar >::SetBand(), TPZFBMatrix< TVar >::TPZFBMatrix(), and TPZFBMatrix< TVar >::Transpose().
|
private |
Definition at line 155 of file pzbndmat.h.
Referenced by TPZFBMatrix< TVar >::Clear(), TPZFBMatrix< TVar >::CopyFrom(), TPZFBMatrix< TVar >::GetBand(), TPZFBMatrix< TVar >::GetBandUpper(), TPZFBMatrix< TVar >::GetVal(), TPZFBMatrix< TVar >::Index(), TPZFBMatrix< TVar >::MultAdd(), TPZFBMatrix< TVar >::operator()(), TPZFBMatrix< TVar >::operator+(), TPZFBMatrix< TVar >::operator+=(), TPZFBMatrix< TVar >::operator-(), TPZFBMatrix< TVar >::operator-=(), TPZFBMatrix< TVar >::operator=(), TPZFBMatrix< TVar >::Put(), TPZFBMatrix< TVar >::PutVal(), TPZFBMatrix< TVar >::Redim(), TPZFBMatrix< TVar >::SetBand(), TPZFBMatrix< TVar >::TPZFBMatrix(), and TPZFBMatrix< TVar >::Transpose().
|
private |
Definition at line 154 of file pzbndmat.h.
Referenced by TPZFBMatrix< TVar >::Clear(), TPZFBMatrix< TVar >::CopyFrom(), TPZFBMatrix< TVar >::GetVal(), TPZFBMatrix< TVar >::operator()(), TPZFBMatrix< TVar >::operator*(), TPZFBMatrix< TVar >::operator*=(), TPZFBMatrix< TVar >::operator+(), TPZFBMatrix< TVar >::operator+=(), TPZFBMatrix< TVar >::operator-(), TPZFBMatrix< TVar >::operator-=(), TPZFBMatrix< TVar >::operator=(), TPZFBMatrix< TVar >::PutVal(), TPZFBMatrix< TVar >::Redim(), TPZFBMatrix< TVar >::SetBand(), TPZFBMatrix< TVar >::Transpose(), and TPZFBMatrix< TVar >::Zero().