NeoPZ
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | List of all members

Computes the contribution over an interface between two discontinuous elements. Computational Element. More...

#include <TPZInterfaceEl.h>

Inheritance diagram for TPZInterfaceElement:
[legend]
Collaboration diagram for TPZInterfaceElement:
[legend]

Public Types

enum  CalcStiffOptions {
  ENone, EStandard, EPenalty, EContDisc,
  EReferred
}
 

Public Member Functions

void GetConnects (TPZCompElSide &elside, TPZVec< TPZConnect *> &connects, TPZVec< int64_t > &connectindex)
 Extract connects from element el. More...
 
void NeighbourSolution (TPZCompElSide &Neighbor, TPZVec< REAL > &qsi, TPZSolVec &sol, TPZGradSolVec &dsol, TPZFMatrix< REAL > &NeighborAxes)
 Compute solution at neighbour element in a given master coordinate qsi. It returns the axes at which respect derivatives are computed. More...
 
bool CheckConsistencyOfMappedQsi (TPZCompElSide &Neighbor, TPZVec< REAL > &qsi, TPZVec< REAL > &NeighIntPoint)
 Check consistency of mapped qsi performed by method TPZInterfaceElement::MapQsi by comparing the X coordinate of qsi and the correspondent NeighIntPoint. More...
 
void ComputeSideTransform (TPZCompElSide &Neighbor, TPZTransform<> &transf)
 
void InitializeElementMatrix (TPZElementMatrix &ef)
 
void InitializeElementMatrix (TPZElementMatrix &ek, TPZElementMatrix &ef)
 
void MapQsi (TPZCompElSide &Neighbor, TPZVec< REAL > &qsi, TPZVec< REAL > &NeighIntPoint)
 Maps qsi coordinate at this master element to qsi coordinate at neighbor master element. More...
 
 TPZInterfaceElement (TPZCompMesh &mesh, TPZGeoEl *geo, int64_t &index, TPZCompElSide &left, TPZCompElSide &right)
 Constuctor to continuous and/or discontinuous neighbours. More...
 
 TPZInterfaceElement (TPZCompMesh &mesh, const TPZInterfaceElement &copy)
 Simple copy constructor. More...
 
 TPZInterfaceElement (TPZCompMesh &mesh, const TPZInterfaceElement &copy, std::map< int64_t, int64_t > &gl2lcConIdx, std::map< int64_t, int64_t > &gl2lcElIdx)
 Clone constructor to a patch mesh. More...
 
 TPZInterfaceElement (TPZCompMesh &mesh, const TPZInterfaceElement &copy, int64_t &index)
 Copy constructor with specified index. More...
 
 TPZInterfaceElement ()
 Empty constructor. More...
 
 TPZInterfaceElement (TPZCompMesh &mesh, TPZGeoEl *geo, int64_t &index)
 Default TPZCompEl constructor. SetLeftRightElements must be called before any computation. More...
 
 ~TPZInterfaceElement ()
 Destructor. More...
 
virtual int IsInterface () override
 
void SetLeftRightElements (TPZCompElSide &left, TPZCompElSide &right)
 Set neighbors. More...
 
virtual TPZCompElClone (TPZCompMesh &mesh) const override
 Makes a clone of this. More...
 
virtual TPZCompElClonePatchEl (TPZCompMesh &mesh, std::map< int64_t, int64_t > &gl2lcConMap, std::map< int64_t, int64_t > &gl2lcElMap) const override
 
TPZCompElCloneInterface (TPZCompMesh &aggmesh, int64_t &index, TPZCompElSide &left, TPZCompElSide &right) const
 Method used in TPZAgglomerateElement::CreateAgglomerateMesh. More...
 
void VolumeEls (TPZCompEl &thirdel)
 Identifies the elements of left and right volume of the interface. More...
 
TPZCompElRightElement () const
 Returns the right element from the element interface. More...
 
TPZCompElLeftElement () const
 Returns the left element from the element interface. More...
 
TPZCompElSideLeftElementSide ()
 Returns left neighbor. More...
 
TPZCompElSideRightElementSide ()
 Returns right neighbor. More...
 
void CenterNormal (TPZVec< REAL > &CenterNormal) const
 Returns the normal of this interface which goes from left to right neighbors. More...
 
void SetCenterNormal (const TPZVec< REAL > &CenterNormal)
 Set the normal to the given vector. More...
 
void Normal (TPZFMatrix< REAL > &axes, TPZVec< REAL > &normal)
 Returns normal based on already computed axes matrix. More...
 
void Normal (TPZVec< REAL > &qsi, TPZVec< REAL > &normal)
 Returns normal at qsi point. More...
 
virtual int NConnects () const override
 Returns the number from connectivities of the element. More...
 
int NRightConnects () const
 Returns the number from connectivities of the element related to right neighbour. More...
 
int NLeftConnects () const
 Returns the number from connectivities of the element related to left neighbour. More...
 
int64_t ConnectIndex (int i) const override
 Its return the connects of the left and right element associates. More...
 
void SetConnectIndex (int node, int64_t index) override
 This function should not be called. More...
 
virtual void BuildCornerConnectList (std::set< int64_t > &connectindexes) const override
 adds the connect indexes associated with base shape functions to the set More...
 
int Dimension () const override
 Returns the dimension from the element interface. More...
 
MElementType Type () override
 Type of the element. More...
 
virtual void LoadSolution () override
 Loads the solution within the internal data structure of the element. More...
 
virtual void CalcStiff (TPZElementMatrix &ek, TPZElementMatrix &ef) override
 CalcStiff computes the element stiffness matrix and right hand side. More...
 
virtual void CalcResidual (TPZElementMatrix &ef) override
 CalcResidual only computes the element residual. More...
 
