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.
boundary_condition.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosSolidMechanicsApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: August 2017 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_BOUNDARY_CONDITION_H_INCLUDED)
11 #define KRATOS_BOUNDARY_CONDITION_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
18 #include "includes/checks.h"
19 #include "includes/condition.h"
20 #include "custom_utilities/solid_mechanics_math_utilities.hpp"
22 #include "custom_utilities/element_utilities.hpp"
23 
24 namespace Kratos
25 {
40 
42 
47 class KRATOS_API(SOLID_MECHANICS_APPLICATION) BoundaryCondition
48  : public Condition
49 {
50 public:
51 
54 
57 
60 
61  // Counted pointer of BoundaryCondition
64 
65 protected:
66 
71  KRATOS_DEFINE_LOCAL_FLAG( COMPUTE_RHS_VECTOR );
72  KRATOS_DEFINE_LOCAL_FLAG( COMPUTE_LHS_MATRIX );
73 
79  {
80  private:
81 
82  //variables including all integration points
84  const Matrix* pNcontainer;
85 
86  public:
87 
88  //for axisymmetric use only
89  double CurrentRadius;
91 
92  //general variables
93  double GeometrySize;
94  double Jacobian;
98 
99  //external boundary values
102 
103  //boundary characteristics
107 
108  //variables including all integration points
111 
116  {
117  pDN_De=&rDN_De;
118  };
119 
120  void SetShapeFunctions(const Matrix& rNcontainer)
121  {
122  pNcontainer=&rNcontainer;
123  };
124 
125 
130  {
131  return *pDN_De;
132  };
133 
135  {
136  return *pNcontainer;
137  };
138 
139  void Initialize( const unsigned int& dimension,
140  const unsigned int& local_dimension,
141  const unsigned int& number_of_nodes )
142  {
143  //doubles
144  //radius
145  CurrentRadius = 0;
146  ReferenceRadius = 0;
147  //jacobians
148  GeometrySize = 1;
149  Jacobian = 1;
150  //external boundary values
151  ExternalScalarValue = 0;
152  ExternalVectorValue.resize(dimension,false);
153  noalias(ExternalVectorValue) = ZeroVector(dimension);
154  //vectors
155  N.resize(number_of_nodes,false);
156  Normal.resize(dimension,false);
157  Tangent1.resize(dimension,false);
158  Tangent2.resize(dimension,false);
159  noalias(N) = ZeroVector(number_of_nodes);
160  noalias(Normal) = ZeroVector(dimension);
161  noalias(Tangent1) = ZeroVector(dimension);
162  noalias(Tangent2) = ZeroVector(dimension);
163  //matrices
164  DN_De.resize(number_of_nodes, local_dimension,false);
165  noalias(DN_De) = ZeroMatrix(number_of_nodes, local_dimension);
166  DeltaPosition.resize(number_of_nodes, dimension,false);
167  noalias(DeltaPosition) = ZeroMatrix(number_of_nodes, dimension);
168  //others
169  J.resize(1,false);
170  j.resize(1,false);
171  J[0].resize(dimension,dimension,false);
172  j[0].resize(dimension,dimension,false);
175  //pointers
176  pDN_De = NULL;
177  pNcontainer = NULL;
178  }
179 
180  };
181 
190  {
191  private:
192 
193  //for calculation local system with compacted LHS and RHS
194  MatrixType *mpLeftHandSideMatrix;
195  VectorType *mpRightHandSideVector;
196 
197  public:
198 
199  //calculation flags
201 
205  void SetLeftHandSideMatrix( MatrixType& rLeftHandSideMatrix ) { mpLeftHandSideMatrix = &rLeftHandSideMatrix; };
206 
207  void SetRightHandSideVector( VectorType& rRightHandSideVector ) { mpRightHandSideVector = &rRightHandSideVector; };
208 
212  MatrixType& GetLeftHandSideMatrix() { return *mpLeftHandSideMatrix; };
213 
214  VectorType& GetRightHandSideVector() { return *mpRightHandSideVector; };
215 
216  };
217 
218 
219 public:
220 
221 
224 
227 
229  BoundaryCondition( IndexType NewId, GeometryType::Pointer pGeometry );
230 
231  BoundaryCondition( IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties );
232 
234  BoundaryCondition( BoundaryCondition const& rOther);
235 
237  ~BoundaryCondition() override;
238 
242 
243 
247 
255  Condition::Pointer Create(IndexType NewId,
256  NodesArrayType const& ThisNodes,
257  PropertiesType::Pointer pProperties ) const override;
258 
259 
267  Condition::Pointer Clone(IndexType NewId,
268  NodesArrayType const& ThisNodes) const override;
269 
270 
271  //************* STARTING - ENDING METHODS
272 
273 
277  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
278 
282  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
283 
287  void InitializeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
288 
289 
290  //************* GETTING METHODS
291 
295  void GetDofList(DofsVectorType& rConditionDofList,
296  const ProcessInfo& rCurrentProcessInfo ) const override;
297 
301  void EquationIdVector(EquationIdVectorType& rResult,
302  const ProcessInfo& rCurrentProcessInfo ) const override;
303 
307  void GetValuesVector(Vector& rValues,
308  int Step = 0 ) const override;
309 
313  void GetFirstDerivativesVector(Vector& rValues,
314  int Step = 0 ) const override;
315 
319  void GetSecondDerivativesVector(Vector& rValues,
320  int Step = 0 ) const override;
321 
322 
323  //************* COMPUTING METHODS
324 
333  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,
334  VectorType& rRightHandSideVector,
335  const ProcessInfo& rCurrentProcessInfo ) override;
336 
343  void CalculateRightHandSide(VectorType& rRightHandSideVector,
344  const ProcessInfo& rCurrentProcessInfo ) override;
345 
352  void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix,
353  const ProcessInfo& rCurrentProcessInfo) override;
354 
361  void CalculateMassMatrix(MatrixType& rMassMatrix,
362  const ProcessInfo& rCurrentProcessInfo ) override;
363 
370  void CalculateDampingMatrix(MatrixType& rDampingMatrix,
371  const ProcessInfo& rCurrentProcessInfo ) override;
372 
373 
383  void AddExplicitContribution(const VectorType& rRHS,
384  const Variable<VectorType>& rRHSVariable,
385  const Variable<array_1d<double,3> >& rDestinationVariable,
386  const ProcessInfo& rCurrentProcessInfo) override;
387 
388 
392  void CalculateOnIntegrationPoints(const Variable<double>& rVariable,
393  std::vector<double>& rOutput,
394  const ProcessInfo& rCurrentProcessInfo) override;
395 
396 
397 
398  //************************************************************************************
399  //************************************************************************************
407  int Check( const ProcessInfo& rCurrentProcessInfo ) const override;
408 
422 
423 protected:
429 
434 
435 
442 
446  void InitializeExplicitContributions();
447 
448 
452  virtual bool HasVariableDof(VariableVectorType& rVariable) const;
453 
457  virtual bool HasVariableDof(VariableScalarType& rVariable) const;
458 
459 
463  virtual unsigned int GetDofsSize() const;
464 
468  virtual void InitializeSystemMatrices(MatrixType& rLeftHandSideMatrix,
469  VectorType& rRightHandSideVector,
470  Flags& rCalculationFlags);
471 
475  virtual void InitializeConditionVariables(ConditionVariables& rVariables,
476  const ProcessInfo& rCurrentProcessInfo);
477 
481  virtual void CalculateKinematics(ConditionVariables& rVariables,
482  const double& rPointNumber);
483 
484 
488  virtual void CalculateConditionSystem(LocalSystemComponents& rLocalSystem,
489  const ProcessInfo& rCurrentProcessInfo);
490 
491 
495  virtual void CalculateAndAddLHS(LocalSystemComponents& rLocalSystem,
496  ConditionVariables& rVariables,
497  double& rIntegrationWeight);
498 
502  virtual void CalculateAndAddRHS(LocalSystemComponents& rLocalSystem,
503  ConditionVariables& rVariables,
504  double& rIntegrationWeight);
505 
506 
510  virtual void CalculateAndAddKuug(MatrixType& rLeftHandSideMatrix,
511  ConditionVariables& rVariables,
512  double& rIntegrationWeight);
513 
514 
518  virtual void CalculateAndAddExternalForces(Vector& rRightHandSideVector,
519  ConditionVariables& rVariables,
520  double& rIntegrationWeight);
521 
525  virtual double& CalculateAndAddExternalEnergy(double& rEnergy,
526  ConditionVariables& rVariables,
527  double& rIntegrationWeight,
528  const ProcessInfo& rCurrentProcessInfo);
532  void GetNodalDeltaMovements(Vector& rValues, const int& rNode);
533 
537  Vector& GetNodalCurrentValue(const Variable<array_1d<double,3> >&rVariable, Vector& rValue, const unsigned int& rNode);
538 
542  Vector& GetNodalPreviousValue(const Variable<array_1d<double,3> >&rVariable, Vector& rValue, const unsigned int& rNode);
543 
544 
548 
549 
553 
554 
558 
559 
561 
562 private:
565 
566 
570 
574 
575 
579 
580 
584 
585 
589 
593 
594  friend class Serializer;
595 
596  void save(Serializer& rSerializer) const override;
597 
598  void load(Serializer& rSerializer) override;
599 
600 
601 }; // class BoundaryCondition.
602 
603 } // namespace Kratos.
604 
605 #endif // KRATOS_BOUNDARY_CONDITION_H_INCLUDED defined
General Boundary Condition base type for 3D and 2D geometries.
Definition: boundary_condition.hpp:49
KRATOS_DEFINE_LOCAL_FLAG(COMPUTE_RHS_VECTOR)
KRATOS_DEFINE_LOCAL_FLAG(COMPUTE_LHS_MATRIX)
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(BoundaryCondition)
Variable< double > VariableScalarType
Definition: boundary_condition.hpp:56
IntegrationMethod mThisIntegrationMethod
Definition: boundary_condition.hpp:433
Variable< array_1d< double, 3 > > VariableVectorType
Definition: boundary_condition.hpp:55
GeometryData::SizeType SizeType
Type for size.
Definition: boundary_condition.hpp:59
Base class for all Conditions.
Definition: condition.h:59
Definition: flags.h:58
IntegrationMethod
Definition: geometry_data.h:76
std::size_t SizeType
Definition: geometry_data.h:173
Definition: amatrix_interface.h:41
void resize(std::size_t NewSize1, std::size_t NewSize2, bool preserve=0)
Definition: amatrix_interface.h:224
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
std::size_t IndexType
The definition of the index type.
Definition: key_hash.h:35
Matrix MatrixType
Definition: geometrical_transformation_utilities.h:55
Modeler::Pointer Create(const std::string &ModelerName, Model &rModel, const Parameters ModelParameters)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:30
pybind11::list CalculateOnIntegrationPoints(TObject &dummy, const Variable< TDataType > &rVariable, const ProcessInfo &rProcessInfo)
Definition: add_mesh_to_python.cpp:142
void InitializeSolutionStep(ConstructionUtility &rThisUtil, std::string ThermalSubModelPartName, std::string MechanicalSubModelPartName, std::string HeatFluxSubModelPartName, std::string HydraulicPressureSubModelPartName, bool thermal_conditions, bool mechanical_conditions, int phase)
Definition: add_custom_utilities_to_python.cpp:45
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
KratosZeroMatrix< double > ZeroMatrix
Definition: amatrix_interface.h:559
ModelPart::NodesContainerType NodesArrayType
Definition: gid_gauss_point_container.h:42
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
int dimension
Definition: isotropic_damage_automatic_differentiation.py:123
def load(f)
Definition: ode_solve.py:307
int j
Definition: quadrature.py:648
J
Definition: sensitivityMatrix.py:58
N
Definition: sensitivityMatrix.py:29
Definition: boundary_condition.hpp:79
double ExternalScalarValue
Definition: boundary_condition.hpp:100
Vector Tangent2
Definition: boundary_condition.hpp:106
Vector Normal
Definition: boundary_condition.hpp:104
GeometryType::JacobiansType J
Definition: boundary_condition.hpp:110
const GeometryType::ShapeFunctionsGradientsType & GetShapeFunctionsGradients()
Definition: boundary_condition.hpp:129
double ReferenceRadius
Definition: boundary_condition.hpp:90
void SetShapeFunctions(const Matrix &rNcontainer)
Definition: boundary_condition.hpp:120
Matrix DeltaPosition
Definition: boundary_condition.hpp:97
Vector Tangent1
Definition: boundary_condition.hpp:105
GeometryType::JacobiansType j
Definition: boundary_condition.hpp:109
void SetShapeFunctionsGradients(const GeometryType::ShapeFunctionsGradientsType &rDN_De)
Definition: boundary_condition.hpp:115
const Matrix & GetShapeFunctions()
Definition: boundary_condition.hpp:134
double CurrentRadius
Definition: boundary_condition.hpp:89
Vector N
Definition: boundary_condition.hpp:95
double Jacobian
Definition: boundary_condition.hpp:94
Matrix DN_De
Definition: boundary_condition.hpp:96
double GeometrySize
Definition: boundary_condition.hpp:93
Vector ExternalVectorValue
Definition: boundary_condition.hpp:101
void Initialize(const unsigned int &dimension, const unsigned int &local_dimension, const unsigned int &number_of_nodes)
Definition: boundary_condition.hpp:139
Definition: boundary_condition.hpp:190
void SetRightHandSideVector(VectorType &rRightHandSideVector)
Definition: boundary_condition.hpp:207
MatrixType & GetLeftHandSideMatrix()
Definition: boundary_condition.hpp:212
Flags CalculationFlags
Definition: boundary_condition.hpp:200
void SetLeftHandSideMatrix(MatrixType &rLeftHandSideMatrix)
Definition: boundary_condition.hpp:205
VectorType & GetRightHandSideVector()
Definition: boundary_condition.hpp:214