NeoPZ
tpzdohrprecond.h
Go to the documentation of this file.
1 
7 #ifndef TPZDOHRPRECOND_H
8 #define TPZDOHRPRECOND_H
9 
10 #include "pzmatrix.h"
11 #include "tpzautopointer.h"
12 #include "tpzdohrsubstruct.h"
13 #include "tpzdohrmatrix.h"
14 #include "tpzdohrassembly.h"
15 #include "tpzdohrassemblelist.h"
16 
17 #include <list>
18 #include <semaphore.h>
19 
20 #include "pz_pthread.h"
21 
31 template <class TVar, class TSubStruct>
32 class TPZDohrPrecond : public TPZMatrix<TVar>
33 {
35  std::list<TPZAutoPointer<TSubStruct> > fGlobal;
38 
40  int64_t fNumCoarse; //n(c)
41 
44 
46 
47  //
48 
49 public:
53  TPZDohrPrecond(const TPZDohrPrecond &copy);
56 
58 
59  // CLONEDEF(TPZDohrPrecond)
60  virtual TPZMatrix<TVar>*Clone() const override { return new TPZDohrPrecond(*this); }
61 
63  std::list<TPZAutoPointer<TSubStruct> > &Global()
64  {
65  return fGlobal;
66  }
67 
70  void Initialize();
71 
76  {
77  fNumThreads = numthreads;
78  }
79 
95  virtual void MultAdd(const TPZFMatrix<TVar> &x,const TPZFMatrix<TVar> &y, TPZFMatrix<TVar> &z, const TVar alpha,const TVar beta,const int opt) const override;
96 
98  virtual void MultAddTBB(const TPZFMatrix<TVar> &x,const TPZFMatrix<TVar> &y, TPZFMatrix<TVar> &z, const TVar alpha,const TVar beta,const int opt) const;
99 
101  void SetSolver(TPZSolver<TVar> &solver);
102 
104  void ComputeV1(const TPZFMatrix<TVar> &x, TPZFMatrix<TVar> &v1) const;
106  void ComputeV2(const TPZFMatrix<TVar> &x, TPZFMatrix<TVar> &v2) const;
107 
109  public:
110 int ClassId() const override;
111 
117  void Read(TPZStream &buf, void *context) override;
123  void Write(TPZStream &buf, int withclassid) const override;
124 
125 };
126 
127 template <class TVar, class TSubStruct>
129  return Hash("TPZDohrPrecond") ^ TPZMatrix<TVar>::ClassId() << 1 ^ TSubStruct().ClassId() << 2;
130 }
131 
132 #include <pthread.h>
133 
137 template <class TVar, class TSubStruct>
139  TPZDohrPrecondThreadV1Data() : fDohrMatrix(0), fInput(0), fOutput(0)
140  {
141  }
143  fInput(&input), fOutput(&output)
144  {
145  }
153  static void *ComputeV1(void *dataptr)
154  {
156  ptr->fDohrMatrix->ComputeV1(*(ptr->fInput),*(ptr->fOutput));
157  return dataptr;
158  }
159 };
160 
164 template <class TVar, class TSubStruct>
166 
167  TPZDohrPrecondV2SubData() : fSubStructure(0), fInput_local(0), fv2_local(0)
168  {
169  }
170 
171  TPZDohrPrecondV2SubData(int subindex, const TPZAutoPointer<TSubStruct> &substruct, TPZAutoPointer<TPZFMatrix<TVar> > res_local) : fSubStructure(substruct),
172  fInput_local(res_local)
173  {
174  fv2_local = new TPZDohrAssembleItem<TVar>(subindex, res_local->Rows());
175  }
177  TPZDohrPrecondV2SubData(const TPZDohrPrecondV2SubData<TVar, TSubStruct> &copy) : fSubStructure(copy.fSubStructure), fInput_local(copy.fInput_local),
178  fv2_local(copy.fv2_local)
179  {
180  }
181 
183  {
184  fSubStructure = copy.fSubStructure;
185  fInput_local = copy.fInput_local;
186  fv2_local = copy.fv2_local;
187  return *this;
188  }
189 
190  bool IsValid()
191  {
192  return fSubStructure;
193  }
194 
196  {
197  }
198 
205 };
206 
207 template<class TVar>
208 struct TPZDohrAssembleList;
209 
213 template <class TVar, class TSubStruct>
215  TPZDohrPrecondV2SubDataList(TPZAutoPointer<TPZDohrAssembleList<TVar> > &assemble) : fAssemblyStructure(assemble)
216  {
217  PZ_PTHREAD_MUTEX_INIT(&fAccessLock, 0, "TPZDohrPrecondV2SubDataList::TPZDohrPrecondV2SubDataList()");
218  }
220  {
221  PZ_PTHREAD_MUTEX_DESTROY(&fAccessLock, "TPZDohrPrecondV2SubDataList::~TPZDohrPrecondV2SubDataList()");
222  }
223 
225  pthread_mutex_t fAccessLock;
226 
228  std::list<TPZDohrPrecondV2SubData<TVar, TSubStruct> > fWork;
229 
232  {
233  PZ_PTHREAD_MUTEX_LOCK(&fAccessLock, "TPZDohrPrecondV2SubDataList::AddItem()");
234  fWork.push_back(data);
235  PZ_PTHREAD_MUTEX_UNLOCK(&fAccessLock, "TPZDohrPrecondV2SubDataList::AddItem()");
236  }
239  {
241  PZ_PTHREAD_MUTEX_LOCK(&fAccessLock, "TPZDohrPrecondV2SubDataList::PopItem()");
242  if (fWork.size()) {
243  result = *fWork.begin();
244  fWork.pop_front();
245  }
246  PZ_PTHREAD_MUTEX_UNLOCK(&fAccessLock, "TPZDohrPrecondV2SubDataList::PopItem()");
247  return result;
248  }
249 
252 
254  static void *ThreadWork(void *voidptr);
255 
256 };
257 
260 #endif
TPZAutoPointer< TPZFMatrix< TVar > > fInput_local
Input matrix.
TPZDohrPrecondThreadV1Data(const TPZDohrPrecond< TVar, TSubStruct > *ptr, const TPZFMatrix< TVar > &input, TPZFMatrix< TVar > &output)
#define PZ_PTHREAD_MUTEX_UNLOCK(mutex, fn)
Definition: pz_pthread.h:79
Assembling using Dohrmann algorithm. Sub structure.
const TPZDohrPrecond< TVar, TSubStruct > * fDohrMatrix
Pointer to the dohr matrix.
Auxiliar structure for v2 vector to compute the preconditioner developed by Dohrmann. Sub Structure.
void ComputeV2(const TPZFMatrix< TVar > &x, TPZFMatrix< TVar > &v2) const
Compute the contribution of the sub domains.
Contains the TPZDohrSubstruct class which implements sub structure matrices using Dohrman algorithm...
pthread_mutex_t fAccessLock
Mutex which will enable the access protection of the list.
Auxiliar structure with list for v2 vector data. Sub Structure.
Defines a abstract class of solvers which will be used by matrix classes. Solver. ...
Definition: pzmatrix.h:32
Defines step solvers class. Solver.
Definition: pzmganalysis.h:17
TPZDohrPrecondV2SubDataList(TPZAutoPointer< TPZDohrAssembleList< TVar > > &assemble)
TPZAutoPointer< TSubStruct > fSubStructure
Pointer to the dohr matrix.
Contains the TPZDohrMatrix class which implements a matrix divided into substructures. Also contains the TPZDohrThreadMultData and TPZDohrThreadMultList structs.
virtual TPZMatrix< TVar > * Clone() const override
int64_t fNumCoarse
Size of the coarse system.
clarg::argString input("-if", "input file", "cube1.txt")
TPZDohrPrecondV2SubData(const TPZDohrPrecondV2SubData< TVar, TSubStruct > &copy)
Implements a matrix divided into substructures. Matrix Sub structure.
Definition: tpzdohrmatrix.h:30
TPZDohrPrecondV2SubData & operator=(const TPZDohrPrecondV2SubData &copy)
Contains the TPZDohrAssembly class which implements assembling using Dohrmann algorithm.
std::list< TPZAutoPointer< TSubStruct > > fGlobal
The matrix class is a placeholder for a list of substructures.
TPZDohrPrecondV2SubData(int subindex, const TPZAutoPointer< TSubStruct > &substruct, TPZAutoPointer< TPZFMatrix< TVar > > res_local)
Contains the TPZDohrAssembleItem and TPZDohrAssembleList structs to assembling using Dohrman algorith...
int fNumThreads
Number of threads used during preconditioning.
std::list< TPZAutoPointer< TSubStruct > > & Global()
The matrix class is a placeholder for a list of substructures.
TPZAutoPointer< TPZDohrAssembleItem< TVar > > fv2_local
The local contribution to the v2 vector.
void AddItem(TPZDohrPrecondV2SubData< TVar, TSubStruct > &data)
Interface to add items in a thread safe way.
Implements a matrix which computes the preconditioner developed by Dohrmann. Sub Structure.
TPZDohrPrecondV2SubData< TVar, TSubStruct > PopItem()
Interface to pop an item in a thread safe way.
#define PZ_PTHREAD_MUTEX_LOCK(mutex, fn)
Definition: pz_pthread.h:76
static void * ComputeV1(void *dataptr)
Compute the contribution of the coarse matrix.
int ClassId() const override
Routines to send and receive messages.
TPZAutoPointer< TPZDohrAssembleList< TVar > > fAssemblyStructure
The local contribution to the v2 vector.
To assembling one item using Dohrmann algorithm. Sub structure.
std::list< TPZDohrPrecondV2SubData< TVar, TSubStruct > > fWork
The list of structures which need to be computed.
List of items to assembling using Dohrmann algorithm.
Contains declaration of the TPZAutoPointer class which has Increment and Decrement actions are mutexe...
#define PZ_PTHREAD_MUTEX_DESTROY(mutex, fn)
Definition: pz_pthread.h:73
Full matrix class. Matrix.
Definition: pzfmatrix.h:32
const TPZFMatrix< TVar > * fInput
Input matrix.
virtual void MultAddTBB(const TPZFMatrix< TVar > &x, const TPZFMatrix< TVar > &y, TPZFMatrix< TVar > &z, const TVar alpha, const TVar beta, const int opt) const
int32_t Hash(std::string str)
Definition: TPZHash.cpp:10
void SetNumThreads(int numthreads)
Contains TPZMatrix<TVar>class, root matrix class.
virtual void MultAdd(const TPZFMatrix< TVar > &x, const TPZFMatrix< TVar > &y, TPZFMatrix< TVar > &z, const TVar alpha, const TVar beta, const int opt) const override
The only method any matrix class needs to implement.
Auxiliar structure with thread to compute the preconditioner developed by Dohrmann. Sub Structure.
TPZDohrPrecond()
Empty constructor for restoring.
void Read(TPZStream &buf, void *context) override
Unpacks the object structure from a stream of bytes.
void Write(TPZStream &buf, int withclassid) const override
Packs the object structure in a stream of bytes.
void ComputeV1(const TPZFMatrix< TVar > &x, TPZFMatrix< TVar > &v1) const
Compute the contribution of the coarse matrix.
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
#define PZ_PTHREAD_MUTEX_INIT(mutex, attr, fn)
Definition: pz_pthread.h:70
TPZStepSolver< TVar > * fCoarse
The global matrix associated with the coarse degrees of freedom.
void Initialize()
Initialize the necessary datastructures.
TPZAutoPointer< TPZDohrAssembly< TVar > > fAssemble
void SetSolver(TPZSolver< TVar > &solver)
Specify the solution process for the coarse matrix.
TPZFMatrix< TVar > * fOutput
Matrix where the coarse solution will be contributed.
int ClassId() const override
Define the class id associated with the class.
Definition: pzmatrix.h:957
Root matrix class (abstract). Matrix.
Definition: pzmatrix.h:60
This class implements a reference counter mechanism to administer a dynamically allocated object...