NeoPZ
pzmgsolver.cpp
Go to the documentation of this file.
1 
6 #include "pzmgsolver.h"
7 #include "pztransfer.h"
9 
10 using namespace std;
11 
12 template <class TVar>
14  TPZAutoPointer<TPZMatrix<TVar> > refmat) :
15 TPZRegisterClassId(&TPZMGSolver::ClassId),TPZMatrixSolver<TVar>(refmat), fTransfer(trf)
16 {
17  this->fCoarse = (TPZMatrixSolver<TVar> *) sol.Clone();
18  this->fNVar = nvar;
19 }
20 
21 template <class TVar>
24 TPZMatrixSolver<TVar>(), fTransfer(trf)
25 {
26  this->fCoarse = (TPZMatrixSolver<TVar> *) sol.Clone();
27  // fTransfer = new TPZMatrixSolver::TPZContainer(trf);
28  this->fNVar = nvar;
29 }
30 
31 template <class TVar>
33  if((!this->Matrix() && residual != 0) || !TransferMatrix()) {
34  cout << "TPZMGSolver::Solve called without a matrix pointer\n";
35  DebugStop();
36  }
38  if(result.Rows() != tr->Cols() || result.Cols() != F.Cols()) {
39  result.Redim(tr->Rows(),F.Cols());
40  }
41 
42  TPZFMatrix<TVar> FCoarse,UCoarse;
43  tr->Multiply(F,FCoarse,0);
44  double norm = Norm(FCoarse);
45  if(norm != 0.)
46  {
47  fCoarse->Solve(FCoarse,UCoarse);
48  }
49  tr->Multiply(UCoarse,result,1);
50  if(residual) this->Matrix()->Residual(F,result,*residual);
51 }
52 
53 template <class TVar>
55 TPZMatrixSolver<TVar>(copy), fTransfer(copy.fTransfer) {
56  fCoarse = (TPZMatrixSolver<TVar> *) copy.fCoarse->Clone();
57  fNVar = copy.fNVar;
58 }
59 
60 template <class TVar>
62  return new TPZMGSolver<TVar>(*this);
63 }
64 
65 template <class TVar>
67  delete fCoarse;
68 }
69 
70 template <class TVar>
73  fTransfer = reset;
74 }
75 
76 template <class TVar>
78  fTransfer = Refmat;
79 }
80 
81 template <class TVar>
82 void TPZMGSolver<TVar>::Write(TPZStream &buf, int withclassid) const
83 {
84  TPZMatrixSolver<TVar>::Write(buf, withclassid);
86  buf.Write(&fNVar);
88 }
89 
90 template <class TVar>
91 void TPZMGSolver<TVar>::Read(TPZStream &buf, void *context)
92 {
93  TPZMatrixSolver<TVar>::Read(buf, context);
95  buf.Read(&fNVar, 1);
96  fTransfer = TPZAutoPointerDynamicCast<TPZMatrix<TVar>>(TPZPersistenceManager::GetAutoPointer(&buf));
97 }
98 
99 template class TPZMGSolver<float>;
100 template class TPZMGSolver<double>;
101 template class TPZMGSolver<long double>;
102 
103 
104 #ifndef BORLAND
105 template class TPZRestoreClass<TPZMGSolver<REAL>>;
106 #endif
Contains the TPZMGSolver class which represents a solution process in three steps.
~TPZMGSolver()
Default destructor.
Definition: pzmgsolver.cpp:66
Defines a abstract class of solvers which will be used by matrix classes. Solver. ...
Definition: pzmatrix.h:32
void Write(TPZStream &buf, int withclassid) const override
Writes this object to the TPZStream buffer. Include the classid if withclassid = true.
Definition: pzsolve.cpp:65
TPZSolver< TVar > * Clone() const override
Clones the current object returning a pointer of type TPZSolver.
Definition: pzmgsolver.cpp:61
Defines a class of matrix solvers. Solver.
Definition: pzanalysis.h:24
TPZAutoPointer< TPZMatrix< TVar > > TransferMatrix()
Gets the transfer matrix.
Definition: pzmgsolver.h:44
static TPZSavable * GetInstance(const int64_t &objId)
Contains the TPZTransfer class which implements a rectangular sparse block matrix.
static TPZAutoPointer< TPZSavable > GetAutoPointer(const int64_t &objId)
virtual void Write(const bool val)
Definition: TPZStream.cpp:8
TVar Norm(const TPZFMatrix< TVar > &A)
Returns the norm of the matrix A.
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
TPZMGSolver()
Default constructor.
Definition: pzmgsolver.h:25
void Read(TPZStream &buf, void *context) override
read objects from the stream
Definition: pzsolve.cpp:106
int64_t Rows() const
Returns number of rows.
Definition: pzmatrix.h:803
Represents a solution process in three steps: transfer of the residual, execute a solver on the coars...
Definition: pzmgsolver.h:21
Full matrix class. Matrix.
Definition: pzfmatrix.h:32
int Redim(const int64_t newRows, const int64_t newCols) override
Redimension a matrix and ZERO your elements.
Definition: pzfmatrix.h:616
void Write(TPZStream &buf, int withclassid) const override
Writes this object to the TPZStream buffer. Include the classid if withclassid = true.
Definition: pzmgsolver.cpp:82
virtual TPZSolver * Clone() const =0
Clones the current object returning a pointer of type TPZSolver.
void Solve(const TPZFMatrix< TVar > &F, TPZFMatrix< TVar > &result, TPZFMatrix< TVar > *residual=0) override
Solves the system of linear equations.
Definition: pzmgsolver.cpp:32
TPZAutoPointer< TPZMatrix< TVar > > Matrix() const
Returns a pointer to TPZMatrix<>
Definition: pzsolve.h:138
int ClassId() const override
Saveable specific methods.
Definition: pzmgsolver.h:69
TPZMatrixSolver< TVar > * fCoarse
Definition: pzmgsolver.h:61
int64_t Cols() const
Returns number of cols.
Definition: pzmatrix.h:809
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
void ResetTransferMatrix()
Clean the transfer matrix.
Definition: pzmgsolver.cpp:71
void SetTransferMatrix(TPZAutoPointer< TPZMatrix< TVar > > Refmat)
Sets the transfer matrix.
Definition: pzmgsolver.cpp:77
static void WritePointer(const TPZSavable *obj, TPZStream *stream)
TPZAutoPointer< TPZMatrix< TVar > > fTransfer
Transfer matrix.
Definition: pzmgsolver.h:64
Implements an interface to register a class id and a restore function. Persistence.
Definition: TPZSavable.h:150
void Read(TPZStream &buf, void *context) override
read objects from the stream
Definition: pzmgsolver.cpp:91
virtual void Read(bool &val)
Definition: TPZStream.cpp:91
Root matrix class (abstract). Matrix.
Definition: pzmatrix.h:60
This class implements a reference counter mechanism to administer a dynamically allocated object...