NeoPZ
tpzdohrmatrix.h
Go to the documentation of this file.
1 
9 #ifndef TPZDOHRMATRIX_H
10 #define TPZDOHRMATRIX_H
11 
12 #include "pzmatrix.h"
13 #include <list>
14 #include <sstream>
15 #include "tpzautopointer.h"
16 #include "tpzdohrsubstruct.h"
18 #include "tpzdohrassembly.h"
19 
20 #include "tpzdohrassemblelist.h"
21 
22 #include "pz_pthread.h"
23 
29 template <class TVar, class TSubStruct>
30 class TPZDohrMatrix : public TPZMatrix<TVar>
31 {
32 public:
34  typedef typename std::list<TPZAutoPointer<TSubStruct> > SubsList;
35 private:
36  SubsList fGlobal;
37 
38  int fNumCoarse; //n(c)
39 
42 
43 public:
44 
46 
48 
50  {
51 
52  }
53 
54  TPZDohrMatrix(const TPZDohrMatrix &cp) : fGlobal(cp.fGlobal), fNumCoarse(cp.fNumCoarse), fNumThreads(cp.fNumThreads),
55  fAssembly(cp.fAssembly)
56  {
57  }
58 
59  // CLONEDEF(TPZDohrMatrix)
60  virtual TPZMatrix<TVar>*Clone() const override { return new TPZDohrMatrix(*this); }
61 
63 
64  const SubsList &SubStructures() const
65  {
66  return fGlobal;
67  }
68 
69  int NumCoarse() const
70  {
71  return fNumCoarse;
72  }
73 
74  int NumThreads() const
75  {
76  return fNumThreads;
77  }
78 
80  {
81  fNumThreads = numthreads;
82  }
83 
86  return (*fGlobal.begin());
87  }
89  void Print(const char *name, std::ostream& out,const MatrixOutputFormat form = EFormatted) const override
90  {
91  out << __PRETTY_FUNCTION__ << std::endl;
92  out << name << std::endl;
93  out << "Number of coarse equations " << fNumCoarse << std::endl;
94  typename SubsList::const_iterator iter;
95  for (iter=fGlobal.begin();iter!=fGlobal.end();iter++) {
96  (*iter)->Print(out);
97  }
98  }
99 
100  void SetNumCornerEqs(int nc)
101  {
102  fNumCoarse = nc;
103  }
105  void Initialize();
108  {
109  fGlobal.push_back(substruct);
110  }
111 
123  virtual void MultAddTBB(const TPZFMatrix<TVar> &x,const TPZFMatrix<TVar> &y, TPZFMatrix<TVar> &z,const TVar alpha,const TVar beta,const int opt) const;
124 
126  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;
127 
130 
132  void AddInternalSolution(TPZFMatrix<TVar> &solution);
133 
144  void Read(TPZStream &buf, void *context) override;
150  void Write(TPZStream &buf, int withclassid) const override;
151 
153  public:
154 int ClassId() const override;
155 
156 };
157 
158 template <class TVar, class TSubStruct>
160  return Hash("TPZDohrMatrix") ^ TPZMatrix<TVar>::ClassId() << 1 ^ TSubStruct().ClassId() << 2;
161 }
162 
167 template <class TSubStruct>
169 {
171  TPZDohrThreadMultData() : fisub(-1), fSub(0)
172  {
173  }
174  TPZDohrThreadMultData(int isub, TPZAutoPointer<TSubStruct> submesh) : fisub(isub), fSub(submesh)
175  {
176  }
178  TPZDohrThreadMultData(const TPZDohrThreadMultData<TSubStruct> &cp) : fisub(cp.fisub), fSub(cp.fSub)
179  {
180  }
183  {
184  fisub = cp.fisub;
185  fSub = cp.fSub;
186  return *this;
187  }
188  int fisub;
190 
191  bool IsValid()
192  {
193  return (fisub >= 0);
194  }
195 };
196 
201 template <class TVar, class TSubStruct>
203 {
207  TVar fAlpha;
209  pthread_mutex_t fAccessLock;
213  std::list<TPZDohrThreadMultData<TSubStruct> > fWork;
216 
217  TPZDohrThreadMultList(const TPZFMatrix<TVar> &x, TVar alpha, TPZAutoPointer<TPZDohrAssembly<TVar> > assembly, TPZAutoPointer<TPZDohrAssembleList<TVar> > &assemblestruct) : fInput(&x), fAlpha(alpha),
218  fAssembly(assembly), fAssemblyStructure(assemblestruct)
219  {
220  PZ_PTHREAD_MUTEX_INIT(&fAccessLock, 0, "TPZDohrThreadMultList::TPZDohrThreadMultList(...)");
221  }
223  {
224  PZ_PTHREAD_MUTEX_DESTROY(&fAccessLock, "TPZDohrThreadMultList::TPZDohrThreadMultList()");
225  }
226 
228  static void *ThreadWork(void *voidptr);
231  {
232  PZ_PTHREAD_MUTEX_LOCK(&fAccessLock, "TPZDohrThreadMultList::AddItem()");
233  fWork.push_back(data);
234  PZ_PTHREAD_MUTEX_UNLOCK(&fAccessLock, "TPZDohrThreadMultList::AddItem()");
235  }
238  {
240  PZ_PTHREAD_MUTEX_LOCK(&fAccessLock, "TPZDohrThreadMultList::PopItem()");
241  if (fWork.size()) {
242  result = *fWork.begin();
243  fWork.pop_front();
244  }
245  PZ_PTHREAD_MUTEX_UNLOCK(&fAccessLock, "TPZDohrThreadMultList::PopItem()");
246  return result;
247  }
248 };
249 
250 #endif
void AdjustResidual(TPZFMatrix< TVar > &res)
Adjust the residual to zero the residual of the internal connects.
#define PZ_PTHREAD_MUTEX_UNLOCK(mutex, fn)
Definition: pz_pthread.h:79
const TPZFMatrix< TVar > * fInput
The vector with which we will multiply.
Assembling using Dohrmann algorithm. Sub structure.
virtual TPZMatrix< TVar > * Clone() const override
Definition: tpzdohrmatrix.h:60
std::list< TPZDohrThreadMultData< TSubStruct > > fWork
The list of data objects which need to treated by the threads.
Contains the TPZDohrSubstruct class which implements sub structure matrices using Dohrman algorithm...
MatrixOutputFormat
Defines output format.
Definition: pzmatrix.h:55
void Read(TPZStream &buf, void *context) override
Unpacks the object structure from a stream of bytes.
TPZAutoPointer< TPZDohrAssembleList< TVar > > fAssemblyStructure
The local contribution to the v2 vector.
void SetNumCornerEqs(int nc)
void AddSubstruct(TPZAutoPointer< TSubStruct > substruct)
It adds a substruct.
TPZDohrThreadMultList(const TPZFMatrix< TVar > &x, TVar alpha, TPZAutoPointer< TPZDohrAssembly< TVar > > assembly, TPZAutoPointer< TPZDohrAssembleList< TVar > > &assemblestruct)
TPZDohrThreadMultData< TSubStruct > & operator=(const TPZDohrThreadMultData< TSubStruct > &cp)
Implement the attribution operator.
Implements a matrix divided into substructures. Matrix Sub structure.
Definition: tpzdohrmatrix.h:30
Contains the TPZDohrAssembly class which implements assembling using Dohrmann algorithm.
virtual void MultAddTBB(const TPZFMatrix< TVar > &x, const TPZFMatrix< TVar > &y, TPZFMatrix< TVar > &z, const TVar alpha, const TVar beta, const int opt) const
It computes z = beta * y + alpha * opt(this)*x but z and x can not overlap in memory.
.. . Sub structure
Contains the TPZDohrAssembleItem and TPZDohrAssembleList structs to assembling using Dohrman algorith...
TPZAutoPointer< TSubStruct > GetFirstSub()
Just a method for tests.
Definition: tpzdohrmatrix.h:85
TPZDohrThreadMultData(int isub, TPZAutoPointer< TSubStruct > submesh)
TVar fAlpha
Scalar multiplication factor.
#define PZ_PTHREAD_MUTEX_LOCK(mutex, fn)
Definition: pz_pthread.h:76
void Print(const char *name, std::ostream &out, const MatrixOutputFormat form=EFormatted) const override
Just a method for tests.
Definition: tpzdohrmatrix.h:89
TPZDohrThreadMultData(const TPZDohrThreadMultData< TSubStruct > &cp)
Copy constructor.
int fNumThreads
Number of threads that will be used during the matrix vector multiplication.
Definition: tpzdohrmatrix.h:41
void AddItem(TPZDohrThreadMultData< TSubStruct > &data)
Interface to add items in a thread safe way.
string res
Definition: test.py:151
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
Contains the TPZDohrSubstructCondense class which condenses matrix divided in sub structures...
void SetNumThreads(int numthreads)
Definition: tpzdohrmatrix.h:79
TPZAutoPointer< TPZDohrAssembly< TVar > > fAssembly
The data structure which defines the assemble destinations.
void AddInternalSolution(TPZFMatrix< TVar > &solution)
Add the solution corresponding to the internal residual.
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
Full matrix class. Matrix.
Definition: pzfmatrix.h:32
TPZDohrMatrix(const TPZDohrMatrix &cp)
Definition: tpzdohrmatrix.h:54
int32_t Hash(std::string str)
Definition: TPZHash.cpp:10
Contains TPZMatrix<TVar>class, root matrix class.
TPZDohrThreadMultData< TSubStruct > PopItem()
Interface to pop an item in a thread safe way.
TPZAutoPointer< TSubStruct > fSub
int NumThreads() const
Definition: tpzdohrmatrix.h:74
TPZDohrThreadMultData()
Default constructor.
std::list< TPZAutoPointer< TSubStruct > > SubsList
The matrix class is a placeholder for a list of substructures.
Definition: tpzdohrmatrix.h:34
const SubsList & SubStructures() const
Definition: tpzdohrmatrix.h:64
void Write(TPZStream &buf, int withclassid) const override
Packs the object structure in a stream of bytes.
int ClassId() const override
Routines to send and receive messages.
.. . Sub structure
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
void Initialize()
Initialize the necessary datastructures.
SubsList fGlobal
Definition: tpzdohrmatrix.h:36
#define PZ_PTHREAD_MUTEX_INIT(mutex, attr, fn)
Definition: pz_pthread.h:70
TPZAutoPointer< TPZDohrAssembly< TVar > > fAssembly
Definition: tpzdohrmatrix.h:45
pthread_mutex_t fAccessLock
Mutex which will enable the access protection of the list.
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...
int NumCoarse() const
Definition: tpzdohrmatrix.h:69