KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
element.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ \.
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Pooyan Dadvand
11 //
12 
13 #if !defined(KRATOS_ELEMENT_H_INCLUDED )
14 #define KRATOS_ELEMENT_H_INCLUDED
15 
16 // System includes
17 
18 // External includes
19 
20 // Project includes
21 #include "includes/properties.h"
22 #include "includes/process_info.h"
27 
28 namespace Kratos
29 {
32 
36 
40 
44 
48 
50 
59 class Element : public GeometricalObject
60 {
61 public:
66 
69 
72 
74  typedef Node NodeType;
75 
81 
84 
87 
88  typedef Vector VectorType;
89 
90  typedef Matrix MatrixType;
91 
92  typedef std::size_t IndexType;
93 
94  typedef std::size_t SizeType;
95 
97 
98  typedef std::vector<std::size_t> EquationIdVectorType;
99 
100  typedef std::vector<DofType::Pointer> DofsVectorType;
101 
103 
106 
109 
112 
121  explicit Element(IndexType NewId = 0)
122  : BaseType(NewId)
123  , mpProperties(nullptr)
124  {
125  }
126 
130  Element(IndexType NewId, const NodesArrayType& ThisNodes)
131  : BaseType(NewId,GeometryType::Pointer(new GeometryType(ThisNodes)))
132  , mpProperties(nullptr)
133  {
134  }
135 
139  Element(IndexType NewId, GeometryType::Pointer pGeometry)
140  : BaseType(NewId,pGeometry)
141  , mpProperties(nullptr)
142  {
143  }
144 
148  Element(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
149  : BaseType(NewId,pGeometry)
150  , mpProperties(pProperties)
151  {
152  }
153 
155 
156  Element(Element const& rOther)
157  : BaseType(rOther)
158  , mpProperties(rOther.mpProperties)
159  {
160  }
161 
163 
164  ~Element() override
165  {
166  }
167 
171 
178 
179  Element & operator=(Element const& rOther)
180  {
182  mpProperties = rOther.mpProperties;
183  return *this;
184  }
185 
189 
202  virtual Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes,
203  PropertiesType::Pointer pProperties) const
204  {
205  KRATOS_TRY
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);
208  KRATOS_CATCH("");
209  }
210 
218  virtual Pointer Create(IndexType NewId,
219  GeometryType::Pointer pGeom,
220  PropertiesType::Pointer pProperties) const
221  {
222  KRATOS_TRY
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);
225  KRATOS_CATCH("");
226  }
227 
235  virtual Pointer Clone (IndexType NewId, NodesArrayType const& ThisNodes) const
236  {
237  KRATOS_TRY
238  KRATOS_WARNING("Element") << " Call base class element Clone " << std::endl;
239  Element::Pointer p_new_elem = Kratos::make_intrusive<Element>(NewId, GetGeometry().Create(ThisNodes), pGetProperties());
240  p_new_elem->SetData(this->GetData());
241  p_new_elem->Set(Flags(*this));
242  return p_new_elem;
243  KRATOS_CATCH("");
244  }
245 
246 
258  virtual void EquationIdVector(EquationIdVectorType& rResult,
259  const ProcessInfo& rCurrentProcessInfo) const
260  {
261  if (rResult.size() != 0) {
262  rResult.resize(0);
263  }
264  }
265 
271  virtual void GetDofList(DofsVectorType& rElementalDofList,
272  const ProcessInfo& rCurrentProcessInfo) const
273  {
274  if (rElementalDofList.size() != 0) {
275  rElementalDofList.resize(0);
276  }
277  }
278 
286  {
287  return pGetGeometry()->GetDefaultIntegrationMethod();
288  }
289 
300  virtual void GetValuesVector(Vector& values, int Step = 0) const
301  {
302  if (values.size() != 0) {
303  values.resize(0, false);
304  }
305  }
306 
310  virtual void GetFirstDerivativesVector(Vector& values, int Step = 0) const
311  {
312  if (values.size() != 0) {
313  values.resize(0, false);
314  }
315  }
316 
320  virtual void GetSecondDerivativesVector(Vector& values, int Step = 0) const
321  {
322  if (values.size() != 0) {
323  values.resize(0, false);
324  }
325  }
326 
341  virtual void Initialize(const ProcessInfo& rCurrentProcessInfo)
342  {
343  }
344 
349  virtual void ResetConstitutiveLaw()
350  {
351  }
352 
365  virtual void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo)
366  {
367  }
368 
372  virtual void InitializeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo)
373  {
374  }
375 
379  virtual void FinalizeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo)
380  {
381  }
382 
386  virtual void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo)
387  {
388  }
389 
405  virtual void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,
406  VectorType& rRightHandSideVector,
407  const ProcessInfo& rCurrentProcessInfo)
408  {
409  if (rLeftHandSideMatrix.size1() != 0) {
410  rLeftHandSideMatrix.resize(0, 0, false);
411  }
412  if (rRightHandSideVector.size() != 0) {
413  rRightHandSideVector.resize(0, false);
414  }
415  }
416 
423  virtual void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix,
424  const ProcessInfo& rCurrentProcessInfo)
425  {
426  if (rLeftHandSideMatrix.size1() != 0) {
427  rLeftHandSideMatrix.resize(0, 0, false);
428  }
429  }
430 
437  virtual void CalculateRightHandSide(VectorType& rRightHandSideVector,
438  const ProcessInfo& rCurrentProcessInfo)
439  {
440  if (rRightHandSideVector.size() != 0) {
441  rRightHandSideVector.resize(0, false);
442  }
443  }
444 
461  virtual void CalculateFirstDerivativesContributions(MatrixType& rLeftHandSideMatrix,
462  VectorType& rRightHandSideVector,
463  const ProcessInfo& rCurrentProcessInfo)
464  {
465  if (rLeftHandSideMatrix.size1() != 0) {
466  rLeftHandSideMatrix.resize(0, 0, false);
467  }
468  if (rRightHandSideVector.size() != 0) {
469  rRightHandSideVector.resize(0, false);
470  }
471  }
472 
479  virtual void CalculateFirstDerivativesLHS(MatrixType& rLeftHandSideMatrix,
480  const ProcessInfo& rCurrentProcessInfo)
481  {
482  if (rLeftHandSideMatrix.size1() != 0) {
483  rLeftHandSideMatrix.resize(0, 0, false);
484  }
485  }
486 
493  virtual void CalculateFirstDerivativesRHS(VectorType& rRightHandSideVector,
494  const ProcessInfo& rCurrentProcessInfo)
495  {
496  if (rRightHandSideVector.size() != 0) {
497  rRightHandSideVector.resize(0, false);
498  }
499  }
500 
518  virtual void CalculateSecondDerivativesContributions(MatrixType& rLeftHandSideMatrix,
519  VectorType& rRightHandSideVector,
520  const ProcessInfo& rCurrentProcessInfo)
521  {
522  if (rLeftHandSideMatrix.size1() != 0) {
523  rLeftHandSideMatrix.resize(0, 0, false);
524  }
525  if (rRightHandSideVector.size() != 0) {
526  rRightHandSideVector.resize(0, false);
527  }
528  }
529 
536  virtual void CalculateSecondDerivativesLHS(MatrixType& rLeftHandSideMatrix,
537  const ProcessInfo& rCurrentProcessInfo)
538  {
539  if (rLeftHandSideMatrix.size1() != 0) {
540  rLeftHandSideMatrix.resize(0, 0, false);
541  }
542  }
543 
550  virtual void CalculateSecondDerivativesRHS(VectorType& rRightHandSideVector,
551  const ProcessInfo& rCurrentProcessInfo)
552  {
553  if (rRightHandSideVector.size() != 0) {
554  rRightHandSideVector.resize(0, false);
555  }
556  }
557 
570  virtual void CalculateMassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo)
571  {
572  if (rMassMatrix.size1() != 0) {
573  rMassMatrix.resize(0, 0, false);
574  }
575  }
576 
583  virtual void CalculateDampingMatrix(MatrixType& rDampingMatrix, const ProcessInfo& rCurrentProcessInfo)
584  {
585  if (rDampingMatrix.size1() != 0) {
586  rDampingMatrix.resize(0, 0, false);
587  }
588  }
589 
597  VectorType& rLumpedMassVector,
598  const ProcessInfo& rCurrentProcessInfo) const
599  {
600  KRATOS_TRY;
601  KRATOS_ERROR << "Calling the CalculateLumpedMassVector() method of the base element. The method must be implemented in the derived element.";
602  KRATOS_CATCH("");
603  }
604 
605 
622  virtual void AddExplicitContribution(const ProcessInfo& rCurrentProcessInfo)
623  {
624  }
625 
635  const VectorType& rRHSVector,
636  const Variable<VectorType>& rRHSVariable,
637  const Variable<double >& rDestinationVariable,
638  const ProcessInfo& rCurrentProcessInfo
639  )
640  {
641  KRATOS_ERROR << "Base element class is not able to assemble rRHS to the desired variable. destination variable is " << rDestinationVariable << std::endl;
642  }
643 
653  const VectorType& rRHSVector,
654  const Variable<VectorType>& rRHSVariable,
655  const Variable<array_1d<double,3> >& rDestinationVariable,
656  const ProcessInfo& rCurrentProcessInfo
657  )
658  {
659  KRATOS_ERROR << "Base element class is not able to assemble rRHS to the desired variable. destination variable is " << rDestinationVariable << std::endl;
660  }
661 
671  const MatrixType& rLHSMatrix,
672  const Variable<MatrixType>& rLHSVariable,
673  const Variable<Matrix>& rDestinationVariable,
674  const ProcessInfo& rCurrentProcessInfo
675  )
676  {
677  KRATOS_ERROR << "Base element class is not able to assemble rLHS to the desired variable. destination variable is " << rDestinationVariable << std::endl;
678  }
679 
686  virtual void Calculate(const Variable<double>& rVariable,
687  double& Output,
688  const ProcessInfo& rCurrentProcessInfo)
689  {
690  }
691 
692  virtual void Calculate(const Variable<array_1d<double, 3 > >& rVariable,
693  array_1d<double, 3 > & Output,
694  const ProcessInfo& rCurrentProcessInfo)
695  {
696  }
697 
698  virtual void Calculate(const Variable<Vector >& rVariable,
699  Vector& Output,
700  const ProcessInfo& rCurrentProcessInfo)
701  {
702  }
703 
704  virtual void Calculate(const Variable<Matrix >& rVariable,
705  Matrix& Output,
706  const ProcessInfo& rCurrentProcessInfo)
707  {
708  }
709 
720  const Variable<bool>& rVariable,
721  std::vector<bool>& rOutput,
722  const ProcessInfo& rCurrentProcessInfo)
723  {
724  }
725 
727  const Variable<int>& rVariable,
728  std::vector<int>& rOutput,
729  const ProcessInfo& rCurrentProcessInfo)
730  {
731  }
732 
734  const Variable<double>& rVariable,
735  std::vector<double>& rOutput,
736  const ProcessInfo& rCurrentProcessInfo)
737  {
738  }
739 
741  const Variable<array_1d<double, 3>>& rVariable,
742  std::vector< array_1d<double, 3>>& rOutput,
743  const ProcessInfo& rCurrentProcessInfo)
744  {
745  }
746 
748  const Variable<array_1d<double, 4>>& rVariable,
749  std::vector< array_1d<double, 4>>& rOutput,
750  const ProcessInfo& rCurrentProcessInfo)
751  {
752  }
753 
755  const Variable<array_1d<double, 6>>& rVariable,
756  std::vector<array_1d<double, 6>>& rOutput,
757  const ProcessInfo& rCurrentProcessInfo)
758  {
759  }
760 
762  const Variable<array_1d<double, 9>>& rVariable,
763  std::vector<array_1d<double, 9>>& rOutput,
764  const ProcessInfo& rCurrentProcessInfo)
765  {
766  }
767 
769  const Variable<Vector>& rVariable,
770  std::vector<Vector>& rOutput,
771  const ProcessInfo& rCurrentProcessInfo)
772  {
773  }
774 
776  const Variable<Matrix>& rVariable,
777  std::vector<Matrix>& rOutput,
778  const ProcessInfo& rCurrentProcessInfo)
779  {
780  }
781 
783  const Variable<ConstitutiveLaw::Pointer>& rVariable,
784  std::vector<ConstitutiveLaw::Pointer>& rOutput,
785  const ProcessInfo& rCurrentProcessInfo)
786  {
787  }
788 
789  // virtual void CalculateOnIntegrationPoints(const Variable<ConstitutiveLaw::StrainVectorType>& rVariable,
790  // std::vector<ConstitutiveLaw::StrainVectorType>& rOutput,
791  // const ProcessInfo& rCurrentProcessInfo)
792  // {
793  // }
794 
795  // virtual void CalculateOnIntegrationPoints(const Variable<ConstitutiveLaw::StressVectorType>& rVariable,
796  // std::vector<ConstitutiveLaw::StressVectorType>& rOutput,
797  // const ProcessInfo& rCurrentProcessInfo)
798  // {
799  // }
800 
801  // virtual void CalculateOnIntegrationPoints(const Variable<ConstitutiveLaw::VoigtSizeMatrixType>& rVariable,
802  // std::vector<ConstitutiveLaw::VoigtSizeMatrixType>& rOutput,
803  // const ProcessInfo& rCurrentProcessInfo)
804  // {
805  // }
806 
807  // virtual void CalculateOnIntegrationPoints(const Variable<ConstitutiveLaw::DeformationGradientMatrixType>& rVariable,
808  // std::vector<ConstitutiveLaw::DeformationGradientMatrixType>& rOutput,
809  // const ProcessInfo& rCurrentProcessInfo)
810  // {
811  // }
812 
824  //SET ON INTEGRATION POINTS - METHODS
826  const Variable<bool>& rVariable,
827  const std::vector<bool>& rValues,
828  const ProcessInfo& rCurrentProcessInfo)
829  {
830  }
832  const Variable<int>& rVariable,
833  const std::vector<int>& rValues,
834  const ProcessInfo& rCurrentProcessInfo)
835  {
836  }
837 
839  const Variable<double>& rVariable,
840  const std::vector<double>& rValues,
841  const ProcessInfo& rCurrentProcessInfo)
842  {
843  }
844 
846  const Variable<array_1d<double, 3>>& rVariable,
847  const std::vector<array_1d<double, 3>>& rValues,
848  const ProcessInfo& rCurrentProcessInfo)
849  {
850  }
851 
853  const Variable<array_1d<double, 4>>& rVariable,
854  const std::vector<array_1d<double, 4>>& rValues,
855  const ProcessInfo& rCurrentProcessInfo)
856  {
857  }
858 
860  const Variable<array_1d<double, 6>>& rVariable,
861  const std::vector<array_1d<double, 6>>& rValues,
862  const ProcessInfo& rCurrentProcessInfo)
863  {
864  }
865 
867  const Variable<array_1d<double, 9>>& rVariable,
868  const std::vector<array_1d<double, 9>>& rValues,
869  const ProcessInfo& rCurrentProcessInfo)
870  {
871  }
872 
874  const Variable<Vector>& rVariable,
875  const std::vector<Vector>& rValues,
876  const ProcessInfo& rCurrentProcessInfo)
877  {
878  }
879 
881  const Variable<Matrix>& rVariable,
882  const std::vector<Matrix>& rValues,
883  const ProcessInfo& rCurrentProcessInfo)
884  {
885  }
886 
888  const Variable<ConstitutiveLaw::Pointer>& rVariable,
889  const std::vector<ConstitutiveLaw::Pointer>& rValues,
890  const ProcessInfo& rCurrentProcessInfo)
891  {
892  }
893 
904  virtual int Check(const ProcessInfo& rCurrentProcessInfo) const
905  {
906  KRATOS_TRY
907 
908  KRATOS_ERROR_IF( this->Id() < 1 ) << "Element found with Id " << this->Id() << std::endl;
909 
910  const double domain_size = this->GetGeometry().DomainSize();
911  KRATOS_ERROR_IF( domain_size <= 0.0 ) << "Element " << this->Id() << " has non-positive size " << domain_size << std::endl;
912 
913  GetGeometry().Check();
914 
915  return 0;
916 
917  KRATOS_CATCH("")
918  }
919 
926  virtual void MassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo)
927  {
928  if (rMassMatrix.size1() != 0) {
929  rMassMatrix.resize(0, 0, false);
930  }
931  }
932 
939  virtual void AddMassMatrix(MatrixType& rLeftHandSideMatrix,
940  double coeff, const ProcessInfo& rCurrentProcessInfo)
941  {
942  }
943 
950  virtual void DampMatrix(MatrixType& rDampMatrix, const ProcessInfo& rCurrentProcessInfo)
951  {
952  if (rDampMatrix.size1() != 0) {
953  rDampMatrix.resize(0, 0, false);
954  }
955  }
956 
961  virtual void AddInertiaForces(VectorType& rRightHandSideVector, double coeff,
962  const ProcessInfo& rCurrentProcessInfo)
963  {
964  }
965 
972  virtual void CalculateLocalVelocityContribution(MatrixType& rDampingMatrix,
973  VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo)
974  {
975  if (rDampingMatrix.size1() != 0) {
976  rDampingMatrix.resize(0, 0, false);
977  }
978  }
979 
983  virtual void CalculateSensitivityMatrix(const Variable<double>& rDesignVariable,
984  Matrix& rOutput,
985  const ProcessInfo& rCurrentProcessInfo)
986  {
987  if (rOutput.size1() != 0)
988  rOutput.resize(0, 0, false);
989  }
990 
994  virtual void CalculateSensitivityMatrix(const Variable<array_1d<double,3> >& rDesignVariable,
995  Matrix& rOutput,
996  const ProcessInfo& rCurrentProcessInfo)
997  {
998  if (rOutput.size1() != 0)
999  rOutput.resize(0, 0, false);
1000  }
1001 
1002  //METHODS TO BE CLEANED: DEPRECATED end
1003 
1007 
1014  PropertiesType::Pointer pGetProperties()
1015  {
1016  return mpProperties;
1017  }
1018 
1019  const PropertiesType::Pointer pGetProperties() const
1020  {
1021  return mpProperties;
1022  }
1023 
1025  {
1026  KRATOS_DEBUG_ERROR_IF(mpProperties == nullptr)
1027  << "Tryining to get the properties of " << Info()
1028  << ", which are uninitialized." << std::endl;
1029  return *mpProperties;
1030  }
1031 
1033  {
1034  KRATOS_DEBUG_ERROR_IF(mpProperties == nullptr)
1035  << "Tryining to get the properties of " << Info()
1036  << ", which are uninitialized." << std::endl;
1037  return *mpProperties;
1038  }
1039 
1040  void SetProperties(PropertiesType::Pointer pProperties)
1041  {
1042  mpProperties = pProperties;
1043  }
1044 
1048 
1050  bool HasProperties() const
1051  {
1052  return mpProperties != nullptr;
1053  }
1054 
1058 
1087  virtual const Parameters GetSpecifications() const
1088  {
1089  const Parameters specifications = Parameters(R"({
1090  "time_integration" : [],
1091  "framework" : "lagrangian",
1092  "symmetric_lhs" : false,
1093  "positive_definite_lhs" : false,
1094  "output" : {
1095  "gauss_point" : [],
1096  "nodal_historical" : [],
1097  "nodal_non_historical" : [],
1098  "entity" : []
1099  },
1100  "required_variables" : [],
1101  "required_dofs" : [],
1102  "flags_used" : [],
1103  "compatible_geometries" : [],
1104  "element_integrates_in_time" : true,
1105  "compatible_constitutive_laws": {
1106  "type" : [],
1107  "dimension" : [],
1108  "strain_size" : []
1109  },
1110  "required_polynomial_degree_of_geometry" : -1,
1111  "documentation" : "This is the base element"
1112 
1113  })");
1114  return specifications;
1115  }
1116 
1118 
1119  std::string Info() const override
1120  {
1121  std::stringstream buffer;
1122  buffer << "Element #" << Id();
1123  return buffer.str();
1124  }
1125 
1127 
1128  void PrintInfo(std::ostream& rOStream) const override
1129  {
1130  rOStream << "Element #" << Id();
1131  }
1132 
1134 
1135  void PrintData(std::ostream& rOStream) const override
1136  {
1137  pGetGeometry()->PrintData(rOStream);
1138  }
1139 
1144 
1145 protected:
1167 
1168 private:
1174 
1178  Properties::Pointer mpProperties;
1179 
1189 
1190  friend class Serializer;
1191 
1192  void save(Serializer& rSerializer) const override
1193  {
1195  rSerializer.save("Properties", mpProperties);
1196  }
1197 
1198  void load(Serializer& rSerializer) override
1199  {
1201  rSerializer.load("Properties", mpProperties);
1202  }
1203 
1214 
1215 }; // Class Element
1216 
1223 
1225 inline std::istream & operator >>(std::istream& rIStream,
1226  Element& rThis);
1227 
1229 
1230 inline std::ostream & operator <<(std::ostream& rOStream,
1231  const Element& rThis)
1232 {
1233  rThis.PrintInfo(rOStream);
1234  rOStream << " : " << std::endl;
1235  rThis.PrintData(rOStream);
1236  return rOStream;
1237 }
1239 
1241 
1242 void KRATOS_API(KRATOS_CORE) AddKratosComponent(std::string const& Name, Element const& ThisComponent);
1243 
1248 #undef KRATOS_EXPORT_MACRO
1249 #define KRATOS_EXPORT_MACRO KRATOS_API
1250 
1252 
1253 #undef KRATOS_EXPORT_MACRO
1254 #define KRATOS_EXPORT_MACRO KRATOS_NO_EXPORT
1255 
1256 } // namespace Kratos.
1257 #endif // KRATOS_ELEMENT_H_INCLUDED defined
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