virtual void ComputeSolution (TPZVec< REAL > &qsi, TPZVec< REAL > &normal, TPZSolVec &leftsol, TPZGradSolVec &dleftsol, TPZFMatrix< REAL > &leftaxes, TPZSolVec &rightsol, TPZGradSolVec &drightsol, TPZFMatrix< REAL > &rightaxes) override
 Computes solution and its derivatives in the local coordinate qsi. More...
 
virtual void ComputeSolution (TPZVec< REAL > &qsi, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphix, const TPZFMatrix< REAL > &axes, TPZSolVec &sol, TPZGradSolVec &dsol) override
 Computes solution and its derivatives in local coordinate qsi. More...
 
virtual void ComputeSolution (TPZVec< REAL > &qsi, TPZSolVec &sol, TPZGradSolVec &dsol, TPZFMatrix< REAL > &axes) override
 Computes solution and its derivatives in the local coordinate qsi. More...
 
virtual void ComputeSolution (TPZVec< REAL > &qsi, TPZMaterialData &data) override
 Computes solution and its derivatives in the local coordinate qsi. More...
 
void VetorialProd (TPZVec< REAL > &ivet, TPZVec< REAL > &jvet, TPZVec< REAL > &kvet)
 
void Print (std::ostream &out=std::cout) const override
 Prints attributes of the object. More...
 
virtual void CreateGraphicalElement (TPZGraphMesh &graphmesh, int dimension) override
 Interface elements does not have graphical representation. More...
 
void CloneInterface (TPZCompMesh *aggmesh)
 Make a clone of the fine mesh into clustered mesh. More...
 
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. More...
 
virtual void ComputeErrorFace (int errorid, TPZVec< STATE > &errorL, TPZVec< STATE > &errorR)
 ComputeError computes the element error estimator. More...
 
virtual void Integrate (int variable, TPZVec< STATE > &value) override
 Integrate a variable over the element. More...
 
void IntegrateInterface (int variable, TPZVec< REAL > &value)
 
void EvaluateInterfaceJump (TPZSolVec &jump, int opt)
 
int ComputeIntegrationOrder () const override
 Returns the unique identifier for reading/writing objects to streams. More...
 
virtual int ClassId () const override
 
virtual void Write (TPZStream &buf, int withclassid) const override
 Saves the element data to a stream. More...
 
virtual void Read (TPZStream &buf, void *context) override
 Reads the element data from a stream. More...
 
- Public Member Functions inherited from TPZCompEl
 TPZCompEl ()
 Simple Constructor. More...
 
virtual ~TPZCompEl ()
 Simple destructor. More...
 
 TPZCompEl (TPZCompMesh &mesh, const TPZCompEl &copy)
 Put a copy of the element in the referred mesh. More...
 
 TPZCompEl (TPZCompMesh &mesh, const TPZCompEl &copy, std::map< int64_t, int64_t > &gl2lcElMap)
 Put a copy of the element in the patch mesh. More...
 
 TPZCompEl (TPZCompMesh &mesh, const TPZCompEl &copy, int64_t &index)
 Copy of the element in the new mesh returning allocated index. More...
 
 TPZCompEl (TPZCompMesh &mesh, TPZGeoEl *gel, int64_t &index)
 Creates a computational element within mesh. Inserts the element within the data structure of the mesh. More...
 
virtual void SetCreateFunctions (TPZCompMesh *mesh)
 Sets create function in TPZCompMesh to create elements of this type. More...
 
virtual REAL VolumeOfEl ()
 Returns the volume of the geometric element associated. More...
 
virtual void LoadElementReference ()
 Loads the geometric element reference. More...
 
virtual REAL CompareElement (int var, char *matname)
 This method computes the norm of the difference of a post processed variable with @ the same post processed variable of the element pointed to by the geometric element. More...
 
virtual void Assemble ()
 Computes the element stifness matrix and right hand side in an internal data structure. Used for initializing condensed element data structures. More...
 
virtual bool HasMaterial (const std::set< int > &materialids) const
 Verifies if the material associated with the element is contained in the set. More...
 
virtual void Divide (int64_t index, TPZVec< int64_t > &subindex, int interpolate=0)
 Divide the computational element. If interpolate = 1, the solution is interpolated to the sub elements. More...
 
virtual void ProjectFlux (TPZElementMatrix &ek, TPZElementMatrix &ef)
 Projects the flux function on the finite element space. More...
 
virtual void ComputeError (int errorid, TPZVec< REAL > &error)
 ComputeError computes the element error estimator. More...
 
virtual void GetMemoryIndices (TPZVec< int64_t > &indices) const
 Get the indices of the vector of element memory associated with the integration points. More...
 
virtual void SetMemoryIndices (TPZVec< int64_t > &indices)
 Set the indices of the vector of element memory associated with the integration points. More...
 
virtual void PrepareIntPtIndices ()
 Prepare the vector of the material withmem with the correct integration point indexes. More...
 
virtual void ForcePrepareIntPtIndices ()
 PrepareIntPtIndices initializes the material damage varibles memory in the proper material class. More...
 
virtual void SetFreeIntPtIndices ()
 Frees the material damage varibles memory in the proper material class. More...
 
virtual int NumberOfCompElementsInsideThisCompEl ()
 Return the size of the elementvec in multiphysics, if it is not multiphysics, just return 1. More...
 
virtual void TransferMultiphysicsElementSolution ()
 
virtual void SetMultiphysicsElementSolution ()
 
