13 #if !defined(KRATOS_ELEMENT_H_INCLUDED )
14 #define KRATOS_ELEMENT_H_INCLUDED
123 , mpProperties(nullptr)
132 , mpProperties(nullptr)
141 , mpProperties(nullptr)
148 Element(
IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
150 , mpProperties(pProperties)
158 , mpProperties(rOther.mpProperties)
182 mpProperties = rOther.mpProperties;
203 PropertiesType::Pointer pProperties)
const
206 KRATOS_ERROR <<
"Please implement the First Create method in your derived Element" <<
Info() << std::endl;
207 return Kratos::make_intrusive<Element>(NewId,
GetGeometry().
Create(ThisNodes), pProperties);
219 GeometryType::Pointer pGeom,
220 PropertiesType::Pointer pProperties)
const
223 KRATOS_ERROR <<
"Please implement the Second Create method in your derived Element" <<
Info() << std::endl;
224 return Kratos::make_intrusive<Element>(NewId, pGeom, pProperties);
238 KRATOS_WARNING(
"Element") <<
" Call base class element Clone " << std::endl;
240 p_new_elem->SetData(this->
GetData());
241 p_new_elem->Set(
Flags(*
this));
261 if (rResult.size() != 0) {
274 if (rElementalDofList.size() != 0) {
275 rElementalDofList.resize(0);
409 if (rLeftHandSideMatrix.size1() != 0) {
410 rLeftHandSideMatrix.
resize(0, 0,
false);
412 if (rRightHandSideVector.size() != 0) {
413 rRightHandSideVector.
resize(0,
false);
426 if (rLeftHandSideMatrix.size1() != 0) {
427 rLeftHandSideMatrix.
resize(0, 0,
false);
440 if (rRightHandSideVector.size() != 0) {
441 rRightHandSideVector.
resize(0,
false);
465 if (rLeftHandSideMatrix.size1() != 0) {
466 rLeftHandSideMatrix.
resize(0, 0,
false);
468 if (rRightHandSideVector.size() != 0) {
469 rRightHandSideVector.
resize(0,
false);
482 if (rLeftHandSideMatrix.size1() != 0) {
483 rLeftHandSideMatrix.
resize(0, 0,
false);
496 if (rRightHandSideVector.size() != 0) {
497 rRightHandSideVector.
resize(0,
false);
522 if (rLeftHandSideMatrix.size1() != 0) {
523 rLeftHandSideMatrix.
resize(0, 0,
false);
525 if (rRightHandSideVector.size() != 0) {
526 rRightHandSideVector.
resize(0,
false);
539 if (rLeftHandSideMatrix.size1() != 0) {
540 rLeftHandSideMatrix.
resize(0, 0,
false);
553 if (rRightHandSideVector.size() != 0) {
554 rRightHandSideVector.
resize(0,
false);
572 if (rMassMatrix.size1() != 0) {
573 rMassMatrix.
resize(0, 0,
false);
585 if (rDampingMatrix.size1() != 0) {
586 rDampingMatrix.
resize(0, 0,
false);
601 KRATOS_ERROR <<
"Calling the CalculateLumpedMassVector() method of the base element. The method must be implemented in the derived element.";
641 KRATOS_ERROR <<
"Base element class is not able to assemble rRHS to the desired variable. destination variable is " << rDestinationVariable << std::endl;
659 KRATOS_ERROR <<
"Base element class is not able to assemble rRHS to the desired variable. destination variable is " << rDestinationVariable << std::endl;
677 KRATOS_ERROR <<
"Base element class is not able to assemble rLHS to the desired variable. destination variable is " << rDestinationVariable << std::endl;
721 std::vector<bool>& rOutput,
728 std::vector<int>& rOutput,
735 std::vector<double>& rOutput,
770 std::vector<Vector>& rOutput,
777 std::vector<Matrix>& rOutput,
784 std::vector<ConstitutiveLaw::Pointer>& rOutput,
827 const std::vector<bool>& rValues,
833 const std::vector<int>& rValues,
840 const std::vector<double>& rValues,
875 const std::vector<Vector>& rValues,
882 const std::vector<Matrix>& rValues,
889 const std::vector<ConstitutiveLaw::Pointer>& rValues,
911 KRATOS_ERROR_IF( domain_size <= 0.0 ) <<
"Element " << this->
Id() <<
" has non-positive size " << domain_size << std::endl;
928 if (rMassMatrix.size1() != 0) {
929 rMassMatrix.
resize(0, 0,
false);
952 if (rDampMatrix.size1() != 0) {
953 rDampMatrix.
resize(0, 0,
false);
975 if (rDampingMatrix.size1() != 0) {
976 rDampingMatrix.
resize(0, 0,
false);
987 if (rOutput.size1() != 0)
988 rOutput.
resize(0, 0,
false);
998 if (rOutput.size1() != 0)
999 rOutput.
resize(0, 0,
false);
1016 return mpProperties;
1021 return mpProperties;
1027 <<
"Tryining to get the properties of " <<
Info()
1028 <<
", which are uninitialized." << std::endl;
1029 return *mpProperties;
1035 <<
"Tryining to get the properties of " <<
Info()
1036 <<
", which are uninitialized." << std::endl;
1037 return *mpProperties;
1042 mpProperties = pProperties;
1052 return mpProperties !=
nullptr;
1090 "time_integration" : [],
1091 "framework" : "lagrangian",
1092 "symmetric_lhs" : false,
1093 "positive_definite_lhs" : false,
1096 "nodal_historical" : [],
1097 "nodal_non_historical" : [],
1100 "required_variables" : [],
1101 "required_dofs" : [],
1103 "compatible_geometries" : [],
1104 "element_integrates_in_time" : true,
1105 "compatible_constitutive_laws": {
1110 "required_polynomial_degree_of_geometry" : -1,
1111 "documentation" : "This is the base element"
1114 return specifications;
1121 std::stringstream buffer;
1122 buffer <<
"Element #" <<
Id();
1123 return buffer.str();
1130 rOStream <<
"Element #" <<
Id();
1178 Properties::Pointer mpProperties;
1192 void save(
Serializer& rSerializer)
const override
1195 rSerializer.
save(
"Properties", mpProperties);
1201 rSerializer.
load(
"Properties", mpProperties);
1234 rOStream <<
" : " << std::endl;
1248 #undef KRATOS_EXPORT_MACRO
1249 #define KRATOS_EXPORT_MACRO KRATOS_API
1253 #undef KRATOS_EXPORT_MACRO
1254 #define KRATOS_EXPORT_MACRO KRATOS_NO_EXPORT
Dof represents a degree of freedom (DoF).
Definition: dof.h:86
Base class for all Elements.
Definition: element.h:60
virtual void CalculateSecondDerivativesLHS(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:536
virtual void CalculateSensitivityMatrix(const Variable< double > &rDesignVariable, Matrix &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:983
virtual void ResetConstitutiveLaw()
Definition: element.h:349
virtual Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const
It creates a new element pointer.
Definition: element.h:218
virtual void CalculateOnIntegrationPoints(const Variable< bool > &rVariable, std::vector< bool > &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:719
Element(IndexType NewId, GeometryType::Pointer pGeometry)
Definition: element.h:139
Element(Element const &rOther)
Copy constructor.
Definition: element.h:156
virtual void CalculateOnIntegrationPoints(const Variable< array_1d< double, 4 >> &rVariable, std::vector< array_1d< double, 4 >> &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:747
virtual void CalculateFirstDerivativesLHS(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:479
virtual void FinalizeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:379
virtual void InitializeSolutionStep(const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:365
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: element.h:105
Element ElementType
definition of element type
Definition: element.h:68
virtual const Parameters GetSpecifications() const
This method provides the specifications/requirements of the element.
Definition: element.h:1087
PropertiesType & GetProperties()
Definition: element.h:1024
virtual void Calculate(const Variable< Matrix > &rVariable, Matrix &Output, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:704
virtual void SetValuesOnIntegrationPoints(const Variable< double > &rVariable, const std::vector< double > &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:838
virtual void AddExplicitContribution(const VectorType &rRHSVector, const Variable< VectorType > &rRHSVariable, const Variable< double > &rDestinationVariable, const ProcessInfo &rCurrentProcessInfo)
This function is designed to make the element to assemble an rRHS vector identified by a variable rRH...
Definition: element.h:634
Element(IndexType NewId, const NodesArrayType &ThisNodes)
Definition: element.h:130
virtual void CalculateSecondDerivativesContributions(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:518
virtual void SetValuesOnIntegrationPoints(const Variable< array_1d< double, 4 >> &rVariable, const std::vector< array_1d< double, 4 >> &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:852
Vector VectorType
Definition: element.h:88
virtual void SetValuesOnIntegrationPoints(const Variable< ConstitutiveLaw::Pointer > &rVariable, const std::vector< ConstitutiveLaw::Pointer > &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:887
virtual void DampMatrix(MatrixType &rDampMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:950
virtual void CalculateSecondDerivativesRHS(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:550
virtual void CalculateOnIntegrationPoints(const Variable< Matrix > &rVariable, std::vector< Matrix > &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:775
PointerVectorSet< DofType > DofsArrayType
Definition: element.h:102
virtual void CalculateSensitivityMatrix(const Variable< array_1d< double, 3 > > &rDesignVariable, Matrix &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:994
virtual void GetFirstDerivativesVector(Vector &values, int Step=0) const
Definition: element.h:310
virtual void CalculateOnIntegrationPoints(const Variable< Vector > &rVariable, std::vector< Vector > &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:768
std::size_t SizeType
Definition: element.h:94
virtual IntegrationMethod GetIntegrationMethod() const
Definition: element.h:285
virtual void AddMassMatrix(MatrixType &rLeftHandSideMatrix, double coeff, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:939
virtual void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:423
const PropertiesType::Pointer pGetProperties() const
Definition: element.h:1019
void SetProperties(PropertiesType::Pointer pProperties)
Definition: element.h:1040
virtual void CalculateOnIntegrationPoints(const Variable< double > &rVariable, std::vector< double > &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:733
virtual void SetValuesOnIntegrationPoints(const Variable< int > &rVariable, const std::vector< int > &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:831
virtual void SetValuesOnIntegrationPoints(const Variable< Matrix > &rVariable, const std::vector< Matrix > &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:880
virtual void SetValuesOnIntegrationPoints(const Variable< Vector > &rVariable, const std::vector< Vector > &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:873
Element(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Definition: element.h:148
virtual void CalculateOnIntegrationPoints(const Variable< int > &rVariable, std::vector< int > &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:726
virtual void CalculateMassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:570
Properties PropertiesType
Definition: element.h:80
Element(IndexType NewId=0)
Definition: element.h:121
virtual void AddExplicitContribution(const VectorType &rRHSVector, const Variable< VectorType > &rRHSVariable, const Variable< array_1d< double, 3 > > &rDestinationVariable, const ProcessInfo &rCurrentProcessInfo)
This function is designed to make the element to assemble an rRHS vector identified by a variable rRH...
Definition: element.h:652
virtual void Calculate(const Variable< Vector > &rVariable, Vector &Output, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:698
virtual void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:437
bool HasProperties() const
Check that the Element has a correctly initialized pointer to a Properties instance.
Definition: element.h:1050
virtual void CalculateDampingMatrix(MatrixType &rDampingMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:583
virtual void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const
Definition: element.h:258
virtual void AddExplicitContribution(const MatrixType &rLHSMatrix, const Variable< MatrixType > &rLHSVariable, const Variable< Matrix > &rDestinationVariable, const ProcessInfo &rCurrentProcessInfo)
This function is designed to make the element to assemble an rRHS vector identified by a variable rRH...
Definition: element.h:670
virtual void SetValuesOnIntegrationPoints(const Variable< array_1d< double, 9 >> &rVariable, const std::vector< array_1d< double, 9 >> &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:866
Element & operator=(Element const &rOther)
Assignment operator.
Definition: element.h:179
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
Geometry< NodeType >::PointsArrayType NodesArrayType
definition of nodes container type, redefined from GeometryType
Definition: element.h:86
virtual void SetValuesOnIntegrationPoints(const Variable< array_1d< double, 3 >> &rVariable, const std::vector< array_1d< double, 3 >> &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:845
virtual void CalculateLocalVelocityContribution(MatrixType &rDampingMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:972
virtual void CalculateFirstDerivativesContributions(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:461
virtual void CalculateLumpedMassVector(VectorType &rLumpedMassVector, const ProcessInfo &rCurrentProcessInfo) const
Definition: element.h:596
virtual Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const
It creates a new element pointer.
Definition: element.h:202
virtual void GetValuesVector(Vector &values, int Step=0) const
Definition: element.h:300
virtual void GetSecondDerivativesVector(Vector &values, int Step=0) const
Definition: element.h:320
Node NodeType
definition of node type (default is: Node)
Definition: element.h:74
~Element() override
Destructor.
Definition: element.h:164
PropertiesType const & GetProperties() const
Definition: element.h:1032
virtual void CalculateOnIntegrationPoints(const Variable< array_1d< double, 3 >> &rVariable, std::vector< array_1d< double, 3 >> &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:740
virtual void AddExplicitContribution(const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:622
std::string Info() const override
Turn back information as a string.
Definition: element.h:1119
virtual void InitializeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:372
virtual Pointer Clone(IndexType NewId, NodesArrayType const &ThisNodes) const
It creates a new element pointer and clones the previous element data.
Definition: element.h:235
virtual void CalculateOnIntegrationPoints(const Variable< array_1d< double, 6 >> &rVariable, std::vector< array_1d< double, 6 >> &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:754
virtual void CalculateOnIntegrationPoints(const Variable< array_1d< double, 9 >> &rVariable, std::vector< array_1d< double, 9 >> &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:761
virtual void Calculate(const Variable< array_1d< double, 3 > > &rVariable, array_1d< double, 3 > &Output, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:692
GeometryData GeometryDataType
Definition: element.h:107
Matrix MatrixType
Definition: element.h:90
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: element.h:1135
virtual void CalculateFirstDerivativesRHS(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:493
virtual void SetValuesOnIntegrationPoints(const Variable< bool > &rVariable, const std::vector< bool > &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:825
virtual void CalculateOnIntegrationPoints(const Variable< ConstitutiveLaw::Pointer > &rVariable, std::vector< ConstitutiveLaw::Pointer > &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:782
Dof< double > DofType
Definition: element.h:96
virtual void AddInertiaForces(VectorType &rRightHandSideVector, double coeff, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:961
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: element.h:1128
virtual void GetDofList(DofsVectorType &rElementalDofList, const ProcessInfo &rCurrentProcessInfo) const
Definition: element.h:271
virtual void SetValuesOnIntegrationPoints(const Variable< array_1d< double, 6 >> &rVariable, const std::vector< array_1d< double, 6 >> &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:859
virtual void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:405
GeometricalObject BaseType
base type: an GeometricalObject that automatically has a unique number
Definition: element.h:71
virtual void Calculate(const Variable< double > &rVariable, double &Output, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:686
virtual int Check(const ProcessInfo &rCurrentProcessInfo) const
Definition: element.h:904
Geometry< NodeType > GeometryType
definition of the geometry type with given NodeType
Definition: element.h:83
virtual void MassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:926
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(Element)
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
PropertiesType::Pointer pGetProperties()
returns the pointer to the property of the element. Does not throw an error, to allow copying of elem...
Definition: element.h:1014
virtual void Initialize(const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:341
virtual void FinalizeSolutionStep(const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:386
std::size_t IndexType
Definition: element.h:92
std::size_t IndexType
Definition: flags.h:74
Flags()
Default constructor.
Definition: flags.h:119
This defines the geometrical object, base definition of the element and condition entities.
Definition: geometrical_object.h:58
GeometryType::Pointer pGetGeometry()
Returns the pointer to the geometry.
Definition: geometrical_object.h:140
GeometryType & GetGeometry()
Returns the reference of the geometry.
Definition: geometrical_object.h:158
GeometricalObject & operator=(GeometricalObject const &rOther)
Assignment operator.
Definition: geometrical_object.h:112
DataValueContainer & GetData()
Definition: geometrical_object.h:212
Definition: geometry_data.h:60
IntegrationMethod
Definition: geometry_data.h:76
Geometry base class.
Definition: geometry.h:71
virtual int Check() const
Definition: geometry.h:3790
virtual double DomainSize() const
This method calculate and return length, area or volume of this geometry depending to it's dimension.
Definition: geometry.h:1371
IndexType Id() const
Definition: indexed_object.h:107
void resize(std::size_t NewSize1, std::size_t NewSize2, bool preserve=0)
Definition: amatrix_interface.h:224
KratosComponents class encapsulates a lookup table for a family of classes in a generic way.
Definition: kratos_components.h:49
This class defines the node.
Definition: node.h:65
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.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
void load(std::string const &rTag, TDataType &rObject)
Definition: serializer.h:207
void save(std::string const &rTag, std::array< TDataType, TDataSize > const &rObject)
Definition: serializer.h:545
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
#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
#define KRATOS_ERROR_IF(conditional)
Definition: exception.h:162
#define KRATOS_DEBUG_ERROR_IF(conditional)
Definition: exception.h:171
#define KRATOS_WARNING(label)
Definition: logger.h:265
#define KRATOS_API_EXTERN
Definition: kratos_export_api.h:57
#define KRATOS_API(...)
Definition: kratos_export_api.h:40
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
void AddKratosComponent(std::string const &Name, ExplicitBuilderType const &ThisComponent)
Definition: register_factories.cpp:23
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
KRATOS_DEFINE_VARIABLE(Vector, BIOT_STRAIN_VECTOR)
list coeff
Definition: bombardelli_test.py:41
list values
Definition: bombardelli_test.py:42
int domain_size
Definition: face_heat.py:4
def load(f)
Definition: ode_solve.py:307