NeoPZ
TPZCompElDisc.h
Go to the documentation of this file.
1 
8 #ifndef ELCOMPDISCHPP
9 #define ELCOMPDISCHPP
10 
11 #include <iostream>
12 #include "pzfmatrix.h"
13 #include "pzvec.h"
14 #include "pzcompel.h"
15 #include "pzinterpolationspace.h"
16 #include "pzgeoel.h"
17 #include "pzreal.h"
18 #include "TPZShapeDisc.h"
19 #include "tpzautopointer.h"
20 #include "pzquad.h"
21 #include "pzfunction.h"
22 #include "pzfmatrix.h"
23 
24 struct TPZElementMatrix;
25 class TPZBndCond;
26 class TPZConnect;
27 class TPZMaterial;
28 class TPZGeoEl;
29 class TPZCompMesh;
30 
36 
37 private:
38 
40 
41 protected:
42 
44 
47 
48 public:
49 
51  static void SetTensorialShape(TPZCompMesh * cmesh);
52 
54  static void SetTotalOrderShape(TPZCompMesh * cmesh);
55 
57  void SetTensorialShape();
58 
60  void SetTotalOrderShape();
61 
64  void SetTensorialShapeFull();
65 
69 
70 protected:
71 
73  int64_t fConnectIndex;
74 
76  REAL fConstC;
77 
79  bool fUseQsiEta;
80 
83 
84 protected:
85 
88 
93  virtual int64_t CreateMidSideConnect();
94 
95 public:
96 
99 
102 
103  int GetMaterial( const TPZGeoElSide& gside );
104 
106  static TPZCompEl *CreateDisc(TPZGeoEl *geo, TPZCompMesh &mesh, int64_t &index);
107 
112  static void SetOrthogonalFunction(void (*orthogonal)(REAL C, REAL x0, REAL x,int degree, TPZFMatrix<REAL> & phi, TPZFMatrix<REAL> & dphi, int n)){
114  }
115 
117  virtual void SetInnerRadius(REAL InnerRadius) {PZError << "TPZCompElDisc::SetInnerRadius - This method should never be called because the inner" << std::endl
118  << "radius is not stored in TPZCompElDisc. It is stored in TPZAgglomerateElement." << std::endl;}
119 
121  virtual void SetNInterfaces(int nfaces) {PZError << "TPZCompElDisc::SetNFaces - This method should never be called because the number of interfaces" << std::endl
122  << "is not stored in TPZCompElDisc. It is only stored by TPZAgglomerateElement." << std::endl;}
123 
125  virtual int NInterfaces();
126 
128  TPZCompElDisc();
130  TPZCompElDisc(TPZCompMesh &mesh,TPZGeoEl *ref,int64_t &index);//original
132  TPZCompElDisc(TPZCompMesh &mesh,int64_t &index);//construtor do aglomerado
134  TPZCompElDisc(TPZCompMesh &mesh, const TPZCompElDisc &copy);
135 
144  const TPZCompElDisc &copy,
145  std::map<int64_t,int64_t> &gl2lcConMap,
146  std::map<int64_t,int64_t> &gl2lcElMap);
147 
148  TPZCompElDisc(TPZCompMesh &mesh, const TPZCompElDisc &copy,int64_t &index);
149 
151  virtual void SetCreateFunctions(TPZCompMesh *mesh) override;
152 
153  virtual TPZCompEl *Clone(TPZCompMesh &mesh) const override {
154  return new TPZCompElDisc(mesh,*this);
155  }
156 
157  virtual TPZCompEl *Clone(TPZCompMesh &mesh,int64_t &index) const {
158  return new TPZCompElDisc(mesh,*this,index);
159  }
160 
165  std::map<int64_t,int64_t> & gl2lcConMap,
166  std::map<int64_t,int64_t> & gl2lcElMap) const override {
167  return new TPZCompElDisc(mesh,*this,gl2lcConMap,gl2lcElMap);
168  }
170  virtual ~TPZCompElDisc();
171 
173  void Divide(int64_t index, TPZVec<int64_t> &subindex, int interpolate = 0) override;
174 
182  virtual void Shape(TPZVec<REAL> &qsi,TPZFMatrix<REAL> &phi,TPZFMatrix<REAL> &dphidxi) override;
183 
184 
185 protected:
186  virtual void ComputeShape(TPZVec<REAL> &intpoint, TPZVec<REAL> &X, TPZFMatrix<REAL> &jacobian, TPZFMatrix<REAL> &axes, REAL &detjac, TPZFMatrix<REAL> &jacinv, TPZFMatrix<REAL> &phi, TPZFMatrix<REAL> &dphi, TPZFMatrix<REAL> &dphidx) override;
187 
188  public:
189 
195  virtual void ComputeShape(TPZVec<REAL> &intpoint,TPZMaterialData &data) override;
196 
199  virtual void ShapeX(TPZVec<REAL> &X, TPZFMatrix<REAL> &phi, TPZFMatrix<REAL> &dphi);
200 
203 
205  virtual const TPZIntPoints &GetIntegrationRule() const override;
206 
208  virtual TPZIntPoints &GetIntegrationRule() override;
209 
211  virtual MElementType Type() override {return EDiscontinuous;}
212 
214  REAL ConstC() const {return fConstC;}
215 
216  void SetConstC(REAL c){fConstC = c;}
217 
219 #ifdef PZDEBUG
220  if(!Reference()) DebugStop();
221 #endif
222  fUseQsiEta = true;
223  fCenterPoint.Fill(0.);
224  fConstC = 1.;
225 
226  TPZGeoEl *gel = Reference();
227  if (gel)
228  {
229  int ns = Reference()->NSides();
230  fCenterPoint.resize(gel->Dimension());
231  this->Reference()->CenterPoint(ns-1, fCenterPoint);
232  }
233  }
234 
236 #ifdef PZDEBUG
237  if(!Reference()) DebugStop();
238 #endif
239  fUseQsiEta = false;
241  Reference()->CenterPoint(Reference()->NSides()-1,csi);
242  Reference()->X(csi,fCenterPoint);
243  fConstC = NormalizeConst();
244  }
245 
247  void InternalPoint(TPZVec<REAL> &point);
248 
250  virtual void Print(std::ostream & out = std::cout) const override;
251 
253  virtual int Degree() const;
254 
256  virtual void SetDegree(int degree);
258  virtual int NConnects() const override;
259 
261  int NCornerConnects() const { return Reference()->NNodes();}
262 
264  virtual void BuildCornerConnectList(std::set<int64_t> &connectindexes) const override;
265 
267  int Dimension() const override { return Reference()->Dimension();}
268 
270  virtual REAL NormalizeConst();
271 
273  int64_t ConnectIndex(int side = 0) const override;
274  void SetConnectIndex(int /*inode*/, int64_t index) override {fConnectIndex = index;}
275 
277  virtual int NSideConnects(int iside) const override
278  {
279  return NConnects();
280  }
281 
287  virtual int SideConnectLocId(int icon,int is) const override
288  {
289 #ifdef PZDEBUG
290  if (icon != 0) {
291  DebugStop();
292  }
293 #endif
294  return 0;
295  }
296 
297 
299  virtual int NShapeF() const override;
300 
302  virtual int MaxOrder() override;
303 
306 
308  virtual int NConnectShapeF(int inod, int order) const override;
309 
310  REAL CenterPoint(int index) const {return fCenterPoint[index];}
311 
312  virtual void CenterPoint(TPZVec<REAL> &center);
313 
314  void SetCenterPoint(int i,REAL x){fCenterPoint[i] = x;}
315 
316  REAL SizeOfElement();
317 
319  virtual REAL VolumeOfEl() override { return Reference()->Volume(); }
320 
327  void CreateGraphicalElement(TPZGraphMesh &grmesh, int dimension) override;
328 
336  void SolutionX(TPZVec<REAL> &x,TPZVec<STATE> &uh);
337 
338  public:
346  virtual void ComputeSolution(TPZVec<REAL> &qsi, TPZSolVec &sol, TPZGradSolVec &dsol,TPZFMatrix<REAL> & axes) override;
347 
348 public:
349 
355  virtual void ComputeSolution(TPZVec<REAL> &qsi, TPZMaterialData &data) override;
356 
366  virtual void ComputeSolution(TPZVec<REAL> &qsi, TPZFMatrix<REAL> &phi, TPZFMatrix<REAL> &dphix,
367  const TPZFMatrix<REAL> &axes, TPZSolVec &sol, TPZGradSolVec &dsol) override;
368 
381  virtual void ComputeSolution(TPZVec<REAL> &qsi,
382  TPZVec<REAL> &normal,
383  TPZSolVec &leftsol, TPZGradSolVec &dleftsol,TPZFMatrix<REAL> &leftaxes,
384  TPZSolVec &rightsol, TPZGradSolVec &drightsol,TPZFMatrix<REAL> &rightaxes) override;
385 
386  virtual void AccumulateIntegrationRule(int degree, TPZStack<REAL> &point, TPZStack<REAL> &weight);
387 
389  virtual void AccumulateVertices(TPZStack<TPZGeoNode *> &nodes);
390 
391  int NSides();
392 
393  void BuildTransferMatrix(TPZCompElDisc &coarsel, TPZTransfer<STATE> &transfer);
394 
396  public:
397 int ClassId() const override;
398 
400  void Write(TPZStream &buf, int withclassid) const override;
401 
403  void Read(TPZStream &buf, void *context) override;
404 
406  virtual void SetPreferredOrder ( int order ) override { SetDegree( order ); }
407 
413  virtual void PRefine ( int order ) override { SetDegree( order ); }
414 
424 
426  static void EvaluateSquareResidual2D(TPZCompMesh &cmesh, TPZVec<REAL> &error, bool verbose = false);
427 
428 
429 };
430 
431 inline TPZCompEl *TPZCompElDisc::CreateDisc(TPZGeoEl *geo, TPZCompMesh &mesh, int64_t &index) {
432  if(!geo->Reference() && geo->NumInterfaces() == 0)
433  return new TPZCompElDisc(mesh,geo,index);
434  return NULL;
435 }
436 //Exemplo do quadrilatero:
437 //acessar com -> TPZGeoElement<TPZShapeQuad,TPZGeoQuad,TPZRefQuad>::SetCreateFunction(TPZCompElDisc::CreateDisc);
438 
439 #endif
Represents a graphical mesh used for post processing purposes. Post processing.
Definition: pzgraphmesh.h:34
int GetMaterial(const TPZGeoElSide &gside)
REAL Volume()
Return the volume of the element.
Definition: pzgeoel.cpp:1052
virtual void AccumulateIntegrationRule(int degree, TPZStack< REAL > &point, TPZStack< REAL > &weight)
Contains the TPZInt1d, TPZIntTriang, TPZIntQuad, TPZIntCube3D, TPZIntTetra3D, TPZIntPyram3D and TPZIn...
REAL CenterPoint(int index) const
void AppendExternalShapeFunctions(TPZVec< REAL > &X, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphi)
Add extenal shape function into already computed phi and dphi discontinuous functions.
virtual int NShapeF() const override
Returns the shapes number of the element.
Represents a set of shape functions associated with a computational element/side. Computational Eleme...
Definition: pzconnect.h:30
void SolutionX(TPZVec< REAL > &x, TPZVec< STATE > &uh)
Computes the solution in function of a point in cartesian space.
virtual void SetCreateFunctions(TPZCompMesh *mesh) override
Set create function in TPZCompMesh to create elements of this type.
virtual TPZCompEl * Clone(TPZCompMesh &mesh, int64_t &index) const
virtual int NConnects() const override
Returns the number of connects.
virtual int SideConnectLocId(int icon, int is) const override
Returns the local node number of icon along is.
virtual void resize(const int64_t newsize)
Definition: pzvec.h:213
clarg::argInt dimension("-d", "Matrices dimension M x M", 1000)
static void SetOrthogonalFunction(void(*orthogonal)(REAL C, REAL x0, REAL x, int degree, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphi, int n))
Sets the orthogonal function which will be used throughout the program.
Contains declaration of TPZCompEl class which defines the interface of a computational element...
Templated vector implementation.
bool fUseQsiEta
Variable to choose the qsi point or the X point in the calculus of the phis and dphis.
Definition: TPZCompElDisc.h:79
void SetFalseUseQsiEta()
void Read(TPZStream &buf, void *context) override
Read the element data from a stream.
virtual int64_t CreateMidSideConnect()
It creates new conect that it associates the degrees of freedom of the element and returns its index...
virtual int NSideConnects(int iside) const override
Returns the number of dof nodes along side iside.
virtual const TPZIntPoints & GetIntegrationRule() const override
Returns a reference to an integration rule suitable for integrating the interior of the element...
virtual void ComputeShape(TPZVec< REAL > &intpoint, TPZVec< REAL > &X, TPZFMatrix< REAL > &jacobian, TPZFMatrix< REAL > &axes, REAL &detjac, TPZFMatrix< REAL > &jacinv, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphi, TPZFMatrix< REAL > &dphidx) override
Compute shape functions based on master element in the classical FEM manner.
Contains declaration of TPZInterpolationSpace class which implements the interface for interpolated c...
void degree(int root, int adj_num, int adj_row[], int adj[], int mask[], int deg[], int *iccsze, int ls[], int node_num)
Definition: rcm.cpp:875
void InternalPoint(TPZVec< REAL > &point)
Returns the center point.
Utility class which represents an element with its side. The Geometric approximation classes Geometry...
Definition: pzgeoelside.h:83
virtual TPZCompEl * ClonePatchEl(TPZCompMesh &mesh, std::map< int64_t, int64_t > &gl2lcConMap, std::map< int64_t, int64_t > &gl2lcElMap) const override
void BuildTransferMatrix(TPZCompElDisc &coarsel, TPZTransfer< STATE > &transfer)
void SetExternalShapeFunction(TPZAutoPointer< TPZFunction< STATE > > externalShapes)
Define external shape functions which are stored in class attribute fExternalShape.
virtual int NSides() const =0
Returns the number of connectivities of the element.
virtual void ComputeSolution(TPZVec< REAL > &qsi, TPZSolVec &sol, TPZGradSolVec &dsol, TPZFMatrix< REAL > &axes) override
Computes solution and its derivatives in the local coordinate qsi.
void CreateGraphicalElement(TPZGraphMesh &grmesh, int dimension) override
Creates corresponding graphical element(s) if the dimension matches graphical elements are used to ge...
Abstract class defining integration rules. Numerical Integration.
Definition: tpzintpoints.h:19
virtual void PRefine(int order) override
Change the preferred order for the element and proceed the adjust of the aproximation space taking in...
void SetConnectIndex(int, int64_t index) override
Set the index i to node inode.
void error(char *string)
Definition: testShape.cc:7
This abstract class defines the behaviour which each derived class needs to implement.
Definition: TPZMaterial.h:39
virtual REAL NormalizeConst()
Calculates the normalizing constant of the bases of the element.
virtual void CenterPoint(int side, TPZVec< REAL > &masscent) const =0
It returns the coordinates from the center of the side of the element in the element coordinate space...
void Write(TPZStream &buf, int withclassid) const override
Save the element data to a stream.
TPZCompElDisc()
Default constructor.
bool HasExternalShapeFunction()
Return whether element has external shape functions set to.
void SetTensorialShape()
Set tensorial shape functions.
void SetTotalOrderShape()
Set total order shape functions.
TPZAutoPointer< TPZIntPoints > CreateIntegrationRule() const
virtual void SetPreferredOrder(int order) override
Define the desired order for entire element.
Defines the behaviour of all geometric elements. GeometryTPZGeoEl is the common denominator for all g...
Definition: pzgeoel.h:43
Contains TPZMatrixclass which implements full matrix (using column major representation).
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
This class defines the boundary condition for TPZMaterial objects.
Definition: pzbndcond.h:29
int ClassId() const override
Returns the unique identifier for reading/writing objects to streams.
virtual void SetNInterfaces(int nfaces)
Sets element&#39;s number of interfaces.
int NumInterfaces()
Returns number of TPZInterfaceElement pointing to this.
Definition: pzgeoel.h:94
TPZAutoPointer< TPZFunction< STATE > > fExternalShape
A pz function to allow the inclusion of extra shape functions which are defined externally.
Definition: TPZCompElDisc.h:82
int64_t fConnectIndex
It preserves index of connect associated to the element.
Definition: TPZCompElDisc.h:73
TPZCompEl * Reference() const
Return a pointer to the element referenced by the geometric element.
Definition: pzgeoel.h:750
Contains declaration of the TPZAutoPointer class which has Increment and Decrement actions are mutexe...
virtual void AccumulateVertices(TPZStack< TPZGeoNode *> &nodes)
Accumulates the vertices of the agglomerated elements.
virtual int NConnectShapeF(int inod, int order) const override
Returns the number of shapefunctions associated with a connect.
int MaxOrderExceptExternalShapes()
Returns the max order of interpolation excluding external shape order.
pzshape::TPZShapeDisc::MShapeType fShapefunctionType
Shape function type used by the element.
Definition: TPZCompElDisc.h:46
int64_t ConnectIndex(int side=0) const override
Returns the connect index from the element.
TPZManVector< REAL, 3 > fCenterPoint
It keeps the interior point coordinations of the element.
Definition: TPZCompElDisc.h:87
virtual TPZCompEl * Clone(TPZCompMesh &mesh) const override
Method for creating a copy of the element.
virtual ~TPZCompElDisc()
Default destructor.
virtual int NNodes() const =0
Returns the number of nodes of the element.
This class associates an element matrix with the coeficients of its contribution in the global stiffn...
Definition: pzelmat.h:30
Implements the interfaces for TPZCompElDisc, TPZInterfaceElement and TPZInterpolatedElement. Computational element.
static TPZCompEl * CreateDisc(TPZGeoEl *geo, TPZCompMesh &mesh, int64_t &index)
Creates discontinuous computational element.
void SetTensorialShapeFull()
Set tensorial shape functions with many derivatives.
virtual int MaxOrder() override
Returns the max order of interpolation.
virtual int Dimension() const =0
Returns the dimension of the element.
virtual void Print(std::ostream &out=std::cout) const override
Prints the features of the element.
virtual void Shape(TPZVec< REAL > &qsi, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphidxi) override
Computes the shape function set at the point x. This method uses the order of interpolation of the el...
TPZGeoEl * Reference() const
Return a pointer to the corresponding geometric element if such exists, return 0 otherwise.
Definition: pzcompel.cpp:1137
virtual void X(TPZVec< REAL > &qsi, TPZVec< REAL > &result) const =0
Return the coordinate in real space of the point coordinate in the master element space...
MElementType
Define the element types.
Definition: pzeltype.h:52
int verbose
Definition: decompose.cpp:67
Implements computational mesh. Computational Mesh.
Definition: pzcmesh.h:47
int NCornerConnects() const
Amount of vertices of the element.
void SetConstC(REAL c)
virtual void SetInnerRadius(REAL InnerRadius)
Sets the inner radius value.
static void(* fOrthogonal)(REAL C, REAL x0, REAL x, int degree, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphi, int n)
Pointer to orthogonal function to use as shape function.
Definition: TPZShapeDisc.h:56
void Fill(const T &copy, const int64_t from=0, const int64_t numelem=-1)
Will fill the elements of the vector with a copy object.
Definition: pzvec.h:460
Contains the declaration of the shape function discontinuous.
virtual void SetDegree(int degree)
Assigns the degree of the element.
virtual int Degree() const
Returns the degree of interpolation of the element.
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
Contains the declaration of TPZFlopCounter class and TPZCounter struct.
void SetCenterPoint(int i, REAL x)
REAL fConstC
Normalizing constant for shape functions.
Definition: TPZCompElDisc.h:76
void Divide(int64_t index, TPZVec< int64_t > &subindex, int interpolate=0) override
Divide the computational element.
void SetTrueUseQsiEta()
This class implements a discontinuous element (for use with discontinuous Galerkin). Computational Element.
Definition: TPZCompElDisc.h:35
virtual REAL VolumeOfEl() override
Returns the volume of the geometric element associated.
virtual void ShapeX(TPZVec< REAL > &X, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphi)
virtual REAL InnerRadius()
Returns the inner radius value.
Implements rectangular matrix which extends a solution vector of the coarse mesh to a solution vector...
Definition: pzcmesh.h:33
virtual int NInterfaces()
Returns the number of interfaces.
int Dimension() const override
Returns dimension from the element.
Defines the interface of a computational element. Computational Element.
Definition: pzcompel.h:59
void SetTotalOrderShapeFull()
Set total order shape functions.
virtual void BuildCornerConnectList(std::set< int64_t > &connectindexes) const override
adds the connect indexes associated with base shape functions to the set
virtual MElementType Type() override
Type of the element.
#define PZError
Defines the output device to error messages and the DebugStop() function.
Definition: pzerror.h:15
static REAL EvaluateSquareResidual2D(TPZInterpolationSpace *cel)
Compute the integral of the square residual over the element domain.
REAL ConstC() const
It returns the constant that normalizes the bases of the element.
TPZAutoPointer< TPZIntPoints > fIntRule
Definition: TPZCompElDisc.h:43