virtual void AddShapeRestraint (TPZOneShapeRestraint restraint)
 Add a shape restraint (meant to fit the pyramid to restraint. More...
 
virtual std::list< TPZOneShapeRestraintGetShapeRestraints () const
 Return a list with the shape restraints. More...
 
virtual void ResetShapeRestraints ()
 Return a list with the shape restraints. More...
 
virtual void Solution (TPZVec< REAL > &qsi, int var, TPZVec< STATE > &sol)
 Calculates the solution - sol - for the variable var at point qsi, where qsi is expressed in terms of the master element coordinates. More...
 
virtual TPZVec< STATE > IntegrateSolution (int var) const
 Compute the integral of a variable. More...
 
virtual TPZVec< STATE > IntegrateSolution (const std::string &varname, const std::set< int > &matids)
 Compute the integral of a variable defined by the string if the material id is included in matids. More...
 
virtual void BuildConnectList (std::set< int64_t > &indepconnectlist, std::set< int64_t > &depconnectlist)
 Builds the list of all connectivities related to the element including the connects pointed to by dependent connects. More...
 
virtual void BuildConnectList (TPZStack< int64_t > &connectlist) const
 Builds the list of all connectivities related to the element including the connects pointed to by dependent connects. More...
 
virtual void BuildConnectList (std::set< int64_t > &connectlist)
 Builds the list of all connectivities related to the element including the connects pointed to by dependent connects. More...
 
virtual int HasDependency ()
 Returns 1 if the element has at least one dependent node. Returns 0 otherwise. More...
 
virtual int PressureConnectIndex () const
 Returns the index of the pressure connect. More...
 
virtual void ReduceInternalNodes ()
 Domain Decomposition.
This method will eliminate the nodes which are internal to the element from the datastructure of the grid
After calling this method, the superelement will statically condense the internal equations. More...
 
virtual void CalcBlockDiagonal (TPZStack< int64_t > &connectlist, TPZBlockDiagonal< STATE > &block)
 Calculates the diagonal block. More...
 
REAL MaximumRadiusOfEl ()
 Will return the maximum distance between the nodes of the reference element. More...
 
REAL LesserEdgeOfEl ()
 Will return the smallest distance between two nodes of the reference element. More...
 
virtual void SetIntegrationRule (TPZIntPoints *intrule)
 Method to set a dynamically allocated integration rule. More...
 
virtual void SetIntegrationRule (int order)
 
virtual const TPZIntPointsGetIntegrationRule () const
 
TPZGeoElReference () const
 Return a pointer to the corresponding geometric element if such exists, return 0 otherwise. More...
 
void SetReference (int64_t referenceindex)
 
virtual bool NeedsComputing (const std::set< int > &materialids)
 return true if the element has a variational statement associated with the material ids More...
 
virtual int NEquations ()
 Returns the number of equations of the element. More...
 
int64_t Index () const
 Returns element index of the mesh fELementVec list. More...
 
void SetIndex (int64_t index)
 Sets element index of the mesh fELementVec list. More...
 
virtual TPZConnectConnect (int i) const
 Returns a pointer to the ith node. More...
 
virtual TPZMaterialMaterial () const
 Identify the material object associated with the element. More...
 
TPZGeoElGetRefElPatch ()
 Returns the reference geometric element patch.
Look for a geometric element which refers to a computational element and is neighbour of the current element AND is larger than the current element. More...
 
void SetMesh (TPZCompMesh *mesh)
 Sets the grid of the element. More...
 
TPZCompMeshMesh () const
 Return a pointer to the grid of the element. More...
 
virtual void PrintSolution (TPZVec< REAL > &point, const char *VarName, std::ostream &out)
 Prints the solution - sol - for the variable "VarName" at point specified in terms of the master element coordinates. More...
 
virtual void PrintCoordinate (TPZVec< REAL > &point, int CoordinateIndex, std::ostream &out)
 Prints one coordinate index corresponding to the point to the output stream. More...
 
virtual void PrintTitle (const char *VarName, std::ostream &out)
 Prints the variables names associated with the element material. More...
 
- Public Member Functions inherited from TPZSavable
 TPZSavable ()
 
virtual ~TPZSavable ()
 
virtual std::list< std::map< std::string, uint64_t > > VersionHistory () const
 
virtual std::pair< std::string, uint64_t > Version () const
 
virtual bool Compare (TPZSavable *copy, bool override=false)
 Compares the object for identity with the object pointed to, eventually copy the object. More...
 
virtual bool Compare (TPZSavable *copy, bool override=false) const
 Compares the object for identity with the object pointed to, eventually copy the object. More...
 
- Public Member Functions inherited from TPZRegisterClassId
template<typename T >
 TPZRegisterClassId (int(T::*)() const)
 
 TPZRegisterClassId ()=default
 

Static Public Member Functions

static int ExistInterfaces (TPZCompElSide &comp)
 Verifies the existence of interfaces associates with the side of an element. More...
 
static int FreeInterface (TPZCompMesh &cmesh)
 
static int main (TPZCompMesh &cmesh)
 
- Static Public Member Functions inherited from TPZCompEl
static int StaticClassId ()
 
static void SetgOrder (int order)
 Sets the value of the default interpolation order. More...
 
static int GetgOrder ()
 Set the default value of the interpolation order. More...
 
static void SetOrthogonalFunction (void(*orthogonal)(REAL x, int num, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphi))
 Sets the orthogonal function which will be used throughout the program by default this function is the Chebyshev function. More...
 
static void Chebyshev (REAL x, int num, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphi)
 Implements of the orthogonal Chebyshev functions. More...
 
- Static Public Member Functions inherited from TPZSavable
static std::set< TPZRestoreClassBase * > & RestoreClassSet ()
 This static function guarantees that the gMap object is available when needed. More...
 
static std::map< int, TPZRestore_t > & ClassIdMap ()
 This static function guarantees that the gMap object is available when needed. More...
 
static std::pair< std::string, uint64_t > NeoPZVersion ()
 
static void Register (TPZRestoreClassBase *restore)
 
static void RegisterClassId (int classid, TPZRestore_t fun)
 
static TPZSavableCreateInstance (const int &classId)
 

Protected Member Functions

void InitMaterialData (TPZMaterialData &data, TPZInterpolationSpace *left)
 Initialize a material data and its attributes based on element dimension, number of state variables and material definitions. More...
 
void InitMaterialData (TPZMaterialData &data)
 Initialize the material data with the geometric data of the interface element. More...
 
void ComputeRequiredData (TPZMaterialData &data, TPZInterpolationSpace *elem, TPZVec< REAL > &IntPoint)
 Compute and fill data with requested attributes for neighbouring element. More...
 
virtual void ComputeRequiredData (TPZMaterialData &data, TPZVec< REAL > &qsi)
 Compute the required geometric data for the interface element. More...
 
virtual void ComputeRequiredData (TPZVec< REAL > &intpointtemp, TPZVec< TPZTransform<> > &trvec, TPZVec< TPZMaterialData > &datavec)
 
void ComputeNormal (TPZVec< REAL > &qsi, TPZVec< REAL > &normal)
 Computes normal at qsi point. More...
 
void ComputeNormal (TPZFMatrix< REAL > &axes, TPZVec< REAL > &normal)
 Computes normal based on already computed axes matrix. More...
 
void ComputeCenterNormal (TPZVec< REAL > &normal)
 Computes normal for linear geometric elements. More...
 

Private Member Functions

virtual void InitializeIntegrationRule () override
 Pointer to the integration rule. More...
 
void DecreaseElConnected ()
 Informs the connects that this element is no longer connected to it. More...
 
void IncrementElConnected ()
 Informs the connects that this element is connected to it. More...
 

Private Attributes

TPZCompElSide fLeftElSide
 Element the left of the normal a interface. More...
 
TPZCompElSide fRightElSide
 Element the right of the normal a interface. More...
 
TPZManVector< REAL, 3 > fCenterNormal
 Normal to the face element. More...
 

Additional Inherited Members

- Protected Attributes inherited from TPZCompEl
TPZCompMeshfMesh
 Computational mesh to which the element belongs. More...
 
int64_t fIndex
 Element index into mesh element vector. More...
 
TPZIntPointsfIntegrationRule
 Integration rule established by the user. More...
 

Detailed Description

Computes the contribution over an interface between two discontinuous elements. Computational Element.

Definition at line 30 of file TPZInterfaceEl.h.

Member Enumeration Documentation

◆ CalcStiffOptions

Enumerator
ENone 
EStandard 
EPenalty 
EContDisc 
EReferred 

Definition at line 127 of file TPZInterfaceEl.h.

Constructor & Destructor Documentation

◆ TPZInterfaceElement() [1/6]

TPZInterfaceElement::TPZInterfaceElement ( TPZCompMesh mesh,
TPZGeoEl geo,
int64_t &  index,
TPZCompElSide left,
TPZCompElSide right 
)

◆ TPZInterfaceElement() [2/6]

TPZInterfaceElement::TPZInterfaceElement ( TPZCompMesh mesh,
const TPZInterfaceElement copy 
)

◆ TPZInterfaceElement() [3/6]

TPZInterfaceElement::TPZInterfaceElement ( TPZCompMesh mesh,
const TPZInterfaceElement copy,
std::map< int64_t, int64_t > &  gl2lcConIdx,
std::map< int64_t, int64_t > &  gl2lcElIdx 
)

Clone constructor to a patch mesh.

Parameters
meshreference to a clone mesh
copyelement to be copied
gl2lcConIdxmap with connects
gl2lcElIdxmap with computational elements

Definition at line 187 of file TPZInterfaceEl.cpp.

References TPZIntPoints::Clone(), DebugStop, TPZCompElSide::Element(), TPZCompMesh::ElementVec(), fCenterNormal, TPZCompEl::fIntegrationRule, fLeftElSide, fRightElSide, IncrementElConnected(), TPZGeoEl::IncrementNumInterfaces(), TPZCompEl::Index(), LOGPZ_ERROR, TPZCompEl::Mesh(), PZError, TPZCompEl::Reference(), TPZCompElSide::SetElement(), TPZCompElSide::SetSide(), and TPZCompElSide::Side().

◆ TPZInterfaceElement() [4/6]

TPZInterfaceElement::TPZInterfaceElement ( TPZCompMesh mesh,
const TPZInterfaceElement copy,
int64_t &  index 
)

◆ TPZInterfaceElement() [5/6]

TPZInterfaceElement::TPZInterfaceElement ( )

Empty constructor.

Definition at line 282 of file TPZInterfaceEl.cpp.

Referenced by Clone(), CloneInterface(), and ClonePatchEl().

◆ TPZInterfaceElement() [6/6]

TPZInterfaceElement::TPZInterfaceElement ( TPZCompMesh mesh,
TPZGeoEl geo,
int64_t &  index 
)

Default TPZCompEl constructor. SetLeftRightElements must be called before any computation.

Definition at line 139 of file TPZInterfaceEl.cpp.

References IncrementElConnected(), TPZGeoEl::IncrementNumInterfaces(), and TPZGeoEl::SetReference().

◆ ~TPZInterfaceElement()

TPZInterfaceElement::~TPZInterfaceElement ( )

Member Function Documentation

◆ BuildCornerConnectList()

void TPZInterfaceElement::BuildCornerConnectList ( std::set< int64_t > &  connectindexes) const
overridevirtual

adds the connect indexes associated with base shape functions to the set

Implements TPZCompEl.

Definition at line 1698 of file TPZInterfaceEl.cpp.

References TPZCompEl::BuildCornerConnectList(), DebugStop, LeftElement(), and RightElement().

Referenced by SetCenterNormal().

◆ CalcResidual()

void TPZInterfaceElement::CalcResidual ( TPZElementMatrix ef)
overridevirtual

◆ CalcStiff()

void TPZInterfaceElement::CalcStiff ( TPZElementMatrix ek,
TPZElementMatrix ef 
)
overridevirtual

◆ CenterNormal()

void TPZInterfaceElement::CenterNormal ( TPZVec< REAL > &  CenterNormal) const

Returns the normal of this interface which goes from left to right neighbors.

Definition at line 756 of file TPZInterfaceEl.cpp.

References fCenterNormal, TPZVec< T >::NElements(), and TPZVec< T >::Resize().

Referenced by Normal(), RightElementSide(), and SetCenterNormal().

◆ CheckConsistencyOfMappedQsi()

bool TPZInterfaceElement::CheckConsistencyOfMappedQsi ( TPZCompElSide Neighbor,
TPZVec< REAL > &  qsi,
TPZVec< REAL > &  NeighIntPoint 
)

Check consistency of mapped qsi performed by method TPZInterfaceElement::MapQsi by comparing the X coordinate of qsi and the correspondent NeighIntPoint.

It return true if everything is ok or false otherwise

Definition at line 1546 of file TPZInterfaceEl.cpp.

References TPZCompElSide::Element(), TPZVec< T >::NElements(), Print(), PZError, TPZCompEl::Reference(), sqrt, pzgeom::tol, TPZGeoEl::X(), and ZeroTolerance().

Referenced by CalcResidual(), CalcStiff(), ComputeErrorFace(), ComputeRequiredData(), and MapQsi().

◆ ClassId()

int TPZInterfaceElement::ClassId ( ) const
overridevirtual

returns the unique identifier for reading/writing objects to streams

Reimplemented from TPZCompEl.

Definition at line 790 of file TPZInterfaceEl.cpp.

References TPZCompEl::ClassId(), and Hash().

Referenced by CreateGraphicalElement().

◆ Clone()

virtual TPZCompEl* TPZInterfaceElement::Clone ( TPZCompMesh mesh) const
inlineoverridevirtual

Makes a clone of this.

Implements TPZCompEl.

Definition at line 166 of file TPZInterfaceEl.h.

References TPZInterfaceElement().

◆ CloneInterface() [1/2]

TPZCompEl * TPZInterfaceElement::CloneInterface ( TPZCompMesh aggmesh,
int64_t &  index,
TPZCompElSide left,
TPZCompElSide right 
) const

◆ CloneInterface() [2/2]

void TPZInterfaceElement::CloneInterface ( TPZCompMesh aggmesh)

Make a clone of the fine mesh into clustered mesh.

◆ ClonePatchEl()

virtual TPZCompEl* TPZInterfaceElement::ClonePatchEl ( TPZCompMesh mesh,
std::map< int64_t, int64_t > &  gl2lcConMap,
std::map< int64_t, int64_t > &  gl2lcElMap 
) const
inlineoverridevirtual
See also
class TPZCompEl

Implements TPZCompEl.

Definition at line 171 of file TPZInterfaceEl.h.

References CloneInterface(), TPZInterfaceElement(), and VolumeEls().

◆ ComputeCenterNormal()

void TPZInterfaceElement::ComputeCenterNormal ( TPZVec< REAL > &  normal)
protected

Computes normal for linear geometric elements.

For linear geometry the normal vector is constant.

Definition at line 611 of file TPZInterfaceEl.cpp.

References TPZGeoEl::CenterPoint(), ComputeNormal(), Dimension(), TPZGeoEl::NSides(), and TPZCompEl::Reference().

Referenced by ComputeRequiredData().

◆ ComputeErrorFace()

void TPZInterfaceElement::ComputeErrorFace ( int  errorid,
TPZVec< STATE > &  errorL,
TPZVec< STATE > &  errorR 
)
virtual

◆ ComputeIntegrationOrder()

int TPZInterfaceElement::ComputeIntegrationOrder ( ) const
overridevirtual

Returns the unique identifier for reading/writing objects to streams.

Reimplemented from TPZCompEl.

Definition at line 403 of file TPZInterfaceEl.cpp.

References LeftElement(), TPZInterpolationSpace::MaxOrder(), and RightElement().

Referenced by CreateGraphicalElement().

◆ ComputeNormal() [1/2]

void TPZInterfaceElement::ComputeNormal ( TPZVec< REAL > &  qsi,
TPZVec< REAL > &  normal 
)
protected

Computes normal at qsi point.

Definition at line 618 of file TPZInterfaceEl.cpp.

References TPZGeoEl::Dimension(), TPZGeoEl::Jacobian(), and TPZCompEl::Reference().

Referenced by ComputeCenterNormal(), ComputeRequiredData(), and Normal().

◆ ComputeNormal() [2/2]

void TPZInterfaceElement::ComputeNormal ( TPZFMatrix< REAL > &  axes,
TPZVec< REAL > &  normal 
)
protected

◆ ComputeRequiredData() [1/3]

void TPZInterfaceElement::ComputeRequiredData ( TPZMaterialData data,
TPZInterpolationSpace elem,
TPZVec< REAL > &  IntPoint 
)
protected

◆ ComputeRequiredData() [2/3]

void TPZInterfaceElement::ComputeRequiredData ( TPZMaterialData data,
TPZVec< REAL > &  qsi 
)
protectedvirtual

◆ ComputeRequiredData() [3/3]

virtual void TPZInterfaceElement::ComputeRequiredData ( TPZVec< REAL > &  intpointtemp,
TPZVec< TPZTransform<> > &  trvec,
TPZVec< TPZMaterialData > &  datavec 
)
inlineprotectedvirtual

◆ ComputeSideTransform()

void TPZInterfaceElement::ComputeSideTransform ( TPZCompElSide Neighbor,
TPZTransform<> &  transf 
)

◆ ComputeSolution() [1/4]

void TPZInterfaceElement::ComputeSolution ( TPZVec< REAL > &  qsi,
TPZVec< REAL > &  normal,
TPZSolVec leftsol,
TPZGradSolVec dleftsol,
TPZFMatrix< REAL > &  leftaxes,
TPZSolVec rightsol,
TPZGradSolVec drightsol,
TPZFMatrix< REAL > &  rightaxes 
)
overridevirtual

Computes solution and its derivatives in the local coordinate qsi.

Parameters
qsi[in] master element coordinate
normalnormal vector
leftsol[out] left finite element solution
rightsol[out] right finite element solution
dleftsol[out] left solution derivatives
drightsol[out] right solution derivatives
leftaxes[out] axes associated with the derivative of the left element
rightaxes[out] axes associated with the derivative of the right element

Reimplemented from TPZCompEl.

Definition at line 1608 of file TPZInterfaceEl.cpp.

References fLeftElSide, fRightElSide, NeighbourSolution(), and Normal().

Referenced by LoadSolution().

◆ ComputeSolution() [2/4]

void TPZInterfaceElement::ComputeSolution ( TPZVec< REAL > &  qsi,
TPZFMatrix< REAL > &  phi,
TPZFMatrix< REAL > &  dphix,
const TPZFMatrix< REAL > &  axes,
TPZSolVec sol,
TPZGradSolVec dsol 
)
overridevirtual

Computes solution and its derivatives in local coordinate qsi.

Parameters
qsimaster element coordinate
phimatrix containing shape functions compute in qsi point
dphixmatrix containing the derivatives of shape functions in the direction of the axes
axesaxes indicating the direction of the derivatives
solfinite element solution
dsolsolution derivatives

Reimplemented from TPZCompEl.

Definition at line 1583 of file TPZInterfaceEl.cpp.

References TPZManVector< T, NumExtAlloc >::Resize().

◆ ComputeSolution() [3/4]

void TPZInterfaceElement::ComputeSolution ( TPZVec< REAL > &  qsi,
TPZSolVec sol,
TPZGradSolVec dsol,
TPZFMatrix< REAL > &  axes 
)
overridevirtual

Computes solution and its derivatives in the local coordinate qsi.

Parameters
qsimaster element coordinate
solfinite element solution
dsolsolution derivatives
axesaxes associated with the derivative of the solution

Reimplemented from TPZCompEl.

Definition at line 1589 of file TPZInterfaceEl.cpp.

References TPZManVector< T, NumExtAlloc >::Resize(), and TPZFMatrix< TVar >::Zero().

◆ ComputeSolution() [4/4]

void TPZInterfaceElement::ComputeSolution ( TPZVec< REAL > &  qsi,
TPZMaterialData data 
)
overridevirtual

Computes solution and its derivatives in the local coordinate qsi.

Parameters
qsimaster element coordinate
datacontains all elements to compute the solution

Reimplemented from TPZCompEl.

Definition at line 1601 of file TPZInterfaceEl.cpp.

References DebugStop.

◆ ConnectIndex()

int64_t TPZInterfaceElement::ConnectIndex ( int  i) const
overridevirtual

◆ CreateGraphicalElement()

virtual void TPZInterfaceElement::CreateGraphicalElement ( TPZGraphMesh graphmesh,
int  dimension 
)
inlineoverridevirtual

Interface elements does not have graphical representation.

See also
Base class for comments

Reimplemented from TPZCompEl.

Definition at line 323 of file TPZInterfaceEl.h.

References ClassId(), CloneInterface(), ComputeErrorFace(), ComputeIntegrationOrder(), EvaluateError(), EvaluateInterfaceJump(), ExistInterfaces(), FreeInterface(), Integrate(), IntegrateInterface(), main(), Read(), val(), and Write().

◆ DecreaseElConnected()

void TPZInterfaceElement::DecreaseElConnected ( )
private

Informs the connects that this element is no longer connected to it.

Definition at line 83 of file TPZInterfaceEl.cpp.

◆ Dimension()

int TPZInterfaceElement::Dimension ( ) const
inlineoverridevirtual

Returns the dimension from the element interface.

Implements TPZCompEl.

Definition at line 237 of file TPZInterfaceEl.h.

References TPZGeoEl::Dimension(), and TPZCompEl::Reference().

Referenced by CalcResidual(), CalcStiff(), ComputeCenterNormal(), ComputeErrorFace(), ComputeRequiredData(), ComputeSideTransform(), InitMaterialData(), and IntegrateInterface().

◆ EvaluateError()

void TPZInterfaceElement::EvaluateError ( std::function< void(const TPZVec< REAL > &loc, TPZVec< STATE > &val, TPZFMatrix< STATE > &deriv)>  func,
TPZVec< REAL > &  errors,
bool  store_error 
)
overridevirtual

Performs an error estimate on the elemen.

Parameters
fpfunction pointer which computes the exact solution
errors[out] the L2 norm of the error of the solution
flux[in] value of the interpolated flux values

Reimplemented from TPZCompEl.

Definition at line 780 of file TPZInterfaceEl.cpp.

References TPZVec< T >::Fill().

Referenced by CreateGraphicalElement().

◆ EvaluateInterfaceJump()

void TPZInterfaceElement::EvaluateInterfaceJump ( TPZSolVec jump,
int  opt 
)

◆ ExistInterfaces()

int TPZInterfaceElement::ExistInterfaces ( TPZCompElSide comp)
static

Verifies the existence of interfaces associates with the side of an element.

Case to interface should exist and exists only a returns 1,
case to interface should not exist and does not exist returns 1,
otherwise returns 0

Definition at line 529 of file TPZInterfaceEl.cpp.

References DebugStop, EInterface, TPZGeoElSide::Element(), TPZCompElSide::Element(), TPZGeoElSide::Exists(), TPZCompElSide::Exists(), TPZCompElSide::HigherLevelElementList(), TPZGeoElSide::Neighbour(), TPZVec< T >::NElements(), PZError, TPZGeoElSide::Reference(), TPZCompElSide::Reference(), TPZManVector< T, NumExtAlloc >::Resize(), and TPZCompEl::Type().

Referenced by CreateGraphicalElement(), and main().

◆ FreeInterface()

int TPZInterfaceElement::FreeInterface ( TPZCompMesh cmesh)
static

◆ GetConnects()

void TPZInterfaceElement::GetConnects ( TPZCompElSide elside,
TPZVec< TPZConnect *> &  connects,
TPZVec< int64_t > &  connectindex 
)

◆ IncrementElConnected()

void TPZInterfaceElement::IncrementElConnected ( )
private

Informs the connects that this element is connected to it.

Definition at line 91 of file TPZInterfaceEl.cpp.

Referenced by TPZInterfaceElement().

◆ InitializeElementMatrix() [1/2]

void TPZInterfaceElement::InitializeElementMatrix ( TPZElementMatrix ef)

◆ InitializeElementMatrix() [2/2]

void TPZInterfaceElement::InitializeElementMatrix ( TPZElementMatrix ek,
TPZElementMatrix ef 
)

◆ InitializeIntegrationRule()

void TPZInterfaceElement::InitializeIntegrationRule ( )
overrideprivatevirtual

◆ InitMaterialData() [1/2]

void TPZInterfaceElement::InitMaterialData ( TPZMaterialData data,
TPZInterpolationSpace left 
)
protected

◆ InitMaterialData() [2/2]

void TPZInterfaceElement::InitMaterialData ( TPZMaterialData data)
protected

◆ Integrate()

void TPZInterfaceElement::Integrate ( int  variable,
TPZVec< STATE > &  value 
)
overridevirtual

Integrate a variable over the element.

Reimplemented from TPZCompEl.

Definition at line 1430 of file TPZInterfaceEl.cpp.

References TPZVec< T >::Fill(), TPZCompEl::Material(), TPZMaterial::NSolutionVariables(), and TPZVec< T >::Resize().

Referenced by CreateGraphicalElement().

◆ IntegrateInterface()

void TPZInterfaceElement::IntegrateInterface ( int  variable,
TPZVec< REAL > &  value 
)

◆ IsInterface()

virtual int TPZInterfaceElement::IsInterface ( )
inlineoverridevirtual

Reimplemented from TPZCompEl.

Definition at line 159 of file TPZInterfaceEl.h.

References SetLeftRightElements().

◆ LeftElement()

TPZCompEl* TPZInterfaceElement::LeftElement ( ) const
inline

◆ LeftElementSide()

TPZCompElSide& TPZInterfaceElement::LeftElementSide ( )
inline

◆ LoadSolution()

virtual void TPZInterfaceElement::LoadSolution ( )
inlineoverridevirtual

Loads the solution within the internal data structure of the element.

Is used to initialize the solution of connect objects with dependency Is also used to load the solution within SuperElements

Reimplemented from TPZCompEl.

Definition at line 251 of file TPZInterfaceEl.h.

References CalcResidual(), CalcStiff(), ComputeSolution(), Print(), and VetorialProd().

◆ main()

int TPZInterfaceElement::main ( TPZCompMesh cmesh)
static

◆ MapQsi()

void TPZInterfaceElement::MapQsi ( TPZCompElSide Neighbor,
TPZVec< REAL > &  qsi,
TPZVec< REAL > &  NeighIntPoint 
)

Maps qsi coordinate at this master element to qsi coordinate at neighbor master element.

Parameters
Neighbor[in] may be this->LeftElementSide() or this->RightElementSide()
qsi[in] is the point at this element master
NeighIntPoint[out] is the point at neighbor element master. X[qsi] is equal to X[NeighIntPoint]

Definition at line 1537 of file TPZInterfaceEl.cpp.

References TPZTransform< T >::Apply(), CheckConsistencyOfMappedQsi(), and ComputeSideTransform().

Referenced by ComputeRequiredData(), and NeighbourSolution().

◆ NConnects()

int TPZInterfaceElement::NConnects ( ) const
overridevirtual

Returns the number from connectivities of the element.

Implements TPZCompEl.

Definition at line 387 of file TPZInterfaceEl.cpp.

References NLeftConnects(), and NRightConnects().

Referenced by EvaluateInterfaceJump(), and SetCenterNormal().

◆ NeighbourSolution()

void TPZInterfaceElement::NeighbourSolution ( TPZCompElSide Neighbor,
TPZVec< REAL > &  qsi,
TPZSolVec sol,
TPZGradSolVec dsol,
TPZFMatrix< REAL > &  NeighborAxes 
)

Compute solution at neighbour element in a given master coordinate qsi. It returns the axes at which respect derivatives are computed.

Parameters
[in]Neighbor
[in]qsi
[out]sol
[out]dsol
[out]NeighborAxes

Definition at line 1573 of file TPZInterfaceEl.cpp.

References TPZCompEl::ComputeSolution(), TPZGeoEl::Dimension(), TPZCompElSide::Element(), MapQsi(), and TPZCompEl::Reference().

Referenced by ComputeRequiredData(), ComputeSolution(), EvaluateInterfaceJump(), and IntegrateInterface().

◆ NLeftConnects()

int TPZInterfaceElement::NLeftConnects ( ) const

Returns the number from connectivities of the element related to left neighbour.

Definition at line 391 of file TPZInterfaceEl.cpp.

References TPZCompElSide::Element(), fLeftElSide, and TPZCompEl::NConnects().

Referenced by ConnectIndex(), NConnects(), and SetCenterNormal().

◆ Normal() [1/2]

void TPZInterfaceElement::Normal ( TPZFMatrix< REAL > &  axes,
TPZVec< REAL > &  normal 
)

Returns normal based on already computed axes matrix.

Axes has been computed in the desired qsi coordinate If geometric element has LinearMapping the CenterNormal is returned

Definition at line 762 of file TPZInterfaceEl.cpp.

References CenterNormal(), ComputeNormal(), TPZGeoEl::IsLinearMapping(), and TPZCompEl::Reference().

Referenced by CalcResidual(), ComputeErrorFace(), ComputeSolution(), IntegrateInterface(), and SetCenterNormal().

◆ Normal() [2/2]

void TPZInterfaceElement::Normal ( TPZVec< REAL > &  qsi,
TPZVec< REAL > &  normal 
)

Returns normal at qsi point.

If geometric element has LinearMapping the CenterNormal is returned

Definition at line 768 of file TPZInterfaceEl.cpp.

References CenterNormal(), ComputeNormal(), TPZGeoEl::IsLinearMapping(), and TPZCompEl::Reference().

◆ NRightConnects()

int TPZInterfaceElement::NRightConnects ( ) const

Returns the number from connectivities of the element related to right neighbour.

Definition at line 397 of file TPZInterfaceEl.cpp.

References TPZCompElSide::Element(), fRightElSide, and TPZCompEl::NConnects().

Referenced by ConnectIndex(), NConnects(), and SetCenterNormal().

◆ Print()

void TPZInterfaceElement::Print ( std::ostream &  out = std::cout) const
overridevirtual

◆ Read()

void TPZInterfaceElement::Read ( TPZStream buf,
void *  context 
)
overridevirtual

Reads the element data from a stream.

Read the element data from a stream

Reimplemented from TPZCompEl.

Definition at line 829 of file TPZInterfaceEl.cpp.

References fCenterNormal, fLeftElSide, fRightElSide, TPZPersistenceManager::GetInstance(), TPZStream::Read(), TPZCompEl::Read(), TPZCompElSide::SetElement(), and TPZCompElSide::SetSide().

Referenced by CreateGraphicalElement().

◆ RightElement()

TPZCompEl* TPZInterfaceElement::RightElement ( ) const
inline

◆ RightElementSide()

TPZCompElSide& TPZInterfaceElement::RightElementSide ( )
inline

◆ SetCenterNormal()

void TPZInterfaceElement::SetCenterNormal ( const TPZVec< REAL > &  CenterNormal)
inline

Set the normal to the given vector.

Definition at line 202 of file TPZInterfaceEl.h.

References BuildCornerConnectList(), CenterNormal(), ConnectIndex(), NConnects(), NLeftConnects(), Normal(), NRightConnects(), and SetConnectIndex().

◆ SetConnectIndex()

void TPZInterfaceElement::SetConnectIndex ( int  node,
int64_t  index 
)
overridevirtual

This function should not be called.

Implements TPZCompEl.

Definition at line 485 of file TPZInterfaceEl.cpp.

References DebugStop.

Referenced by SetCenterNormal().

◆ SetLeftRightElements()

void TPZInterfaceElement::SetLeftRightElements ( TPZCompElSide left,
TPZCompElSide right 
)

◆ Type()

MElementType TPZInterfaceElement::Type ( )
inlineoverridevirtual

Type of the element.

Reimplemented from TPZCompEl.

Definition at line 242 of file TPZInterfaceEl.h.

References EInterface.

◆ VetorialProd()

void TPZInterfaceElement::VetorialProd ( TPZVec< REAL > &  ivet,
TPZVec< REAL > &  jvet,
TPZVec< REAL > &  kvet 
)

Definition at line 748 of file TPZInterfaceEl.cpp.

References TPZVec< T >::Resize().

Referenced by ComputeNormal(), and LoadSolution().

◆ VolumeEls()

void TPZInterfaceElement::VolumeEls ( TPZCompEl thirdel)

Identifies the elements of left and right volume of the interface.

Referenced by ClonePatchEl().

◆ Write()

void TPZInterfaceElement::Write ( TPZStream buf,
int  withclassid 
) const
overridevirtual

Saves the element data to a stream.

Save the element data to a stream

Reimplemented from TPZCompEl.

Definition at line 802 of file TPZInterfaceEl.cpp.

References DebugStop, TPZCompElSide::Element(), fCenterNormal, fLeftElSide, fRightElSide, TPZCompEl::Index(), PZError, TPZCompElSide::Side(), TPZStream::Write(), TPZCompEl::Write(), and TPZPersistenceManager::WritePointer().

Referenced by CreateGraphicalElement().

Member Data Documentation

◆ fCenterNormal

TPZManVector<REAL,3> TPZInterfaceElement::fCenterNormal
private

Normal to the face element.

Definition at line 41 of file TPZInterfaceEl.h.

Referenced by CenterNormal(), InitMaterialData(), Print(), Read(), TPZInterfaceElement(), and Write().

◆ fLeftElSide

TPZCompElSide TPZInterfaceElement::fLeftElSide
private

Element the left of the normal a interface.

Definition at line 35 of file TPZInterfaceEl.h.

Referenced by ComputeSolution(), ConnectIndex(), LeftElementSide(), NLeftConnects(), Read(), TPZInterfaceElement(), and Write().

◆ fRightElSide

TPZCompElSide TPZInterfaceElement::fRightElSide
private

Element the right of the normal a interface.

Definition at line 38 of file TPZInterfaceEl.h.

Referenced by ComputeSolution(), ConnectIndex(), NRightConnects(), Read(), RightElementSide(), TPZInterfaceElement(), and Write().


The documentation for this class was generated from the following files: