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.
condition.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_CONDITION_H_INCLUDED )
14 #define KRATOS_CONDITION_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"
26 
27 namespace Kratos
28 {
31 
35 
39 
43 
47 
49 
59 {
60 public:
63 
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 
110 
114 
123  explicit Condition(IndexType NewId = 0)
124  : BaseType(NewId)
125  , mpProperties(nullptr)
126  {
127  }
128 
132  Condition(IndexType NewId, const NodesArrayType& ThisNodes)
133  : BaseType(NewId,GeometryType::Pointer(new GeometryType(ThisNodes)))
134  , mpProperties(nullptr)
135  {
136  }
137 
141  Condition(IndexType NewId, GeometryType::Pointer pGeometry)
142  : BaseType(NewId,pGeometry)
143  , mpProperties(nullptr)
144  {
145  }
146 
150  Condition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
151  : BaseType(NewId,pGeometry)
152  , mpProperties(pProperties)
153  {
154  }
155 
157 
158  Condition(Condition const& rOther)
159  : BaseType(rOther)
160  , mpProperties(rOther.mpProperties)
161  {
162  }
163 
165 
166  ~Condition() override
167  {
168  }
169 
173 
180 
181  Condition & operator=(Condition const& rOther)
182  {
183  BaseType::operator=(rOther);
184  mpProperties = rOther.mpProperties;
185 
186  return *this;
187  }
188 
192 
205  virtual Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes,
206  PropertiesType::Pointer pProperties) const
207  {
208  KRATOS_TRY
209  KRATOS_ERROR << "Please implement the First Create method in your derived Condition" << Info() << std::endl;
210  return Kratos::make_intrusive<Condition>(NewId, GetGeometry().Create(ThisNodes), pProperties);
211  KRATOS_CATCH("");
212  }
213 
221  virtual Pointer Create(IndexType NewId,
222  GeometryType::Pointer pGeom,
223  PropertiesType::Pointer pProperties) const
224  {
225  KRATOS_TRY
226  KRATOS_ERROR << "Please implement the Second Create method in your derived Condition" << Info() << std::endl;
227  return Kratos::make_intrusive<Condition>(NewId, pGeom, pProperties);
228  KRATOS_CATCH("");
229  }
230 
238  virtual Pointer Clone (IndexType NewId, NodesArrayType const& ThisNodes) const
239  {
240  KRATOS_TRY
241  KRATOS_WARNING("Condition") << " Call base class condition Clone " << std::endl;
242  Condition::Pointer p_new_cond = Kratos::make_intrusive<Condition>(NewId, GetGeometry().Create(ThisNodes), pGetProperties());
243  p_new_cond->SetData(this->GetData());
244  p_new_cond->Set(Flags(*this));
245  return p_new_cond;
246  KRATOS_CATCH("");
247  }
248 
260  virtual void EquationIdVector(EquationIdVectorType& rResult,
261  const ProcessInfo& rCurrentProcessInfo) const
262  {
263  if (rResult.size() != 0) {
264  rResult.resize(0);
265  }
266  }
267 
273  virtual void GetDofList(DofsVectorType& rElementalDofList,
274  const ProcessInfo& rCurrentProcessInfo) const
275  {
276  if (rElementalDofList.size() != 0) {
277  rElementalDofList.resize(0);
278  }
279  }
280 
289  {
290  return pGetGeometry()->GetDefaultIntegrationMethod();
291  }
292 
303  virtual void GetValuesVector(Vector& values, int Step = 0) const
304  {
305  if (values.size() != 0) {
306  values.resize(0, false);
307  }
308  }
309 
313  virtual void GetFirstDerivativesVector(Vector& values, int Step = 0) const
314  {
315  if (values.size() != 0) {
316  values.resize(0, false);
317  }
318  }
319 
323  virtual void GetSecondDerivativesVector(Vector& values, int Step = 0) const
324  {
325  if (values.size() != 0) {
326  values.resize(0, false);
327  }
328  }
329 
344  virtual void Initialize(const ProcessInfo& rCurrentProcessInfo)
345  {
346  }
347 
352  virtual void ResetConstitutiveLaw()
353  {
354  }
355 
368  virtual void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo)
369  {
370  }
371 
375  virtual void InitializeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo)
376  {
377  }
378 
382  virtual void FinalizeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo)
383  {
384  }
385 
389  virtual void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo)
390  {
391  }
392 
408  virtual void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,
409  VectorType& rRightHandSideVector,
410  const ProcessInfo& rCurrentProcessInfo)
411  {
412  if (rLeftHandSideMatrix.size1() != 0) {
413  rLeftHandSideMatrix.resize(0, 0, false);
414  }
415  if (rRightHandSideVector.size() != 0) {
416  rRightHandSideVector.resize(0, false);
417  }
418  }
419 
426  virtual void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix,
427  const ProcessInfo& rCurrentProcessInfo)
428  {
429  if (rLeftHandSideMatrix.size1() != 0) {
430  rLeftHandSideMatrix.resize(0, 0, false);
431  }
432  }
433 
440  virtual void CalculateRightHandSide(VectorType& rRightHandSideVector,
441  const ProcessInfo& rCurrentProcessInfo)
442  {
443  if (rRightHandSideVector.size() != 0) {
444  rRightHandSideVector.resize(0, false);
445  }
446  }
447 
465  virtual void CalculateFirstDerivativesContributions(MatrixType& rLeftHandSideMatrix,
466  VectorType& rRightHandSideVector,
467  const ProcessInfo& rCurrentProcessInfo)
468  {
469  if (rLeftHandSideMatrix.size1() != 0) {
470  rLeftHandSideMatrix.resize(0, 0, false);
471  }
472  if (rRightHandSideVector.size() != 0) {
473  rRightHandSideVector.resize(0, false);
474  }
475  }
476 
483  virtual void CalculateFirstDerivativesLHS(MatrixType& rLeftHandSideMatrix,
484  const ProcessInfo& rCurrentProcessInfo)
485  {
486  if (rLeftHandSideMatrix.size1() != 0) {
487  rLeftHandSideMatrix.resize(0, 0, false);
488  }
489  }
490 
497  virtual void CalculateFirstDerivativesRHS(VectorType& rRightHandSideVector,
498  const ProcessInfo& rCurrentProcessInfo)
499  {
500  if (rRightHandSideVector.size() != 0) {
501  rRightHandSideVector.resize(0, false);
502  }
503  }
504 
522  virtual void CalculateSecondDerivativesContributions(MatrixType& rLeftHandSideMatrix,
523  VectorType& rRightHandSideVector,
524  const ProcessInfo& rCurrentProcessInfo)
525  {
526  if (rLeftHandSideMatrix.size1() != 0) {
527  rLeftHandSideMatrix.resize(0, 0, false);
528  }
529  if (rRightHandSideVector.size() != 0) {
530  rRightHandSideVector.resize(0, false);
531  }
532  }
533 
540  virtual void CalculateSecondDerivativesLHS(MatrixType& rLeftHandSideMatrix,
541  const ProcessInfo& rCurrentProcessInfo)
542  {
543  if (rLeftHandSideMatrix.size1() != 0) {
544  rLeftHandSideMatrix.resize(0, 0, false);
545  }
546  }
547 
554  virtual void CalculateSecondDerivativesRHS(VectorType& rRightHandSideVector,
555  const ProcessInfo& rCurrentProcessInfo)
556  {
557  if (rRightHandSideVector.size() != 0) {
558  rRightHandSideVector.resize(0, false);
559  }
560  }
561 
574  virtual void CalculateMassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo)
575  {
576  if (rMassMatrix.size1() != 0) {
577  rMassMatrix.resize(0, 0, false);
578  }
579  }
586  virtual void CalculateDampingMatrix(MatrixType& rDampingMatrix, const ProcessInfo& rCurrentProcessInfo)
587  {
588  if (rDampingMatrix.size1() != 0) {
589  rDampingMatrix.resize(0, 0, false);
590  }
591  }
592 
609  virtual void AddExplicitContribution(const ProcessInfo& rCurrentProcessInfo)
610  {
611  }
612 
622  const VectorType& rRHSVector,
623  const Variable<VectorType>& rRHSVariable,
624  const Variable<double >& rDestinationVariable,
625  const ProcessInfo& rCurrentProcessInfo
626  )
627  {
628  KRATOS_ERROR << "Base condition class is not able to assemble rRHS to the desired variable. destination variable is " << rDestinationVariable << std::endl;
629  }
630 
640  const VectorType& rRHSVector,
641  const Variable<VectorType>& rRHSVariable,
642  const Variable<array_1d<double,3> >& rDestinationVariable,
643  const ProcessInfo& rCurrentProcessInfo
644  )
645  {
646  KRATOS_ERROR << "Base condition class is not able to assemble rRHS to the desired variable. destination variable is " << rDestinationVariable << std::endl;
647  }
648 
658  const MatrixType& rLHSMatrix,
659  const Variable<MatrixType>& rLHSVariable,
660  const Variable<Matrix>& rDestinationVariable,
661  const ProcessInfo& rCurrentProcessInfo
662  )
663  {
664  KRATOS_ERROR << "Base condition class is not able to assemble rLHS to the desired variable. destination variable is " << rDestinationVariable << std::endl;
665  }
666 
673  virtual void Calculate(const Variable<double >& rVariable,
674  double& Output,
675  const ProcessInfo& rCurrentProcessInfo)
676  {
677  }
678 
679  virtual void Calculate(const Variable< array_1d<double,3> >& rVariable,
680  array_1d<double,3>& Output,
681  const ProcessInfo& rCurrentProcessInfo)
682  {
683  }
684 
685  virtual void Calculate(const Variable<Vector >& rVariable,
686  Vector& Output,
687  const ProcessInfo& rCurrentProcessInfo)
688  {
689  }
690 
691  virtual void Calculate(const Variable<Matrix >& rVariable,
692  Matrix& Output,
693  const ProcessInfo& rCurrentProcessInfo)
694  {
695  }
696 
707  const Variable<bool>& rVariable,
708  std::vector<bool>& rOutput,
709  const ProcessInfo& rCurrentProcessInfo)
710  {
711  }
712 
714  const Variable<int>& rVariable,
715  std::vector<int>& rOutput,
716  const ProcessInfo& rCurrentProcessInfo)
717  {
718  }
719 
721  const Variable<double>& rVariable,
722  std::vector<double>& rOutput,
723  const ProcessInfo& rCurrentProcessInfo)
724  {
725  }
726 
728  const Variable<array_1d<double, 3>>& rVariable,
729  std::vector< array_1d<double, 3>>& rOutput,
730  const ProcessInfo& rCurrentProcessInfo)
731  {
732  }
733 
735  const Variable<array_1d<double, 4>>& rVariable,
736  std::vector< array_1d<double, 4>>& rOutput,
737  const ProcessInfo& rCurrentProcessInfo)
738  {
739  }
740 
742  const Variable<array_1d<double, 6>>& rVariable,
743  std::vector<array_1d<double, 6>>& rOutput,
744  const ProcessInfo& rCurrentProcessInfo)
745  {
746  }
747 
749  const Variable<array_1d<double, 9>>& rVariable,
750  std::vector<array_1d<double, 9>>& rOutput,
751  const ProcessInfo& rCurrentProcessInfo)
752  {
753  }
754 
756  const Variable<Vector>& rVariable,
757  std::vector<Vector>& rOutput,
758  const ProcessInfo& rCurrentProcessInfo)
759  {
760  }
761 
763  const Variable<Matrix>& rVariable,
764  std::vector<Matrix>& rOutput,
765  const ProcessInfo& rCurrentProcessInfo)
766  {
767  }
768 
780  //SET ON INTEGRATION POINTS - METHODS
781 
783  const Variable<bool>& rVariable,
784  const std::vector<bool>& rValues,
785  const ProcessInfo& rCurrentProcessInfo)
786  {
787  }
788 
790  const Variable<int>& rVariable,
791  const std::vector<int>& rValues,
792  const ProcessInfo& rCurrentProcessInfo)
793  {
794  }
795 
797  const Variable<double>& rVariable,
798  const std::vector<double>& rValues,
799  const ProcessInfo& rCurrentProcessInfo)
800  {
801  }
802 
804  const Variable<array_1d<double, 3>>& rVariable,
805  const std::vector<array_1d<double, 3>>& rValues,
806  const ProcessInfo& rCurrentProcessInfo)
807  {
808  }
809 
811  const Variable<array_1d<double, 4>>& rVariable,
812  const std::vector<array_1d<double, 4>>& rValues,
813  const ProcessInfo& rCurrentProcessInfo)
814  {
815  }
816 
818  const Variable<array_1d<double, 6>>& rVariable,
819  const std::vector<array_1d<double, 6>>& rValues,
820  const ProcessInfo& rCurrentProcessInfo)
821  {
822  }
823 
825  const Variable<array_1d<double, 9>>& rVariable,
826  const std::vector<array_1d<double, 9>>& rValues,
827  const ProcessInfo& rCurrentProcessInfo)
828  {
829  }
830 
832  const Variable<Vector>& rVariable,
833  const std::vector<Vector>& rValues,
834  const ProcessInfo& rCurrentProcessInfo)
835  {
836  }
837 
839  const Variable<Matrix>& rVariable,
840  const std::vector<Matrix>& rValues,
841  const ProcessInfo& rCurrentProcessInfo)
842  {
843  }
844 
854  virtual int Check(const ProcessInfo& rCurrentProcessInfo) const
855  {
856  KRATOS_TRY
857 
858  KRATOS_ERROR_IF( this->Id() < 1 ) << "Condition found with Id " << this->Id() << std::endl;
859 
860  const double domain_size = this->GetGeometry().DomainSize();
861  KRATOS_ERROR_IF( domain_size < 0.0 ) << "Condition " << this->Id() << " has negative size " << domain_size << std::endl;
862 
863  GetGeometry().Check();
864 
865  return 0;
866 
867  KRATOS_CATCH("")
868  }
869 
876  virtual void MassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo)
877  {
878  if (rMassMatrix.size1() != 0) {
879  rMassMatrix.resize(0, 0, false);
880  }
881  }
882 
889  virtual void AddMassMatrix(MatrixType& rLeftHandSideMatrix,
890  double coeff, const ProcessInfo& rCurrentProcessInfo)
891  {
892  }
893 
900  virtual void DampMatrix(MatrixType& rDampMatrix, const ProcessInfo& rCurrentProcessInfo)
901  {
902  if (rDampMatrix.size1() != 0) {
903  rDampMatrix.resize(0, 0, false);
904  }
905  }
906 
911  virtual void AddInertiaForces(VectorType& rRightHandSideVector, double coeff,
912  const ProcessInfo& rCurrentProcessInfo)
913  {
914  }
915 
922  virtual void CalculateLocalVelocityContribution(MatrixType& rDampingMatrix,
923  VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo)
924  {
925  if (rDampingMatrix.size1() != 0) {
926  rDampingMatrix.resize(0, 0, false);
927  }
928  }
929 
933  virtual void CalculateSensitivityMatrix(const Variable<double>& rDesignVariable,
934  Matrix& rOutput,
935  const ProcessInfo& rCurrentProcessInfo)
936  {
937  if (rOutput.size1() != 0)
938  rOutput.resize(0, 0, false);
939  }
940 
944  virtual void CalculateSensitivityMatrix(const Variable<array_1d<double,3> >& rDesignVariable,
945  Matrix& rOutput,
946  const ProcessInfo& rCurrentProcessInfo)
947  {
948  if (rOutput.size1() != 0)
949  rOutput.resize(0, 0, false);
950  }
951 
952  //METHODS TO BE CLEANED: DEPRECATED end
953 
957 
964  PropertiesType::Pointer pGetProperties()
965  {
966  return mpProperties;
967  }
968 
969  const PropertiesType::Pointer pGetProperties() const
970  {
971  return mpProperties;
972  }
973 
975  {
976  KRATOS_DEBUG_ERROR_IF(mpProperties == nullptr)
977  << "Tryining to get the properties of " << Info()
978  << ", which are uninitialized." << std::endl;
979  return *mpProperties;
980  }
981 
983  {
984  KRATOS_DEBUG_ERROR_IF(mpProperties == nullptr)
985  << "Tryining to get the properties of " << Info()
986  << ", which are uninitialized." << std::endl;
987  return *mpProperties;
988  }
989 
990  void SetProperties(PropertiesType::Pointer pProperties)
991  {
992  mpProperties = pProperties;
993  }
994 
998 
1000  bool HasProperties() const
1001  {
1002  return mpProperties != nullptr;
1003  }
1004 
1008 
1037  virtual const Parameters GetSpecifications() const
1038  {
1039  const Parameters specifications = Parameters(R"({
1040  "time_integration" : [],
1041  "framework" : "lagrangian",
1042  "symmetric_lhs" : false,
1043  "positive_definite_lhs" : false,
1044  "output" : {
1045  "gauss_point" : [],
1046  "nodal_historical" : [],
1047  "nodal_non_historical" : [],
1048  "entity" : []
1049  },
1050  "required_variables" : [],
1051  "required_dofs" : [],
1052  "flags_used" : [],
1053  "compatible_geometries" : [],
1054  "element_integrates_in_time" : true,
1055  "compatible_constitutive_laws": {
1056  "type" : [],
1057  "dimension" : [],
1058  "strain_size" : []
1059  },
1060  "required_polynomial_degree_of_geometry" : -1,
1061  "documentation" : "This is the base condition"
1062 
1063  })");
1064  return specifications;
1065  }
1066 
1068 
1069  std::string Info() const override
1070  {
1071  std::stringstream buffer;
1072  buffer << "Condition #" << Id();
1073  return buffer.str();
1074  }
1075 
1077 
1078  void PrintInfo(std::ostream& rOStream) const override
1079  {
1080  rOStream << "Condition #" << Id();
1081  }
1082 
1084 
1085  void PrintData(std::ostream& rOStream) const override
1086  {
1087  pGetGeometry()->PrintData(rOStream);
1088  }
1089 
1094 
1095 protected:
1117 
1118 private:
1124 
1128  Properties::Pointer mpProperties;
1129 
1139 
1140  friend class Serializer;
1141 
1142  void save(Serializer& rSerializer) const override
1143  {
1145  rSerializer.save("Properties", mpProperties);
1146  }
1147 
1148  void load(Serializer& rSerializer) override
1149  {
1151  rSerializer.load("Properties", mpProperties);
1152  }
1153 
1164 
1165 }; // Class Condition
1166 
1173 
1175 inline std::istream & operator >>(std::istream& rIStream,
1176  Condition& rThis);
1177 
1179 
1180 inline std::ostream & operator <<(std::ostream& rOStream,
1181  const Condition& rThis)
1182 {
1183  rThis.PrintInfo(rOStream);
1184  rOStream << " : " << std::endl;
1185  rThis.PrintData(rOStream);
1186 
1187  return rOStream;
1188 }
1190 
1192 
1193 void KRATOS_API(KRATOS_CORE) AddKratosComponent(std::string const& Name, Condition const& ThisComponent);
1194 
1199 #undef KRATOS_EXPORT_MACRO
1200 #define KRATOS_EXPORT_MACRO KRATOS_API
1201 
1203 
1204 #undef KRATOS_EXPORT_MACRO
1205 #define KRATOS_EXPORT_MACRO KRATOS_NO_EXPORT
1206 
1207 } // namespace Kratos.
1208 #endif // KRATOS_CONDITION_H_INCLUDED defined
1209 
Base class for all Conditions.
Definition: condition.h:59
virtual void GetDofList(DofsVectorType &rElementalDofList, const ProcessInfo &rCurrentProcessInfo) const
Definition: condition.h:273
virtual void CalculateSensitivityMatrix(const Variable< double > &rDesignVariable, Matrix &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:933
virtual void CalculateOnIntegrationPoints(const Variable< bool > &rVariable, std::vector< bool > &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:706
virtual Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const
It creates a new condition pointer.
Definition: condition.h:205
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: condition.h:1085
std::size_t SizeType
Definition: condition.h:94
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: condition.h:105
virtual const Parameters GetSpecifications() const
This method provides the specifications/requirements of the element.
Definition: condition.h:1037
virtual void CalculateOnIntegrationPoints(const Variable< array_1d< double, 3 >> &rVariable, std::vector< array_1d< double, 3 >> &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:727
virtual void Calculate(const Variable< array_1d< double, 3 > > &rVariable, array_1d< double, 3 > &Output, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:679
virtual void SetValuesOnIntegrationPoints(const Variable< array_1d< double, 3 >> &rVariable, const std::vector< array_1d< double, 3 >> &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:803
virtual void CalculateLocalVelocityContribution(MatrixType &rDampingMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:922
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(Condition)
Pointer definition of Condition.
GeometricalObject BaseType
base type: an GeometricalObject that automatically has a unique number
Definition: condition.h:71
virtual void ResetConstitutiveLaw()
Definition: condition.h:352
Condition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Definition: condition.h:150
virtual void DampMatrix(MatrixType &rDampMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:900
PointerVectorSet< DofType > DofsArrayType
Definition: condition.h:102
virtual void CalculateSecondDerivativesLHS(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:540
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 condition to assemble an rRHS vector identified by a variable r...
Definition: condition.h:639
Dof< double > DofType
Definition: condition.h:96
virtual void FinalizeSolutionStep(const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:389
std::size_t IndexType
Definition: condition.h:92
Node NodeType
definition of node type (default is: Node)
Definition: condition.h:74
virtual Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const
It creates a new condition pointer.
Definition: condition.h:221
Geometry< NodeType >::PointsArrayType NodesArrayType
definition of nodes container type, redefined from GeometryType
Definition: condition.h:86
std::vector< std::size_t > EquationIdVectorType
Definition: condition.h:98
virtual void GetSecondDerivativesVector(Vector &values, int Step=0) const
Definition: condition.h:323
virtual void CalculateFirstDerivativesLHS(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:483
virtual void Calculate(const Variable< Vector > &rVariable, Vector &Output, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:685
virtual void CalculateOnIntegrationPoints(const Variable< int > &rVariable, std::vector< int > &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:713
virtual void AddMassMatrix(MatrixType &rLeftHandSideMatrix, double coeff, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:889
Matrix MatrixType
Definition: condition.h:90
virtual void Calculate(const Variable< Matrix > &rVariable, Matrix &Output, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:691
virtual void AddInertiaForces(VectorType &rRightHandSideVector, double coeff, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:911
virtual Pointer Clone(IndexType NewId, NodesArrayType const &ThisNodes) const
It creates a new condition pointer and clones the previous condition data.
Definition: condition.h:238
std::string Info() const override
Turn back information as a string.
Definition: condition.h:1069
Condition(IndexType NewId, GeometryType::Pointer pGeometry)
Definition: condition.h:141
PropertiesType & GetProperties()
Definition: condition.h:974
virtual void CalculateOnIntegrationPoints(const Variable< double > &rVariable, std::vector< double > &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:720
virtual void SetValuesOnIntegrationPoints(const Variable< bool > &rVariable, const std::vector< bool > &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:782
virtual void CalculateFirstDerivativesRHS(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:497
virtual void InitializeSolutionStep(const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:368
bool HasProperties() const
Check that the Condition has a correctly initialized pointer to a Properties instance.
Definition: condition.h:1000
std::vector< DofType::Pointer > DofsVectorType
Definition: condition.h:100
virtual void SetValuesOnIntegrationPoints(const Variable< array_1d< double, 9 >> &rVariable, const std::vector< array_1d< double, 9 >> &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:824
void SetProperties(PropertiesType::Pointer pProperties)
Definition: condition.h:990
Vector VectorType
Definition: condition.h:88
virtual void GetFirstDerivativesVector(Vector &values, int Step=0) const
Definition: condition.h:313
virtual void SetValuesOnIntegrationPoints(const Variable< int > &rVariable, const std::vector< int > &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:789
virtual void SetValuesOnIntegrationPoints(const Variable< Vector > &rVariable, const std::vector< Vector > &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:831
virtual void CalculateOnIntegrationPoints(const Variable< array_1d< double, 6 >> &rVariable, std::vector< array_1d< double, 6 >> &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:741
virtual void SetValuesOnIntegrationPoints(const Variable< array_1d< double, 4 >> &rVariable, const std::vector< array_1d< double, 4 >> &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:810
Properties PropertiesType
Definition: condition.h:80
GeometryData GeometryDataType
Definition: condition.h:107
virtual void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:426
Condition ConditionType
definition of condition type
Definition: condition.h:68
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: condition.h:1078
virtual void CalculateMassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:574
virtual void CalculateFirstDerivativesContributions(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:465
virtual void CalculateSensitivityMatrix(const Variable< array_1d< double, 3 > > &rDesignVariable, Matrix &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:944
virtual void MassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:876
virtual void SetValuesOnIntegrationPoints(const Variable< array_1d< double, 6 >> &rVariable, const std::vector< array_1d< double, 6 >> &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:817
virtual int Check(const ProcessInfo &rCurrentProcessInfo) const
Definition: condition.h:854
const PropertiesType::Pointer pGetProperties() const
Definition: condition.h:969
virtual void AddExplicitContribution(const MatrixType &rLHSMatrix, const Variable< MatrixType > &rLHSVariable, const Variable< Matrix > &rDestinationVariable, const ProcessInfo &rCurrentProcessInfo)
This function is designed to make the condition to assemble an rRHS vector identified by a variable r...
Definition: condition.h:657
Condition & operator=(Condition const &rOther)
Assignment operator.
Definition: condition.h:181
virtual void SetValuesOnIntegrationPoints(const Variable< double > &rVariable, const std::vector< double > &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:796
virtual void Initialize(const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:344
~Condition() override
Destructor.
Definition: condition.h:166
virtual void SetValuesOnIntegrationPoints(const Variable< Matrix > &rVariable, const std::vector< Matrix > &rValues, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:838
Condition(Condition const &rOther)
Copy constructor.
Definition: condition.h:158
PropertiesType const & GetProperties() const
Definition: condition.h:982
virtual void InitializeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:375
virtual void CalculateOnIntegrationPoints(const Variable< array_1d< double, 9 >> &rVariable, std::vector< array_1d< double, 9 >> &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:748
virtual void CalculateDampingMatrix(MatrixType &rDampingMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:586
PropertiesType::Pointer pGetProperties()
returns the pointer to the property of the condition. Does not throw an error, to allow copying of co...
Definition: condition.h:964
virtual IntegrationMethod GetIntegrationMethod() const
Definition: condition.h:288
virtual void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const
Definition: condition.h:260
virtual void CalculateSecondDerivativesRHS(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:554
virtual void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:440
virtual void CalculateOnIntegrationPoints(const Variable< Matrix > &rVariable, std::vector< Matrix > &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:762
Condition(IndexType NewId=0)
Definition: condition.h:123
virtual void Calculate(const Variable< double > &rVariable, double &Output, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:673
virtual void FinalizeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:382
Geometry< NodeType > GeometryType
definition of the geometry type with given NodeType
Definition: condition.h:83
virtual void CalculateOnIntegrationPoints(const Variable< array_1d< double, 4 >> &rVariable, std::vector< array_1d< double, 4 >> &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:734
virtual void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:408
virtual void AddExplicitContribution(const VectorType &rRHSVector, const Variable< VectorType > &rRHSVariable, const Variable< double > &rDestinationVariable, const ProcessInfo &rCurrentProcessInfo)
This function is designed to make the condition to assemble an rRHS vector identified by a variable r...
Definition: condition.h:621
Condition(IndexType NewId, const NodesArrayType &ThisNodes)
Definition: condition.h:132
virtual void GetValuesVector(Vector &values, int Step=0) const
Definition: condition.h:303
virtual void AddExplicitContribution(const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:609
virtual void CalculateOnIntegrationPoints(const Variable< Vector > &rVariable, std::vector< Vector > &rOutput, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:755
virtual void CalculateSecondDerivativesContributions(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:522
Dof represents a degree of freedom (DoF).
Definition: dof.h:86
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
This class is a vector which stores global pointers.
Definition: global_pointers_vector.h:61
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