NeoPZ
pzanalysis.h
Go to the documentation of this file.
1 
6 #ifndef ANALYSISH
7 #define ANALYSISH
8 
9 #include <iostream> // for string, cout, ostream
10 #include <set> // for set
11 #include "TPZGuiInterface.h" // for TPZGuiInterface
12 #include "pzerror.h" // for DebugStop
13 #include "pzmatrix.h" // for TPZFMatrix, TPZMatrix
14 #include "pzreal.h" // for STATE, REAL
15 #include "TPZRenumbering.h" // for TPZRenumbering
16 #include "pzstrmatrix.h" // for TPZStructMatrix
17 #include "pzvec.h" // for TPZVec
18 #include "tpzautopointer.h" // for TPZAutoPointer
19 class TPZCompEl;
20 class TPZCompMesh;
21 class TPZConnect;
22 class TPZGeoMesh;
23 class TPZGraphMesh;
24 template <class TVar> class TPZMatrixSolver;
25 
32 class TPZAnalysis : public TPZSavable {
33 
34 public:
35 
38 
39 
40 protected:
60  int fStep;
62  REAL fTime;
63 
66 
69 
72 
75 
77  class TTablePostProcess : public TPZSavable {
78  public :
86 
87  int ClassId() const override;
88 
89  void Write(TPZStream &buf, int withclassid) const override;
90 
91  void Read(TPZStream &buf, void *context) override;
92  };
93 
95 
96  public :
97 
99  std::function<void (const TPZVec<REAL> &loc, TPZVec<STATE> &result, TPZFMatrix<STATE> &deriv)> fExact;
100 
102  TPZAnalysis(TPZCompMesh *mesh, bool mustOptimizeBandwidth = true, std::ostream &out = std::cout);
103 
105  TPZAnalysis(TPZAutoPointer<TPZCompMesh> mesh, bool mustOptimizeBandwidth = true, std::ostream &out = std::cout);
106 
109  fGuiInterface = gui;
110  }
111 
114  return fGuiInterface;
115  }
116 
118  bool AmIKilled(){
119  if(fGuiInterface){
120  return fGuiInterface->AmIKilled();
121  }
122  else return false;
123  }
124 
126  virtual void SetCompMesh(TPZCompMesh * mesh, bool mustOptimizeBandwidth);
127 
129  TPZAnalysis();
130 
131  void Write(TPZStream &buf, int withclassid) const override;
132 
133  void Read(TPZStream &buf, void *context) override;
134 
136  virtual ~TPZAnalysis(void);
137 
139  void CleanUp();
140 
143  {
144  fRenumber = renumber;
145  }
146 
148  void OptimizeBandwidth();
149 
151  int HighestDimension();
152 
154  void Resequence(int firstel = -1);
155 
162 
163 
165  virtual void Assemble();
166 
168  virtual void AssembleResidual();
169 
171  virtual void Solve();
172 
174  TPZFMatrix<STATE> &Rhs() { return fRhs;}
175 
178 
180  TPZCompMesh *Mesh()const { return fCompMesh;}
183  if(!fStructMatrix)
184  {
185  DebugStop();
186  }
187  return fStructMatrix;
188  }
189 
192  TPZMatrixSolver<STATE> *BuildPreconditioner(EPrecond preconditioner, bool overlap);
193 
195  void SetStep(int step)
196  {
197  fStep = step;
198  }
199 
201  {
202  fNthreadsError = nthreads;
203  }
204 
205  int GetStep()
206  {
207  return fStep;
208  }
209 
211  void SetTime(REAL time);
213  REAL GetTime();
214 
215 private:
216 
218  TPZMatrixSolver<STATE> *BuildSequenceSolver(TPZVec<int64_t> &graph, TPZVec<int64_t> &graphindex, int64_t neq, int numcolors, TPZVec<int> &colors);
219 
220 public:
222  void ShowShape(const std::string &plotfile, TPZVec<int64_t> &equationindices);
224  void ShowShape(const std::string &plotfile, TPZVec<int64_t> &equationindices, int matid, const TPZVec<std::string> &varname);
226  void LoadShape(double dx,double dy, int64_t numelem,TPZConnect* nod);
227 
229  virtual void Run(std::ostream &out = std::cout);
231  virtual void DefineGraphMesh(int dimension, const TPZVec<std::string> &scalnames, const TPZVec<std::string> &vecnames, const std::string &plotfile);
233  virtual void DefineGraphMesh(int dimension, const TPZVec<std::string> &scalnames, const TPZVec<std::string> &vecnames, const TPZVec<std::string> &tensnames, const std::string &plotfile);
235  virtual void DefineGraphMesh(int dimension, const std::set<int> & matids ,const TPZVec<std::string> &scalnames, const TPZVec<std::string> &vecnames, const std::string &plotfile);
237  virtual void DefineGraphMesh(int dimension, const std::set<int> & matids, const TPZVec<std::string> &scalnames, const TPZVec<std::string> &vecnames, const TPZVec<std::string> &tensnames, const std::string &plotfile);
238 
240  virtual void CloseGraphMesh();
241 
243  TPZGraphMesh *GraphMesh(int dimension) {
244  return fGraphMesh[dimension-1];
245  }
247  virtual void PostProcess(int resolution);
249  virtual void PostProcess(int resolution, int dimension);
250 
252  void IdentifyPostProcessingMatIds(int dimension, std::set<int> & matids);
253 
254 
261  virtual void DefineElementTable(int dimension, TPZVec<int64_t> &GeoElIds, TPZVec<REAL> &points);
263  virtual void SetTableVariableNames(TPZVec<std::string> varnames);
265  virtual void PrePostProcessTable(std::ostream &out_file);
267  virtual void PostProcessTable(std::ostream &out_file);
268 
269 
271  void PostProcessTable( TPZFMatrix<REAL> &pos,std::ostream &out= std::cout );
272 
276  virtual void LoadSolution();
278  virtual void LoadSolution(const TPZFMatrix<STATE> &sol){
279  fSolution = sol;
280  this->LoadSolution();
281  }
282 
284  TPZVec<STATE> Integrate(const std::string &varname, const std::set<int> &matids);
285 
287  void SetExact(std::function<void (const TPZVec<REAL> &loc, TPZVec<STATE> &result, TPZFMatrix<STATE> &deriv)> f);
289  virtual void PostProcess(TPZVec<REAL> &loc, std::ostream &out = std::cout);
290 
295  virtual void PostProcessError(TPZVec<REAL> &, bool store_error = true, std::ostream &out = std::cout);
296 
297  virtual void PostProcessErrorSerial(TPZVec<REAL> &, bool store_error = true, std::ostream &out = std::cout);
298 
299  virtual void PostProcessErrorParallel(TPZVec<REAL> &, bool store_error = true, std::ostream &out = std::cout);
300 
302 
304  void Print( const std::string &name , std::ostream &out );
305 
307  void PrintVectorByElement(std::ostream &out, TPZFMatrix<STATE> &vec, REAL tol = 1.e-10);
308 
312  void AnimateRun(int64_t num_iter, int steps,
313  TPZVec<std::string> &scalnames, TPZVec<std::string> &vecnames, const std::string &plotfile);
315  void SetSolver(TPZMatrixSolver<STATE> &solver);
319  void SetStructuralMatrix(TPZStructMatrix &strmatrix);
320 
321  public:
322 int ClassId() const override;
323 
324  struct ThreadData{
325 
327 
328  ThreadData(TPZAdmChunkVector<TPZCompEl *> &elvec, bool store_error, std::function<void (const TPZVec<REAL> &loc, TPZVec<STATE> &result, TPZFMatrix<STATE> &deriv)> f);
329 
330  ~ThreadData();
331 
332  static void *ThreadWork(void *threaddata);
333 
334  std::function<void (const TPZVec<REAL> &loc, TPZVec<STATE> &result, TPZFMatrix<STATE> &deriv)> fExact;
335 
336  int64_t fNextElement;
337 
338  int ftid;
339 
340  bool fStoreError = false;
341 
342  // Vector with errors. Assuming no more than a 100 threads
344 
346  pthread_mutex_t fAccessElement;
347 
349  pthread_mutex_t fGetUniqueId;
350 
351  };
352 
353  friend struct ThreadData;
354 
355 };
356 
357 
358 inline void
359 
360 TPZAnalysis::SetExact(std::function<void (const TPZVec<REAL> &loc, TPZVec<STATE> &result, TPZFMatrix<STATE> &deriv)> f)
361 {
362  fExact=f;
363 }
364 
365 inline TPZMatrixSolver<STATE> &
366 
368  return (*fSolver);
369 }
370 
371 inline void TPZAnalysis::SetTime(REAL time){
372  this->fTime = time;
373 }
374 
375 inline REAL TPZAnalysis::GetTime(){
376  return this->fTime;
377 }
378 
379 #endif
pthread_mutex_t fGetUniqueId
Mutexes (to sum error)
Definition: pzanalysis.h:349
TPZVec< TPZCompEl * > fCompElPtr
Definition: pzanalysis.h:80
TPZMatrixSolver< STATE > * fSolver
Type of solver to be applied.
Definition: pzanalysis.h:52
void LoadShape(double dx, double dy, int64_t numelem, TPZConnect *nod)
Make assembling and clean the load and solution vectors.
Definition: pzanalysis.cpp:908
Represents a graphical mesh used for post processing purposes. Post processing.
Definition: pzgraphmesh.h:34
REAL GetTime()
Gets time used in dx files.
Definition: pzanalysis.h:375
Represents a set of shape functions associated with a computational element/side. Computational Eleme...
Definition: pzconnect.h:30
virtual void Solve()
Invert the stiffness matrix.
Definition: pzanalysis.cpp:352
void SetStructuralMatrix(TPZAutoPointer< TPZStructMatrix > strmatrix)
Set structural matrix as auto pointer for analysis.
Definition: pzanalysis.cpp:82
Implements a vector class which allows to use external storage provided by the user. Utility.
Definition: pzquad.h:16
void Print(const std::string &name, std::ostream &out)
Print connect and solution information.
Definition: pzanalysis.cpp:447
void OptimizeBandwidth()
Sets the computer connection block number from the graphical connections block number otimization...
Definition: pzanalysis.cpp:195
TPZFMatrix< STATE > fSolution
Solution vector.
Definition: pzanalysis.h:50
virtual void PostProcessTable(std::ostream &out_file)
Print the solution related with the computational element vector in post process. ...
clarg::argInt dimension("-d", "Matrices dimension M x M", 1000)
void Write(TPZStream &buf, int withclassid) const override
Writes this object to the TPZStream buffer. Include the classid if withclassid = true.
int fNthreadsError
Number of threads to be used for post-processing error.
Definition: pzanalysis.h:65
Templated vector implementation.
Contains declaration of TPZGuiInterface class.
void SetSolver(TPZMatrixSolver< STATE > &solver)
Set solver matrix.
Contains the TPZRenumbering class which defines the behavior to implementing node sequence numbering ...
TPZCompMesh * fCompMesh
Computational mesh.
Definition: pzanalysis.h:44
Defines PZError.
Defines a class of matrix solvers. Solver.
Definition: pzanalysis.h:24
virtual void Assemble()
Assemble the stiffness matrix and load vector.
Definition: pzanalysis.cpp:304
void PrintVectorByElement(std::ostream &out, TPZFMatrix< STATE > &vec, REAL tol=1.e-10)
Print the residual vector for those elements with entry above a given tolerance.
TPZCompMesh * Mesh() const
Returns the pointer to the computational mesh.
Definition: pzanalysis.h:180
TPZVec< STATE > Integrate(const std::string &varname, const std::set< int > &matids)
Integrate the postprocessed variable name over the elements included in the set matids.
virtual void PrePostProcessTable(std::ostream &out_file)
Prepare data to print post processing and print coordinates.
virtual void DefineGraphMesh(int dimension, const TPZVec< std::string > &scalnames, const TPZVec< std::string > &vecnames, const std::string &plotfile)
Define GrapMesh as V3D, DX, MV or VTK depending on extension of the file.
Definition: pzanalysis.cpp:952
Implements a chunk vector with free store administration. Utility.
Definition: TPZStream.h:39
TPZVec< std::string > fVariableNames
Definition: pzanalysis.h:83
Contains the TPZStructMatrixOR class which responsible for a interface among Matrix and Finite Elemen...
void SetGuiInterface(TPZAutoPointer< TPZGuiInterface > gui)
Defines gui interface object.
Definition: pzanalysis.h:108
EPrecond
Preconditioners which can be created by objects of this class.
Definition: pzanalysis.h:37
TPZMatrixSolver< STATE > * BuildPreconditioner(EPrecond preconditioner, bool overlap)
Define the type of preconditioner used.
int fStep
Time step.
Definition: pzanalysis.h:60
std::function< void(const TPZVec< REAL > &loc, TPZVec< STATE > &result, TPZFMatrix< STATE > &deriv)> fExact
Pointer to Exact solution function, it is necessary to calculating errors.
Definition: pzanalysis.h:99
TPZManVector< TPZManVector< REAL, 10 >, 100 > fvalues
Definition: pzanalysis.h:343
TPZMatrixSolver< STATE > & Solver()
Get the solver matrix.
Definition: pzanalysis.h:367
void SetExact(std::function< void(const TPZVec< REAL > &loc, TPZVec< STATE > &result, TPZFMatrix< STATE > &deriv)> f)
Sets the pointer of the exact solution function.
Definition: pzanalysis.h:360
virtual ~TPZAnalysis(void)
Destructor: deletes all protected dynamic allocated objects.
Definition: pzanalysis.cpp:166
int ClassId() const override
Define the class id associated with the class.
TPZFMatrix< STATE > & Solution()
Returns the solution matrix.
Definition: pzanalysis.h:177
Refines geometrical mesh (all the elements) num times.
Definition: pzstrmatrix.h:35
void CleanUp()
deletes all data structures
Definition: pzanalysis.cpp:171
TPZAutoPointer< TPZGuiInterface > fGuiInterface
Pointer for gui interface object.
Definition: pzanalysis.h:74
TPZAdmChunkVector< TPZCompEl * > fElvec
Definition: pzanalysis.h:326
TPZAutoPointer< TPZStructMatrix > StructMatrix()
Returns a reference to the structural matrix.
Definition: pzanalysis.h:182
TPZAnalysis()
Create an empty TPZAnalysis object.
Definition: pzanalysis.cpp:85
void SetThreadsForError(int nthreads)
Definition: pzanalysis.h:200
void ShowShape(const std::string &plotfile, TPZVec< int64_t > &equationindices)
Graphic of the solution as V3DGrap visualization.
Definition: pzanalysis.cpp:834
void SetTime(REAL time)
Sets time will be used in dx files.
Definition: pzanalysis.h:371
friend struct ThreadData
Definition: pzanalysis.h:353
TPZGraphMesh * GraphMesh(int dimension)
Defines the postprocessing parameters for the graphical grid.
Definition: pzanalysis.h:243
static const double tol
Definition: pzgeoprism.cpp:23
Implements the sequence of actions to perform a finite element analysis. Analysis.
Definition: pzanalysis.h:32
f
Definition: test.py:287
expr_ dx(i) *cos(expr_.val())
int nthreads
Definition: numatst.cpp:315
virtual void LoadSolution()
Load the solution into the computable grid.
Definition: pzanalysis.cpp:441
TPZAutoPointer< TPZRenumbering > fRenumber
Renumbering scheme.
Definition: pzanalysis.h:71
virtual void PostProcessErrorParallel(TPZVec< REAL > &, bool store_error=true, std::ostream &out=std::cout)
Definition: pzanalysis.cpp:652
TPZGeoMesh * fGeoMesh
Geometric Mesh.
Definition: pzanalysis.h:42
virtual void CloseGraphMesh()
Clean the GrapMesh vector.
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
TPZVec< std::string > fScalarNames[3]
Scalar variables names - to post process.
Definition: pzanalysis.h:54
void Read(TPZStream &buf, void *context) override
read objects from the stream
TPZVec< int64_t > fGeoElId
Definition: pzanalysis.h:79
TPZMatrixSolver< STATE > * BuildSequenceSolver(TPZVec< int64_t > &graph, TPZVec< int64_t > &graphindex, int64_t neq, int numcolors, TPZVec< int > &colors)
Build a sequence solver based on the block graph and its colors.
Datastructure which defines postprocessing for one dimensional meshes.
Definition: pzanalysis.h:77
TTablePostProcess fTable
Definition: pzanalysis.h:94
~TTablePostProcess()
TTablePostProcess()
Contains declaration of the TPZAutoPointer class which has Increment and Decrement actions are mutexe...
TPZVec< std::string > fVectorNames[3]
Vectorial variables names - to post process.
Definition: pzanalysis.h:56
void IdentifyPostProcessingMatIds(int dimension, std::set< int > &matids)
Fill mat ids with materials with provided dimension wich are not boundary conditinos or interface...
pthread_mutex_t fAccessElement
Mutexes (to choose which element is next)
Definition: pzanalysis.h:346
void Resequence(int firstel=-1)
Recompute the node sequence.
virtual void PostProcess(int resolution)
Draw solution over mesh for all dimensions.
int fDimension
Definition: pzanalysis.h:81
virtual void Run(std::ostream &out=std::cout)
Calls the appropriate sequence of methods to build a solution or a time stepping sequence.
Definition: pzanalysis.cpp:920
Contains TPZMatrix<TVar>class, root matrix class.
bool AmIKilled()
Returns if the process was canceled through gui interface.
Definition: pzanalysis.h:118
virtual void DefineElementTable(int dimension, TPZVec< int64_t > &GeoElIds, TPZVec< REAL > &points)
Fill the computational element vector to post processing depending over geometric mesh defined...
void SetStep(int step)
ste the step for post processing
Definition: pzanalysis.h:195
TPZGraphMesh * fGraphMesh[3]
Graphical mesh.
Definition: pzanalysis.h:46
TPZAutoPointer< TPZStructMatrix > fStructMatrix
Structural matrix.
Definition: pzanalysis.h:68
virtual void AssembleResidual()
Assemble the load vector.
Definition: pzanalysis.cpp:288
int GetStep()
Definition: pzanalysis.h:205
This class implements a geometric mesh for the pz environment. Geometry.
Definition: pzgmesh.h:48
Implements computational mesh. Computational Mesh.
Definition: pzcmesh.h:47
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
Contains the declaration of TPZFlopCounter class and TPZCounter struct.
virtual void LoadSolution(const TPZFMatrix< STATE > &sol)
Load the solution into the computable mesh considering sol as Solution vector of the analysis...
Definition: pzanalysis.h:278
TPZFMatrix< STATE > & Rhs()
Returns the load vector.
Definition: pzanalysis.h:174
virtual void PostProcessErrorSerial(TPZVec< REAL > &, bool store_error=true, std::ostream &out=std::cout)
Definition: pzanalysis.cpp:747
void CreateListOfCompElsToComputeError(TPZAdmChunkVector< TPZCompEl *> &elvec)
Definition: pzanalysis.cpp:586
TPZFMatrix< STATE > fRhs
Load vector.
Definition: pzanalysis.h:48
TPZVec< std::string > fTensorNames[3]
Tensorial variables names - to post process.
Definition: pzanalysis.h:58
int HighestDimension()
Returns the dimension of the material which has the highest dimension.
This class defines the interface to save and restore objects from TPZStream objects. Persistency.
Definition: TPZSavable.h:67
void AnimateRun(int64_t num_iter, int steps, TPZVec< std::string > &scalnames, TPZVec< std::string > &vecnames, const std::string &plotfile)
Run and print the solution step by step.
Defines the interface of a computational element. Computational Element.
Definition: pzcompel.h:59
virtual void PostProcessError(TPZVec< REAL > &, bool store_error=true, std::ostream &out=std::cout)
Compute the local error over all elements and global errors in several norms and print out without ca...
Definition: pzanalysis.cpp:573
TPZVec< REAL > fLocations
Definition: pzanalysis.h:82
void SetRenumber(TPZAutoPointer< TPZRenumbering > renumber)
Change the renumbering scheme.
Definition: pzanalysis.h:142
virtual void SetCompMesh(TPZCompMesh *mesh, bool mustOptimizeBandwidth)
Set the computational mesh of the analysis.
Definition: pzanalysis.cpp:115
TPZAutoPointer< TPZGuiInterface > GetGuiInterface() const
Gets gui interface object.
Definition: pzanalysis.h:113
REAL fTime
Time variable which is used in dx output.
Definition: pzanalysis.h:62
virtual void SetTableVariableNames(TPZVec< std::string > varnames)
Sets the names of the variables into the data structure for post processing.
int ComputeNumberofLoadCases()
Determine the number of load cases from the material objects and return its value.
Definition: pzanalysis.cpp:263
This class implements a reference counter mechanism to administer a dynamically allocated object...