NeoPZ
pzmultiphysicscompel.h
Go to the documentation of this file.
1 
6 #ifndef PZMULTIPHYSICCOMPELH
7 #define PZMULTIPHYSICCOMPELH
8 
9 #include <iostream>
10 
11 #include "pzmultiphysicselement.h"
12 #include "pzmaterialdata.h"
13 
14 #include "pzelctemp.h"
15 #include "pzreducedspace.h"
16 
17 template<class T>
18 class TPZTransform;
19 
21 template <class TGeometry>
23 
24 protected:
25 
28 
31 
33  typename TGeometry::IntruleType fIntRule;
34 
35 
36 
37 public:
44  TPZMultiphysicsCompEl(TPZCompMesh &mesh, TPZGeoEl *gel, int64_t &index);
47 
53  std::map<int64_t,int64_t> & gl2lcConMap,
54  std::map<int64_t,int64_t> & gl2lcElMap);
55 
57  virtual ~TPZMultiphysicsCompEl();
58 
60  virtual TPZManVector<TPZCompElSide,5> &ElementVec() override { return fElementVec; }
61 
66  virtual void AffineTransform(TPZVec<TPZTransform<> > &trVec) const override;
67 
74  virtual void EvaluateError(std::function<void(const TPZVec<REAL> &loc,TPZVec<STATE> &val,TPZFMatrix<STATE> &deriv)> func,
75  TPZVec<REAL> &errors, bool store_error ) override;
76 
77  virtual void EvaluateError(TPZFunction<STATE> &func,
78  TPZVec<STATE> &errors, bool store_error) override;
79 
80 
86  void GetReferenceIndexVec(TPZManVector<TPZCompMesh *> cmeshVec, std::set<int64_t> &refIndexVec);
87 
89  virtual TPZCompEl *Clone(TPZCompMesh &mesh) const override;
90 
102  virtual TPZCompEl *ClonePatchEl(TPZCompMesh &mesh,
103  std::map<int64_t,int64_t> & gl2lcConMap,
104  std::map<int64_t,int64_t> & gl2lcElMap) const override;
105 
107  virtual int NConnects() const override;
108 
110  virtual int64_t NMeshes() override
111  {
112  return fElementVec.size();
113  }
114 
119  virtual int64_t ConnectIndex(int i) const override;
120 
121  virtual int64_t ConnectIndex(int elem, int connect) const ;
122 
124  virtual int Dimension() const override;
125 
127  virtual void PolynomialOrder(TPZVec<int> &order) const override;
128 
139  virtual void Integrate(int variable, TPZVec<STATE> & value) override;
140 
142  virtual void ComputeRequiredData(TPZVec<REAL> &point, TPZVec<TPZTransform<> > &trvec, TPZVec<TPZMaterialData> &datavec);
143 
144  virtual void ComputeRequiredData(TPZMaterialData &data, TPZVec<REAL> &point) override
145  {
147  }
148 
159  virtual void Solution(TPZVec<REAL> &qsi,int var,TPZVec<STATE> &sol) override;
160 
164  virtual TPZVec<STATE> IntegrateSolution(int var) const override;
165 
166 
174  virtual void ComputeSolution(TPZVec<REAL> &qsi,
175  TPZSolVec &sol, TPZGradSolVec &dsol,TPZFMatrix<REAL> &axes) override;
176 
189  virtual void ComputeSolution(TPZVec<REAL> &qsi,
190  TPZVec<REAL> &normal,
191  TPZSolVec &leftsol, TPZGradSolVec &dleftsol,TPZFMatrix<REAL> &leftaxes,
192  TPZSolVec &rightsol, TPZGradSolVec &drightsol,TPZFMatrix<REAL> &rightaxes) override;
193 
203  virtual void ComputeSolution(TPZVec<REAL> &qsi, TPZFMatrix<REAL> &phi, TPZFMatrix<REAL> &dphix,
204  const TPZFMatrix<REAL> &axes, TPZSolVec &sol, TPZGradSolVec &dsol) override;
205 
211  virtual void SetConnectIndex(int inode, int64_t index) override;
212 
213 
215  virtual void SetCreateFunctions(TPZCompMesh *mesh) override;
216 
218  virtual void AddElement(TPZCompEl *cel, int64_t meshindex) override
219  {
220  if (fElementVec.size() <= meshindex)
221  {
222  fElementVec.resize(meshindex+1);
223  fActiveApproxSpace.Resize(meshindex+1, 1);
224  }
225  if (cel)
226  {
227  TPZGeoEl *gel = cel->Reference();
228  TPZCompElSide celside(cel,gel->NSides()-1);
229  fElementVec[meshindex] = celside;
230  }
231  else
232  {
233  fElementVec[meshindex] = TPZCompElSide();
234  }
235  }
236 
238  virtual void AddElement(const TPZCompElSide &celside, int64_t meshindex) override
239  {
240  if (fElementVec.size() <= meshindex)
241  {
242  fElementVec.resize(meshindex+1);
243  fActiveApproxSpace.Resize(meshindex+1, 1);
244  }
245  fElementVec[meshindex] = celside;
246  }
247 
248  virtual TPZCompEl *Element(int64_t elindex) override
249  {
250  return fElementVec[elindex].Element();
251  }
252 
254  virtual TPZCompEl *ReferredElement(int64_t mesh) override
255  {
256 
257 #ifdef PZDEBUG
258  if (fElementVec.size() <= mesh) {
259  PZError << "Error at " << __PRETTY_FUNCTION__ << " index does not exist!\n";
260  DebugStop();
261  };
262 #endif
263 
264  return fElementVec[mesh].Element();
265  }
266 
267 
272  virtual void SetConnectIndexes(TPZVec<int64_t> &indexes) override
273  {
274  fConnectIndexes = indexes;
275  }
276 
281  virtual void Print(std::ostream &out = std::cout) const override;
282 
288  virtual void CalcStiff(TPZElementMatrix &ek, TPZElementMatrix &ef) override;
289 
295  virtual void CalcResidual(TPZElementMatrix &ef) override;
296 
299 
302 
307  void InitMaterialData(TPZVec<TPZMaterialData > &dataVec, TPZVec<int64_t> *indices = 0) override;
308 
309  virtual void CreateGraphicalElement(TPZGraphMesh &grmesh, int dimension) override;
310 
311  //virtual void CreateGraphicalElement(TPZGraphMesh &grmesh, std::set<int> dimension, std::set<int> MaterialID);
312 
313  virtual void SetIntegrationRule(int int_order) override;
314 
316  virtual void InitializeIntegrationRule() override;
317 
319  virtual const TPZIntPoints &GetIntegrationRule() const override;
320 
322  virtual TPZIntPoints &GetIntegrationRule() override;
323 
325  virtual int NumberOfCompElementsInsideThisCompEl() override {
326  return fElementVec.NElements();
327  }
328  public:
329 virtual int ClassId() const override;
330 
331 };
332 
333 
335 TPZCompEl *CreateMultiphysicsPointEl(TPZGeoEl *gel,TPZCompMesh &mesh,int64_t &index);
336 
338 TPZCompEl *CreateMultiphysicsLinearEl(TPZGeoEl *gel,TPZCompMesh &mesh,int64_t &index);
339 
341 TPZCompEl *CreateMultiphysicsQuadEl(TPZGeoEl *gel,TPZCompMesh &mesh,int64_t &index);
342 
344 TPZCompEl *CreateMultiphysicsTriangleEl(TPZGeoEl *gel,TPZCompMesh &mesh,int64_t &index);
345 
347 TPZCompEl *CreateMultiphysicsCubeEl(TPZGeoEl *gel,TPZCompMesh &mesh,int64_t &index);
348 
350 TPZCompEl *CreateMultiphysicsPrismEl(TPZGeoEl *gel,TPZCompMesh &mesh,int64_t &index);
351 
353 TPZCompEl *CreateMultiphysicsPyramEl(TPZGeoEl *gel,TPZCompMesh &mesh,int64_t &index);
354 
356 TPZCompEl *CreateMultiphysicsTetraEl(TPZGeoEl *gel,TPZCompMesh &mesh,int64_t &index);
357 
358 //--------------------- WITH MEMORY ----------------------
359 
362 
365 
368 
371 
374 
377 
380 
383 
384 #include "pzcmesh.h"
385 
386 template<class TGeometry>
389 }
390 
391 #endif
Represents a graphical mesh used for post processing purposes. Post processing.
Definition: pzgraphmesh.h:34
TPZCompEl * CreateMultiphysicsQuadElWithMem(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational quadrilateral element for Multiphysics approximate space.
TPZManVector< TPZCompElSide,5 > fElementVec
List of pointers to computational elements.
virtual void ComputeRequiredData(TPZVec< REAL > &point, TPZVec< TPZTransform<> > &trvec, TPZVec< TPZMaterialData > &datavec, TPZVec< int64_t > indices)
TPZCompEl * CreateMultiphysicsPrismElWithMem(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational prismal element for Multiphysics approximate space.
TPZCompEl * CreateMultiphysicsPyramEl(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational pyramidal element for Multiphysics approximate space.
Contains declaration of TPZIntelGen class which implements a generic computational element...
virtual void SetConnectIndex(int inode, int64_t index) override
Set the index i to node inode.
Implements computational element and a side. Computational Element.
Definition: pzcompel.h:632
Implements a vector class which allows to use external storage provided by the user. Utility.
Definition: pzquad.h:16
virtual int64_t NMeshes() override
Return the number of meshes associated with the element.
void GetReferenceIndexVec(TPZManVector< TPZCompMesh *> cmeshVec, std::set< int64_t > &refIndexVec)
Method to obtain an reference index set of multiphysics computational elements.
TPZCompEl * CreateMultiphysicsLinearElWithMem(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational linear element for Multiphysics approximate space.
virtual void resize(const int64_t newsize)
Definition: pzvec.h:213
virtual void AffineTransform(TPZVec< TPZTransform<> > &trVec) const override
Compute the map of a paramenter point in the multiphysic element to a parameter point in the super el...
clarg::argInt dimension("-d", "Matrices dimension M x M", 1000)
virtual TPZCompEl * ReferredElement(int64_t mesh) override
Returns referred element of this.
Contains the declaration of the TPZMultiphysicsElement class. This class is abstract.
TPZCompEl * CreateMultiphysicsTriangleElWithMem(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational triangular element for Multiphysics approximate space.
REAL val(STATE &number)
Returns value of the variable.
Definition: pzartdiff.h:23
virtual void InitializeIntegrationRule() override
After adding the elements initialize the integration rule.
virtual int NSides() const =0
Returns the number of connectivities of the element.
TPZCompEl * CreateMultiphysicsTriangleEl(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational triangular element for Multiphysics approximate space.
virtual TPZCompEl * ClonePatchEl(TPZCompMesh &mesh, std::map< int64_t, int64_t > &gl2lcConMap, std::map< int64_t, int64_t > &gl2lcElMap) const override
Method for creating a copy of the element in a patch mesh.
TPZMultiphysicsCompEl()
Default constructor.
virtual void Solution(TPZVec< REAL > &qsi, int var, TPZVec< STATE > &sol) override
Post processing method which computes the solution for the var post processed variable.
virtual int NumberOfCompElementsInsideThisCompEl() override
Return the size of the elementvec in multiphysics, if it is not multiphysics, just return 1...
virtual void Resize(const int64_t newsize, const T &object)
Resizes the vector object.
Definition: pzmanvector.h:426
Abstract class defining integration rules. Numerical Integration.
Definition: tpzintpoints.h:19
TPZManVector< int, 5 > fActiveApproxSpace
List of active approximation spaces.
virtual void ComputeRequiredData(TPZMaterialData &data, TPZVec< REAL > &point) override
Compute and fill data with requested attributes.
virtual void Print(std::ostream &out=std::cout) const override
Prints element data.
int64_t size() const
Returns the number of elements of the vector.
Definition: pzvec.h:196
TPZCompEl * CreateMultiphysicsPyramElWithMem(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational pyramidal element for Multiphysics approximate space.
Defines the behaviour of all geometric elements. GeometryTPZGeoEl is the common denominator for all g...
Definition: pzgeoel.h:43
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
virtual void SetCreateFunctions(TPZCompMesh *mesh) override
Sets create function in TPZCompMesh to create elements of this type.
void InitializeElementMatrix(TPZElementMatrix &ek, TPZElementMatrix &ef)
Initialize element matrix in which is computed CalcStiff.
virtual int Dimension() const override
Dimension of the element.
TPZCompEl * CreateMultiphysicsTetraEl(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational tetrahedral element for Multiphysics approximate space.
virtual ~TPZMultiphysicsCompEl()
Default destructor.
virtual TPZVec< STATE > IntegrateSolution(int var) const override
Compute the integral of a variable.
virtual void PolynomialOrder(TPZVec< int > &order) const override
virtual void CalcStiff(TPZElementMatrix &ek, TPZElementMatrix &ef) override
Computes the element stiffness matrix and right hand side.
Contains the TPZMaterialData class which implements an interface between TPZCompEl::CalcStiff and TPZ...
Contains declaration of TPZCompMesh class which is a repository for computational elements...
virtual TPZCompEl * Clone(TPZCompMesh &mesh) const override
Method for creating a copy of the element.
void InitMaterialData(TPZVec< TPZMaterialData > &dataVec, TPZVec< int64_t > *indices=0) override
Initialize a material data vector and its attributes based on element dimension, number of state vari...
TPZCompEl * CreateMultiphysicsTetraElWithMem(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational tetrahedral element for Multiphysics approximate space.
This class associates an element matrix with the coeficients of its contribution in the global stiffn...
Definition: pzelmat.h:30
virtual int ClassId() const override
Define the class id associated with the class.
class to create a compute element multiphysics
virtual const TPZIntPoints & GetIntegrationRule() const override
Returns a reference to an integration rule suitable for integrating the interior of the element...
TPZCompEl * CreateMultiphysicsPointElWithMem(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational point element for Multiphysics approximate space.
virtual void CalcResidual(TPZElementMatrix &ef) override
Computes the element stiffness matrix and right hand side.
TPZGeoEl * Reference() const
Return a pointer to the corresponding geometric element if such exists, return 0 otherwise.
Definition: pzcompel.cpp:1137
TPZCompEl * CreateMultiphysicsPointEl(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational point element for Multiphysics approximate space.
TPZCompEl * CreateMultiphysicsCubeElWithMem(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational cube element for Multiphysics approximate space.
virtual void ComputeRequiredData(TPZVec< REAL > &point, TPZVec< TPZTransform<> > &trvec, TPZVec< TPZMaterialData > &datavec)
Compute and fill data with requested attributes for each of the compels in fElementVec.
TPZCompEl * CreateMultiphysicsPrismEl(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational prismal element for Multiphysics approximate space.
virtual int64_t ConnectIndex(int i) const override
Returns the index of the ith connectivity of the element.
Implements computational mesh. Computational Mesh.
Definition: pzcmesh.h:47
virtual void CreateGraphicalElement(TPZGraphMesh &grmesh, int dimension) override
Creates corresponding graphical element(s) if the dimension matches graphical elements are used to ge...
TPZVec< int64_t > fConnectIndexes
Indexes of the connects of the element.
virtual void EvaluateError(std::function< void(const TPZVec< REAL > &loc, TPZVec< STATE > &val, TPZFMatrix< STATE > &deriv)> func, TPZVec< REAL > &errors, bool store_error) override
Performs an error estimate on the elemen.
virtual TPZManVector< TPZCompElSide, 5 > & ElementVec() override
Returns a reference to the element pointers vector.
virtual void ComputeSolution(TPZVec< REAL > &qsi, TPZSolVec &sol, TPZGradSolVec &dsol, TPZFMatrix< REAL > &axes) override
Computes solution and its derivatives in the local coordinate qsi.
virtual TPZCompEl * Element(int64_t elindex) override
int64_t NElements() const
Returns the number of elements of the vector.
Definition: pzvec.h:190
virtual void SetConnectIndexes(TPZVec< int64_t > &indexes) override
Sets indexes of the connects of the element.
Implements an affine transformation between points in parameter space. Topology Utility.
Definition: pzmganalysis.h:14
Contains the declaration of the Reduced Space class.
TPZCompEl * CreateMultiphysicsLinearEl(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational linear element for Multiphysics approximate space.
Defines the interface of a computational element. Computational Element.
Definition: pzcompel.h:59
void SetAllCreateFunctionsMultiphysicElem()
Definition: pzcmesh.h:542
virtual int NConnects() const override
Returns the number of nodes of the element.
virtual void SetIntegrationRule(int int_order) override
TPZCompEl * CreateMultiphysicsQuadEl(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational quadrilateral element for Multiphysics approximate space.
TPZCompEl * CreateMultiphysicsCubeEl(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index)
Creates computational cube element for Multiphysics approximate space.
#define PZError
Defines the output device to error messages and the DebugStop() function.
Definition: pzerror.h:15
TGeometry::IntruleType fIntRule
Integration rule associated with the element.
virtual void AddElement(TPZCompEl *cel, int64_t meshindex) override
add an element to the datastructure
virtual void Integrate(int variable, TPZVec< STATE > &value) override
Post processing method which computes the solution for the var post processed variable.
virtual void AddElement(const TPZCompElSide &celside, int64_t meshindex) override
add an element to the datastructure