13 #if !defined(KRATOS_LAPLACE_ELEMENT_H_INCLUDED)
14 #define KRATOS_LAPLACE_ELEMENT_H_INCLUDED
31 template <
unsigned int TDim,
unsigned int TNumNodes>
99 GeometryType::Pointer pGeometry)
109 GeometryType::Pointer pGeometry,
110 PropertiesType::Pointer pProperties)
111 :
Element(NewId, pGeometry, pProperties)
148 PropertiesType::Pointer pProperties)
const override
152 "LaplaceElement instances."
166 GeometryType::Pointer pGeom,
167 PropertiesType::Pointer pProperties)
const override
171 "LaplaceElement instances."
189 "LaplaceElement instances."
199 "GetVariable method. Please implement it in the "
213 const ProcessInfo& CurrentProcessInfo)
const override;
222 const ProcessInfo& CurrentProcessInfo)
const override;
226 int Step = 0)
const override;
295 std::string
Info()
const override
297 std::stringstream buffer;
298 buffer <<
"LaplaceElement #" << this->
Id();
304 rOStream <<
"LaplaceElement #" << this->
Id();
332 void save(
Serializer& rSerializer)
const override
357 template <
unsigned int TDim,
unsigned int TNumNodes>
362 template <
unsigned int TDim,
unsigned int TNumNodes>
367 rOStream <<
" : " << std::endl;
Base class for all Elements.
Definition: element.h:60
Vector VectorType
Definition: element.h:88
Element(IndexType NewId=0)
Definition: element.h:121
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
Node NodeType
definition of node type (default is: Node)
Definition: element.h:74
Matrix MatrixType
Definition: element.h:90
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: element.h:1135
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
std::size_t IndexType
Definition: flags.h:74
This defines the geometrical object, base definition of the element and condition entities.
Definition: geometrical_object.h:58
IntegrationMethod
Definition: geometry_data.h:76
Geometry base class.
Definition: geometry.h:71
PointerVector< TPointType > PointsArrayType
Definition: geometry.h:118
GeometryData::ShapeFunctionsGradientsType ShapeFunctionsGradientsType
Definition: geometry.h:189
IndexType Id() const
Definition: indexed_object.h:107
Definition: laplace_element.h:33
BaseType::MatrixType MatrixType
Matrix type for local contributions to the linear system.
Definition: laplace_element.h:55
LaplaceElement(IndexType NewId, GeometryType::Pointer pGeometry)
Definition: laplace_element.h:97
virtual void CalculateGeometryData(Vector &rGaussWeights, Matrix &rNContainer, ShapeFunctionDerivativesArrayType &rDN_DX) const
Calculates shape function data for this element.
Definition: laplace_element.cpp:218
LaplaceElement(IndexType NewId=0)
Definition: laplace_element.h:78
GeometryType::ShapeFunctionsGradientsType ShapeFunctionDerivativesArrayType
Type for an array of shape function gradient matrices.
Definition: laplace_element.h:64
LaplaceElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Definition: laplace_element.h:107
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Definition: laplace_element.h:145
std::string Info() const override
Turn back information as a string.
Definition: laplace_element.h:295
virtual void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: laplace_element.cpp:99
virtual const Variable< double > & GetVariable() const
Definition: laplace_element.h:194
LaplaceElement(LaplaceElement const &rOther)
Definition: laplace_element.h:118
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: laplace_element.h:302
void CalculateBoundedLeftHandSide(BoundedMatrix< double, TNumNodes, TNumNodes > &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) const
Definition: laplace_element.cpp:121
virtual void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: laplace_element.cpp:174
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &CurrentProcessInfo) const override
Definition: laplace_element.cpp:31
Element::Pointer Clone(IndexType NewId, NodesArrayType const &ThisNodes) const override
Definition: laplace_element.h:183
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(LaplaceElement)
BaseType::DofsVectorType DofsVectorType
Definition: laplace_element.h:59
void GetValuesVector(VectorType &rValues, int Step=0) const override
Definition: laplace_element.cpp:63
Element::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override
Definition: laplace_element.h:164
void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: laplace_element.cpp:155
void GetValuesArray(BoundedVector< double, TNumNodes > &rValues, int Step=0) const
Definition: laplace_element.cpp:78
void GetDofList(DofsVectorType &rElementalDofList, const ProcessInfo &CurrentProcessInfo) const override
Definition: laplace_element.cpp:47
BaseType::EquationIdVectorType EquationIdVectorType
Definition: laplace_element.h:57
GeometryData::IntegrationMethod GetIntegrationMethod() const override
Definition: laplace_element.cpp:93
LaplaceElement(IndexType NewId, const NodesArrayType &ThisNodes)
Definition: laplace_element.h:87
int Check(const ProcessInfo &rCurrentProcessInfo) const override
Definition: laplace_element.cpp:197
~LaplaceElement() override=default
This class defines the node.
Definition: node.h:65
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
Properties encapsulates data shared by different Elements or Conditions. It can store any type of dat...
Definition: properties.h:69
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
#define KRATOS_ERROR
Definition: exception.h:161
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
def load(f)
Definition: ode_solve.py:307