NeoPZ
TPZParFrontMatrix.cpp
Go to the documentation of this file.
1 
6 #include "TPZParFrontMatrix.h"
7 #include "TPZFrontMatrix.h"
8 #include "pzsfulmat.h"
9 #include "TPZFront.h"
10 #include "pzstack.h"
11 #include "pzreal.h"
12 #include <math.h>
13 #include "pz_pthread.h"
14 
15 #include "tpzeqnarray.h"
16 #include <iostream>
17 #include <stdio.h>
18 #include <stdlib.h>
19 #include <fstream>
20 
21 using namespace std;
22 
23 #include "pzlog.h"
24 
25 #ifdef LOG4CXX
26 
27 static LoggerPtr logger(Logger::getLogger("pz.strmatrix.frontstructmatrix"));
28 static LoggerPtr loggerfw(Logger::getLogger("pz.frontal.frontmatrix.fw"));
29 
30 #endif
31 
33 pthread_mutex_t mutex_write = PTHREAD_MUTEX_INITIALIZER;
35 pthread_cond_t conda_write = PTHREAD_COND_INITIALIZER;
36 
37 // At the class constructor creates a thread
38 // this thread will be active while ParFrontMatrix is active
39 // It will check if a stack contains some equations to be writen to the disk
40 
41 template<class TVar, class store, class front>
43 TPZRegisterClassId(&TPZParFrontMatrix::ClassId), fFinish(0)
44 {
45  fEqnStack.Resize(0);
46  pthread_mutex_t mlocal = PTHREAD_MUTEX_INITIALIZER;
47  fwritelock = mlocal;
48  pthread_cond_t clocal = PTHREAD_COND_INITIALIZER;
49  fwritecond = clocal;
50  /* fFront.Reset();
51  fStorage.Reset();
52  fNumElConnected.Resize(0);
53  fLastDecomposed = -1;
54  fNumEq=0;
55  */
56 }
57 
58 template<class TVar, class store, class front>
60 TPZRegisterClassId(&TPZParFrontMatrix::ClassId), TPZFrontMatrix<TVar, store, front>(globalsize),
61 fFinish(0)
62 {
63  fEqnStack.Resize(0);
64  pthread_mutex_t mlocal = PTHREAD_MUTEX_INITIALIZER;
65  fwritelock = mlocal;
66  pthread_cond_t clocal = PTHREAD_COND_INITIALIZER;
67  fwritecond = clocal;
68 }
69 
70 template<class TVar, class store, class front>
72 }
73 
74 template<class TVar, class store, class front>
76 {
77 
78  // message #1.3 to fFront:TPZFront
79  this->fFront.AddKel(elmat, destinationindex);
80 #ifdef LOG4CXX
81  {
82  std::stringstream sout;
83  sout << "Frondwidth after AddKel "<< this->fFront.FrontSize();
84  LOGPZ_INFO(loggerfw,sout.str())
85  }
86 #endif
87 
88  int64_t mineq, maxeq;
89  this->EquationsToDecompose(destinationindex, mineq, maxeq);
90  if(maxeq >= mineq) {
91 
93 
94  this->fFront.DecomposeEquations(mineq,maxeq,*AuxEqn);
95  this->CheckCompress();
96  PZ_PTHREAD_MUTEX_LOCK(&fwritelock,"TPZParFrontMatrix<...>::AddKel()");
97  fEqnStack.Push(AuxEqn);
98  if(maxeq == this->Rows()-1){
99  cout << "Decomposition finished" << endl;
100  cout.flush();
101  FinishWriting();
102  //fStorage.ReOpen();
103  }
104  PZ_PTHREAD_MUTEX_UNLOCK(&fwritelock,"TPZParFrontMatrix<...>::AddKel()");
105  PZ_PTHREAD_COND_SIGNAL(&fwritecond,"TPZParFrontMatrix<...>::AddKel()");
106  }
107  this->fDecomposed = this->fFront.GetDecomposeType();
108 }
109 template<class TVar, class store, class front>
111 {
112  this->fFront.AddKel(elmat, sourceindex, destinationindex);
113 #ifdef LOG4CXX
114  if (loggerfw->isDebugEnabled())
115  {
116  std::stringstream sout;
117  sout << "Frondwidth after AddKel "<< this->fFront.FrontSize();
118  LOGPZ_DEBUG(loggerfw,sout.str())
119  }
120 #endif
121  int64_t mineq, maxeq;
122  this->EquationsToDecompose(destinationindex, mineq, maxeq);
123 
124  if(maxeq >= mineq) {
125  TPZEqnArray<TVar> *AuxEqn = new TPZEqnArray<TVar>;
126 
127  this->fFront.DecomposeEquations(mineq,maxeq,*AuxEqn);
128  this->CheckCompress();
129  PZ_PTHREAD_MUTEX_LOCK(&fwritelock,"TPZParFrontMatrix<...>::AddKel()");
130  fEqnStack.Push(AuxEqn);
131  if(maxeq == this->Rows()-1){
132  //check if writeing is over and closes file
133  cout << endl << "Decomposition finished" << endl;
134  cout.flush();
135  FinishWriting();
136  this->fFront.Reset(0);
137  //fStorage.ReOpen();
138  }
139  PZ_PTHREAD_MUTEX_UNLOCK(&fwritelock,"TPZParFrontMatrix<...>::AddKel()");
140  PZ_PTHREAD_COND_SIGNAL(&fwritecond,"TPZParFrontMatrix<...>::AddKel()");
141  }
142  this->fDecomposed = this->fFront.GetDecomposeType();
143 }
144 
145 template<class TVar, class store, class front>
147  // FinishWriting already has a lock
148  //PZ_PTHREAD_MUTEX_LOCK(&fwritelock);
149  cout << endl << "FinishWriting" << endl;
150  cout.flush();
151  fFinish = 1;
152  // FinishWriting already has a lock
153  //pthread_mutex_unlock(&fwritelock);
154  PZ_PTHREAD_COND_SIGNAL(&fwritecond,"TPZParFrontMatrix<...>::FinishWriting()");
155 }
156 
157 template<class TVar, class store, class front>
160 #ifdef LOG4CXX
161  if (logger->isDebugEnabled())
162  {
163  std::stringstream sout;
164  sout << "Entering WriteFile thread execution";
165  LOGPZ_DEBUG(logger,sout.str())
166  }
167 #endif
168  cout << endl << "Entering Decomposition" << endl;
169  cout.flush();
170  while(1){
172  PZ_PTHREAD_MUTEX_LOCK(&parfront->fwritelock,"TPZParFrontMatrix<...>::WriteFile()");
173 #ifdef LOG4CXX
174  if (logger->isDebugEnabled())
175  {
176  std::stringstream sout;
177  sout << "Acquired writelock";
178  LOGPZ_DEBUG(logger,sout.str())
179  }
180 #endif
181  if(parfront->fEqnStack.NElements() == 0){
182  if(parfront->fFinish == 1) {
183 #ifdef LOG4CXX
184  if (logger->isDebugEnabled())
185  {
186  std::stringstream sout;
187  sout << "Terminating WriteFile thread execution";
188  LOGPZ_DEBUG(logger,sout.str())
189  }
190 #endif
191  cout << "Leaving WHILE" << endl;
192  cout.flush();
193  break;
194  }
195 #ifdef LOG4CXX
196  if (logger->isDebugEnabled())
197  {
198  std::stringstream sout;
199  sout << "Entering cond_wait on fwritecond variable";
200  LOGPZ_DEBUG(logger,sout.str())
201  }
202 #endif
203  PZ_PTHREAD_COND_WAIT(&parfront->fwritecond, &parfront->fwritelock,"TPZParFrontMatrix<...>::WriteFile()");
204  }
205 
206  local = parfront->fEqnStack;
207  parfront->fEqnStack.Resize(0);
208 #ifdef LOG4CXX
209  if (logger->isDebugEnabled())
210  {
211  std::stringstream sout;
212  sout << "Copied the equation stack releasing the writelock";
213  LOGPZ_DEBUG(logger,sout.str())
214  }
215 #endif
216 
217  PZ_PTHREAD_MUTEX_UNLOCK(&parfront->fwritelock,"TPZParFrontMatrix<...>::WriteFile()");
218  int64_t neqn = local.NElements();
219 
220  int64_t eq;
221  for(eq=0; eq<neqn; eq++) {
222  parfront->fStorage.AddEqnArray(local[eq]);
223  delete local[eq];
224  }
225  }
226  parfront->fStorage.FinishWriting();
227  parfront->fStorage.ReOpen();
228  parfront->fFinish = 0;
229 #ifdef LOG4CXX
230  if (logger->isDebugEnabled())
231  {
232  std::stringstream sout;
233  sout << "Releasing writelock";
234  LOGPZ_DEBUG(logger,sout.str())
235  }
236 #endif
237 
238  PZ_PTHREAD_MUTEX_UNLOCK(&parfront->fwritelock,"TPZParFrontMatrix<...>::WriteFile()");
239 #ifdef LOG4CXX
240  if (logger->isDebugEnabled())
241  {
242  std::stringstream sout;
243  sout << "Falling through on the write thread";
244  LOGPZ_DEBUG(logger,sout.str())
245  }
246 #endif
247  std::cout << "Terminating write thread\n";
248  return (0);
249 }
250 
251 template<class TVar>
252 class TPZStackEqnStorage;
253 template<class TVar>
254 class TPZFileEqnStorage;
255 template<class TVar>
256 class TPZFrontSym;
257 template<class TVar>
258 class TPZFrontNonSym;
259 
260 
262 template class TPZParFrontMatrix<float, TPZFileEqnStorage<float>, TPZFrontSym<float> >;
264 template class TPZParFrontMatrix<float, TPZFileEqnStorage<float>, TPZFrontNonSym<float> >;
265 
267 template class TPZParFrontMatrix<double, TPZFileEqnStorage<double>, TPZFrontSym<double> >;
269 template class TPZParFrontMatrix<double, TPZFileEqnStorage<double>, TPZFrontNonSym<double> >;
270 
272 template class TPZParFrontMatrix<long double, TPZFileEqnStorage<long double>, TPZFrontSym<long double> >;
274 template class TPZParFrontMatrix<long double, TPZFileEqnStorage<long double>, TPZFrontNonSym<long double> >;
275 
277 template class TPZParFrontMatrix<std::complex<float>, TPZFileEqnStorage<std::complex<float> >, TPZFrontSym<std::complex<float> > >;
278 template class TPZParFrontMatrix<std::complex<float>, TPZStackEqnStorage<std::complex<float> >, TPZFrontNonSym<std::complex<float> > >;
279 template class TPZParFrontMatrix<std::complex<float>, TPZFileEqnStorage<std::complex<float> >, TPZFrontNonSym<std::complex<float> > >;
280 
282 template class TPZParFrontMatrix<std::complex<double>, TPZFileEqnStorage<std::complex<double> >, TPZFrontSym<std::complex<double> > >;
283 template class TPZParFrontMatrix<std::complex<double>, TPZStackEqnStorage<std::complex<double> >, TPZFrontNonSym<std::complex<double> > >;
284 template class TPZParFrontMatrix<std::complex<double>, TPZFileEqnStorage<std::complex<double> >, TPZFrontNonSym<std::complex<double> > >;
285 
287 template class TPZParFrontMatrix<std::complex<long double>, TPZFileEqnStorage<std::complex<long double> >, TPZFrontSym<std::complex<long double> > >;
288 template class TPZParFrontMatrix<std::complex<long double>, TPZStackEqnStorage<std::complex<long double> >, TPZFrontNonSym<std::complex<long double> > >;
289 template class TPZParFrontMatrix<std::complex<long double>, TPZFileEqnStorage<std::complex<long double> >, TPZFrontNonSym<std::complex<long double> > >;
290 
#define PZ_PTHREAD_MUTEX_UNLOCK(mutex, fn)
Definition: pz_pthread.h:79
Has the same porpouse of EqnStack but stores the EqnArrays in a different form (binary files)...
TPZParFrontMatrix()
Simple Constructor.
#define PZ_PTHREAD_COND_SIGNAL(cond, fn)
Definition: pz_pthread.h:130
void CheckCompress()
Checks if FrontMatrix needs a compression, if so calls Compress method.
Contains definitions to LOGPZ_DEBUG, LOGPZ_INFO, LOGPZ_WARN, LOGPZ_ERROR and LOGPZ_FATAL, and the implementation of the inline InitializePZLOG(string) function using log4cxx library or not. It must to be called out of "#ifdef LOG4CXX" scope.
FrontMatrix with parallel techniques included. Frontal.
Responsible for the frontal method as a whole. Frontal.
pthread_mutex_t mutex_write
Initializing semaphore.
int fFinish
Boolean responsibility. Assumes values 0 and 1.
It is an equation array, generally in its decomposed form. Frontal.
Definition: tpzeqnarray.h:36
Contains the TPZFront class which implements decomposition process of the frontal matrix...
Responsible for storing arrays of equations (mostly in a decomposed form). Frontal.
store fStorage
Indicates storage schema. Assumes values TPZFileEqnStorage for binary file and TPZStackEqnStorage for...
static void * WriteFile(void *t)
Used in an independent thread to write decomposed equations to a binary file.
char fDecomposed
Decomposition type used to decompose the current matrix.
Definition: pzmatrix.h:783
Contains the TPZParFrontMatrix class which implements FrontMatrix with parallel techniques.
Abstract class implements storage and decomposition process of the frontal matrix involving simmetry ...
#define PZ_PTHREAD_COND_WAIT(cond, mutex, fn)
Definition: pz_pthread.h:127
TPZStack< TPZEqnArray< TVar > * > fEqnStack
Buffer of pointers to decomposed equations. Stored in a Stack.
pthread_cond_t conda_write
Initializing condition.
virtual void Resize(const int64_t newsize, const T &object)
Resizes the vector object.
Definition: pzmanvector.h:426
Contains the TPZEqnArray class which implements an equation array.
void Push(const T object)
Pushes a copy of the object on the stack.
Definition: pzstack.h:80
#define LOGPZ_INFO(A, B)
Define log for informations.
Definition: pzlog.h:89
#define PZ_PTHREAD_MUTEX_LOCK(mutex, fn)
Definition: pz_pthread.h:76
pthread_cond_t fwritecond
Condition variable used in management of writeing to disk and decomposition.
#define LOGPZ_DEBUG(A, B)
Define log for debug info.
Definition: pzlog.h:87
int64_t Rows() const
Returns number of rows.
Definition: pzmatrix.h:803
front fFront
Indicates Front matrix type. Assumes values TPZFrontSym for symmetric front and TPZFrontNonSym for no...
Full matrix class. Matrix.
Definition: pzfmatrix.h:32
virtual void AddKel(TPZFMatrix< TVar > &elmat, TPZVec< int64_t > &sourceindex, TPZVec< int64_t > &destinationindex) override
Add a contribution of a stiffness matrix.
A simple stack.
~TPZParFrontMatrix()
Simple Destructor.
void EquationsToDecompose(TPZVec< int64_t > &destinationindex, int64_t &lower_eq, int64_t &upper_eq)
Sends a message to decompose equations from lower_eq to upper_eq, according to destination index...
This class implements a stack object. Utility.
Definition: pzcheckmesh.h:14
Contains the TPZAbstractFrontMatrix class which implements a matrix stored in a frontal decomposition...
int64_t NElements() const
Returns the number of elements of the vector.
Definition: pzvec.h:190
Contains the declaration of TPZFlopCounter class and TPZCounter struct.
void FinishWriting()
Sets the flag fFinish to its true value.
pthread_mutex_t fwritelock
Mutual exclusion locks used in management of writeing to disk and decomposition.
Abstract class implements storage and decomposition process of the frontal matrix involving non-simet...
Contains TPZSFMatrix class which implements a symmetric full matrix.
int ClassId() const override
Define the class id associated with the class.