NeoPZ
pzinterpolationspace.h
Go to the documentation of this file.
1 
6 #ifndef PZINTERPOLATIONSPACE_H
7 #define PZINTERPOLATIONSPACE_H
8 
9 #include "pzcompel.h"
10 class TPZMaterialData;
11 
18 {
19 public:
20 
21  public:
22 virtual int ClassId() const override;
23 
24 
27 
29  virtual ~TPZInterpolationSpace();
30 
33 
35  TPZInterpolationSpace(TPZCompMesh &mesh, const TPZInterpolationSpace &copy, std::map<int64_t,int64_t> &gl2lcElMap);
36 
38  TPZInterpolationSpace(TPZCompMesh &mesh, const TPZInterpolationSpace &copy, int64_t &index);
39 
47  TPZInterpolationSpace(TPZCompMesh &mesh, TPZGeoEl *gel, int64_t &index);
48 
56  virtual void Print(std::ostream &out = std::cout) const override;
57 
59  int NSideShapeF(int side) const
60  {
61  return NShapeF();
62  }
63 
65  virtual int NSideConnects(int iside) const = 0;
66 
72  virtual int SideConnectLocId(int icon,int is) const = 0;
73 // {
74 // return icon;
75 // }
76 
78  int64_t SideConnectIndex(int icon,int is) const
79  {
80  int locid = SideConnectLocId(icon, is);
81  return ConnectIndex(locid);
82  }
83 
85  TPZConnect &SideConnect(int icon,int is) const
86  {
87  return Connect(SideConnectLocId(icon,is));
88  }
89 
91  virtual int NShapeF() const = 0;
92 
94  virtual int NConnectShapeF(int icon, int order) const = 0;
95 
97  virtual int MaxOrder();
98 
100  virtual void AdjustIntegrationRule();
101 
103  virtual int ComputeIntegrationOrder() const override;
104 
105 
106  virtual void SetIntegrationRule(int order) override{
107  std::cout << "TPZInterpolationSpace::SetIntegrationRule called\n";
108  }
109 
110 
121  virtual void Shape(TPZVec<REAL> &qsi,TPZFMatrix<REAL> &phi,TPZFMatrix<REAL> &dphidxi) = 0;
122 
123 
124 
125 
133  virtual void ComputeShape(TPZVec<REAL> &intpoint, TPZVec<REAL> &X, TPZFMatrix<REAL> &jacobian, TPZFMatrix<REAL> &axes,
134  REAL &detjac, TPZFMatrix<REAL> &jacinv, TPZFMatrix<REAL> &phi, TPZFMatrix<REAL> &dphi, TPZFMatrix<REAL> &dphidx);
135 
137  static void Convert2Axes(const TPZFMatrix<REAL> &dphi, const TPZFMatrix<REAL> &jacinv, TPZFMatrix<REAL> &dphidx);
138 
139 
141  virtual void SideShapeFunction(int side, TPZVec<REAL> &point, TPZFMatrix<REAL> &phi, TPZFMatrix<REAL> &dphi)
142  {
143  TPZGeoEl *gel = Reference();
144  TPZManVector<REAL,3> ptout(gel->Dimension());
145  Reference()->ProjectPoint(side, point, gel->NSides()-1, ptout);
146  Shape(ptout, phi, dphi);
147  }
148 
152  public:
153 
154 
160  virtual void ComputeShape(TPZVec<REAL> &intpoint, TPZMaterialData &data);
161 
166  virtual void InitMaterialData(TPZMaterialData &data);
167 
169  virtual void ComputeRequiredData(TPZMaterialData &data,
170  TPZVec<REAL> &qsi);
171 
173  virtual void ComputeRequiredData(TPZVec<REAL> &intpointtemp, TPZVec<TPZTransform<REAL> > &trvec, TPZVec<TPZMaterialData> &datavec)
174  {
175  PZError << "This Should never be called in this class, only in its children" << std::endl;
176  DebugStop();
177  }
178 
179 
181  virtual void ComputeNormal(TPZMaterialData & data);
182 
184  void VectorialProd(TPZVec<REAL> & ivec, TPZVec<REAL> & jvec, TPZVec<REAL> & kvec, bool unitary = false);
185 
191  virtual void CalcStiff(TPZElementMatrix &ek, TPZElementMatrix &ef) override;
192 
197  virtual void CalcResidual(TPZElementMatrix &ef) override;
198 
201 
203  virtual void InitializeElementMatrix(TPZElementMatrix &ef);
204 
211 
213  virtual const TPZIntPoints &GetIntegrationRule() const override = 0;
214 
216  virtual TPZIntPoints &GetIntegrationRule() = 0;
217 
219  virtual REAL InnerRadius();
220 
231  virtual void Solution(TPZVec<REAL> &qsi,int var,TPZVec<STATE> &sol) override;
232 
238 
244  void CreateInterfaces(bool BetweenContinuous = false);
245 
251  TPZInterfaceElement * CreateInterface(int side, bool BetweenContinuous = false);
252 
254  int ExistsInterface(TPZGeoElSide geosd);
255 
257  void RemoveInterfaces();
258 
260  void RemoveInterface(int side);
261 
268  virtual void EvaluateError( std::function<void(const TPZVec<REAL> &loc,TPZVec<STATE> &val,TPZFMatrix<STATE> &deriv)> func,
269  TPZVec<REAL> &errors, bool store_error ) override;
270 
272  virtual void ComputeError(int errorid, TPZVec<REAL> &error) override;
273 
275  virtual TPZVec<STATE> IntegrateSolution(int variable) const override;
276 
277  virtual void Integrate(int variable, TPZVec<STATE> & value) override;//AQUIFRAN
278 
280 // virtual void IntegrateSolution(TPZVec<STATE> & value);
281 
291  void ProjectFlux(TPZElementMatrix &ek, TPZElementMatrix &ef) override;
292 
293 protected:
294 
297 
298 public:
299 
301  virtual void SetPreferredOrder ( int order ) = 0;
302 
304  virtual int GetPreferredOrder () { return fPreferredOrder; }
305 
311  virtual void PRefine ( int order ) = 0;
312 
318  virtual int GetSideOrient(int side);
319 
325  virtual void SetSideOrient(int side, int sideorient);
326 
327 public:
328 
330  virtual void Write(TPZStream &buf, int withclassid) const override;
331 
333  virtual void Read(TPZStream &buf, void *context) override;
334 
346 
347 protected:
348 
362  void ExpandShapeFunctions(TPZVec<int64_t> &connectlist, TPZVec<int> &dependencyorder, TPZVec<int> &blocksizes, TPZFMatrix<REAL> &phi, TPZFMatrix<REAL> &dphi);
363 
364 };
365 
366 #endif
virtual ~TPZInterpolationSpace()
Default destructor.
virtual void ProjectPoint(int sidefrom, TPZVec< REAL > &ptin, int sideto, TPZVec< REAL > &ptout)
Project the point from one side to another. The dimension of the points needs to be configured proper...
Definition: pzgeoel.h:453
virtual int ClassId() const override
Define the class id associated with the class.
void CreateInterfaces(bool BetweenContinuous=false)
Create interfaces between this and its neighbours.
virtual int GetSideOrient(int side)
It returns the normal orientation of the reference element by the side. Only side that has dimension ...
Represents a set of shape functions associated with a computational element/side. Computational Eleme...
Definition: pzconnect.h:30
virtual void InitMaterialData(TPZMaterialData &data)
Initialize a material data and its attributes based on element dimension, number of state variables a...
int ExistsInterface(TPZGeoElSide geosd)
Verify existence of interface.
virtual int NConnectShapeF(int icon, int order) const =0
Returns the number of shapefunctions associated with a connect.
virtual void Integrate(int variable, TPZVec< STATE > &value) override
Integrates a variable over the element.
Contains declaration of TPZCompEl class which defines the interface of a computational element...
virtual int ComputeIntegrationOrder() const override
Compute integration order according to ... .
Utility class which represents an element with its side. The Geometric approximation classes Geometry...
Definition: pzgeoelside.h:83
virtual void CalcResidual(TPZElementMatrix &ef) override
Only computes the element residual.
REAL val(STATE &number)
Returns value of the variable.
Definition: pzartdiff.h:23
virtual int NSides() const =0
Returns the number of connectivities of the element.
virtual const TPZIntPoints & GetIntegrationRule() const override=0
Returns a reference to an integration rule suitable for integrating the interior of the element...
virtual void SideShapeFunction(int side, TPZVec< REAL > &point, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphi)
Compute the values of the shape function along the side.
virtual void Read(TPZStream &buf, void *context) override
Reads the element data from a stream.
TPZInterpolationSpace()
Default constructor.
Abstract class defining integration rules. Numerical Integration.
Definition: tpzintpoints.h:19
virtual void SetIntegrationRule(int order) override
void error(char *string)
Definition: testShape.cc:7
void BuildTransferMatrix(TPZInterpolationSpace &coarsel, TPZTransform<> &t, TPZTransfer< STATE > &transfer)
Accumulates the transfer coefficients between the current element and the coarse element into the t...
virtual void Shape(TPZVec< REAL > &qsi, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphidxi)=0
Computes the shape function set at the point x.
void ProjectFlux(TPZElementMatrix &ek, TPZElementMatrix &ef) override
Integrate the solution over the element.
void RemoveInterfaces()
Remove interfaces connected to this 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)
Compute shape functions based on master element in the classical FEM manner.
Defines the behaviour of all geometric elements. GeometryTPZGeoEl is the common denominator for all g...
Definition: pzgeoel.h:43
virtual int NSideConnects(int iside) const =0
Returns the number of dof nodes along side iside.
int64_t SideConnectIndex(int icon, int is) const
Returns the index of the c th connect object along side is.
TPZInterfaceElement * CreateInterface(int side, bool BetweenContinuous=false)
Create an interface between this and the neighbour by side side.
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
static void Convert2Axes(const TPZFMatrix< REAL > &dphi, const TPZFMatrix< REAL > &jacinv, TPZFMatrix< REAL > &dphidx)
convert a shapefunction derivative in xi-eta to a function derivative in axes
virtual void ComputeRequiredData(TPZVec< REAL > &intpointtemp, TPZVec< TPZTransform< REAL > > &trvec, TPZVec< TPZMaterialData > &datavec)
Compute and fill data with requested attributes for each of the compels in fElementVec.
virtual void ComputeError(int errorid, TPZVec< REAL > &error) override
Computes the element error estimator.
virtual void PRefine(int order)=0
Change the preferred order for the element and proceed the adjust of the aproximation space taking ...
TPZConnect & SideConnect(int icon, int is) const
Returns a pointer to the icon th connect object along side is.
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 void InitializeElementMatrix(TPZElementMatrix &ek, TPZElementMatrix &ef)
Initialize element matrix in which is computed CalcStiff.
virtual TPZConnect & Connect(int i) const
Returns a pointer to the ith node.
Definition: pzcompel.cpp:298
void VectorialProd(TPZVec< REAL > &ivec, TPZVec< REAL > &jvec, TPZVec< REAL > &kvec, bool unitary=false)
Computes the vectorial product of two vectors and normalize the result if unitary is set to true...
void InterpolateSolution(TPZInterpolationSpace &coarsel)
Interpolates the solution into the degrees of freedom nodes from the degrees of freedom nodes from th...
virtual int GetPreferredOrder()
Returns the prefered order for the element.
virtual int64_t ConnectIndex(int i) const =0
Returns the index of the ith connectivity of the element.
Computes the contribution over an interface between two discontinuous elements. Computational 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.
virtual void SetPreferredOrder(int order)=0
Defines the desired order for entire element.
virtual void SetSideOrient(int side, int sideorient)
It set the normal orientation of the element by the side. Only side that has dimension equal to my di...
virtual void Write(TPZStream &buf, int withclassid) const override
Saves the element data to a stream.
virtual int SideConnectLocId(int icon, int is) const =0
Returns the local node number of icon along is.
virtual TPZVec< STATE > IntegrateSolution(int variable) const override
Integrate a variable over the element.
virtual int Dimension() const =0
Returns the dimension of the element.
virtual void CalcStiff(TPZElementMatrix &ek, 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
virtual int NShapeF() const =0
It returns the shapes number of the element.
void ExpandShapeFunctions(TPZVec< int64_t > &connectlist, TPZVec< int > &dependencyorder, TPZVec< int > &blocksizes, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphi)
Auxiliary method to expand a vector of shapefunctions and their derivatives to acount for constraints...
Implements computational mesh. Computational Mesh.
Definition: pzcmesh.h:47
virtual int MaxOrder()
Returns the max order of interpolation.
virtual void ComputeNormal(TPZMaterialData &data)
Computes the proper normal vector towards the neighbour element.
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
void MinMaxSolutionValues(TPZVec< STATE > &min, TPZVec< STATE > &max)
Returns minimum and maximum values for each state variable.
Implements an affine transformation between points in parameter space. Topology Utility.
Definition: pzmganalysis.h:14
virtual void ComputeRequiredData(TPZMaterialData &data, TPZVec< REAL > &qsi)
Compute and fill data with requested attributes.
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
int NSideShapeF(int side) const
Returns the number of shape functions on a side.
int fPreferredOrder
Preferred polynomial order.
virtual void AdjustIntegrationRule()
Adjust the integration rule according to the polynomial order of shape functions. ...
Defines the interface of a computational element. Computational Element.
Definition: pzcompel.h:59
#define PZError
Defines the output device to error messages and the DebugStop() function.
Definition: pzerror.h:15
void RemoveInterface(int side)
Remove interface which is neighbour from side side.
virtual void Print(std::ostream &out=std::cout) const override
Prints the relevant data of the element to the output stream.
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.