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.
adjoint_monolithic_wall_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: Suneth Warnakulasuriya
11 //
12 
13 #if !defined(KRATOS_ADJOINT_MONOLITHIC_WALL_CONDITION_H_INCLUDED)
14 #define KRATOS_ADJOINT_MONOLITHIC_WALL_CONDITION_H_INCLUDED
15 
16 // System includes
17 #include <string>
18 
19 // External includes
20 
21 // Project includes
22 #include "includes/condition.h"
23 #include "includes/define.h"
24 #include "includes/process_info.h"
25 #include "includes/serializer.h"
28 
29 // Application includes
30 
31 namespace Kratos
32 {
35 
38 
39 template <unsigned int TDim, unsigned int TNumNodes = TDim>
40 class KRATOS_API(FLUID_DYNAMICS_APPLICATION) AdjointMonolithicWallCondition : public Condition
41 {
42  class ThisExtensions : public AdjointExtensions
43  {
44  Condition* mpCondition;
45 
46  public:
47  explicit ThisExtensions(Condition* pCondition)
48  : mpCondition{pCondition}
49  {
50  }
51 
52  void GetFirstDerivativesVector(
53  std::size_t NodeId,
54  std::vector<IndirectScalar<double>>& rVector,
55  std::size_t Step) override;
56 
57  void GetSecondDerivativesVector(
58  std::size_t NodeId,
59  std::vector<IndirectScalar<double>>& rVector,
60  std::size_t Step) override;
61 
62  void GetAuxiliaryVector(
63  std::size_t NodeId,
64  std::vector<IndirectScalar<double>>& rVector,
65  std::size_t Step) override;
66 
67  void GetFirstDerivativesVariables(
68  std::vector<VariableData const*>& rVariables) const override;
69 
70  void GetSecondDerivativesVariables(
71  std::vector<VariableData const*>& rVariables) const override;
72 
73  void GetAuxiliaryVariables(
74  std::vector<VariableData const*>& rVariables) const override;
75  };
76 
77 public:
80 
82 
83  using NodeType = Node;
84 
86 
88 
90 
91  using VectorType = Vector;
92 
93  using MatrixType = Matrix;
94 
95  using IndexType = std::size_t;
96 
97  using SizeType = std::size_t;
98 
99  using EquationIdVectorType = std::vector<std::size_t>;
100 
101  using DofsVectorType = std::vector<Dof<double>::Pointer>;
102 
104 
105  static constexpr IndexType TBlockSize = (TDim + 1);
106 
107  static constexpr IndexType TFluidLocalSize = TBlockSize * TNumNodes;
108 
109  static constexpr IndexType TCoordsLocalSize = TDim * TNumNodes;
110 
114 
116  IndexType NewId = 0) : Condition(NewId)
117  {
118  }
119 
121  IndexType NewId,
122  const NodesArrayType& ThisNodes)
123  : Condition(NewId, ThisNodes)
124  {
125  }
126 
128  IndexType NewId,
129  GeometryType::Pointer pGeometry)
130  : Condition(NewId, pGeometry)
131  {
132  }
133 
135  IndexType NewId,
136  GeometryType::Pointer pGeometry,
137  PropertiesType::Pointer pProperties)
138  : Condition(NewId, pGeometry, pProperties)
139  {
140  }
141 
143  : Condition(rOther)
144  {
145  }
146 
147  ~AdjointMonolithicWallCondition() override = default;
148 
152 
155  {
156  Condition::operator=(rOther);
157  return *this;
158  }
159 
163 
164  Condition::Pointer Create(
165  IndexType NewId,
166  NodesArrayType const& ThisNodes,
167  PropertiesType::Pointer pProperties) const override
168  {
169  return Kratos::make_intrusive<AdjointMonolithicWallCondition>(
170  NewId, GetGeometry().Create(ThisNodes), pProperties);
171  }
172 
173  Condition::Pointer Create(
174  IndexType NewId,
175  GeometryType::Pointer pGeom,
176  PropertiesType::Pointer pProperties) const override
177  {
178  return Kratos::make_intrusive<AdjointMonolithicWallCondition>(NewId, pGeom, pProperties);
179  }
180 
181  Condition::Pointer Clone(
182  IndexType NewId,
183  NodesArrayType const& rThisNodes) const override
184  {
185  Condition::Pointer pNewCondition =
186  Create(NewId, GetGeometry().Create(rThisNodes), pGetProperties());
187 
188  pNewCondition->SetData(this->GetData());
189  pNewCondition->SetFlags(this->GetFlags());
190 
191  return pNewCondition;
192  }
193 
194  void Initialize(const ProcessInfo& rCurrentProcessInfo) override
195  {
196  this->SetValue(ADJOINT_EXTENSIONS, Kratos::make_shared<ThisExtensions>(this));
197  }
198 
199  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
200 
202  EquationIdVectorType& rResult,
203  const ProcessInfo& rCurrentProcessInfo) const override;
204 
206  DofsVectorType& ConditionDofList,
207  const ProcessInfo& CurrentProcessInfo) const override;
208 
210  Vector& Values,
211  int Step = 0) const override;
212 
213  void GetFirstDerivativesVector(
214  Vector& Values,
215  int Step = 0) const override;
216 
218  Vector& Values,
219  int Step = 0) const override;
220 
232  void CalculateLocalSystem(
233  MatrixType& rLeftHandSideMatrix,
234  VectorType& rRightHandSideVector,
235  const ProcessInfo& rCurrentProcessInfo) override;
236 
248  void CalculateLocalVelocityContribution(
249  MatrixType& rDampingMatrix,
250  VectorType& rRightHandSideVector,
251  const ProcessInfo& rCurrentProcessInfo) override;
252 
259  void CalculateLeftHandSide(
260  Matrix& rLeftHandSideMatrix,
261  const ProcessInfo& rCurrentProcessInfo) override;
262 
269  void CalculateFirstDerivativesLHS(
270  Matrix& rLeftHandSideMatrix,
271  const ProcessInfo& rCurrentProcessInfo) override;
272 
279  void CalculateSecondDerivativesLHS(
280  Matrix& rLeftHandSideMatrix,
281  const ProcessInfo& rCurrentProcessInfo) override;
282 
290  void CalculateSensitivityMatrix(
291  const Variable<array_1d<double, 3>>& rDesignVariable,
292  Matrix& rOutput,
293  const ProcessInfo& rCurrentProcessInfo) override;
294 
298 
300  std::string Info() const override
301  {
302  std::stringstream buffer;
303  buffer << "AdjointMonolithicWallCondition" << TDim << "D";
304  return buffer.str();
305  }
306 
308  void PrintInfo(std::ostream& rOStream) const override
309  {
310  rOStream << "AdjointMonolithicWallCondition";
311  }
312 
314  void PrintData(std::ostream& rOStream) const override
315  {
316  }
317 
319 
320 protected:
323 
324  virtual void ApplyNeumannCondition(
325  MatrixType& rLocalMatrix,
326  VectorType& rLocalVector,
327  const ProcessInfo& rCurrentProcessInfo) const;
328 
329  virtual void ApplyWallLaw(
330  MatrixType& rLocalMatrix,
331  VectorType& rLocalVector,
332  const ProcessInfo& rCurrentProcessInfo) const;
333 
334  virtual void ApplyNeumannConditionShapeDerivatives(
335  MatrixType& rLocalMatrix,
336  const ProcessInfo& rCurrentProcessInfo) const;
337 
338  virtual void ApplyWallLawStateDerivatives(
339  MatrixType& rLocalMatrix,
340  const ProcessInfo& rCurrentProcessInfo) const;
341 
342  virtual void ApplyWallLawShapeDerivatives(
343  MatrixType& rLocalMatrix,
344  const ProcessInfo& rCurrentProcessInfo) const;
345 
346  void CalculateData(
347  double& rArea,
348  double& rDetJ,
349  array_1d<double, TDim>& rUnitNormal) const;
350 
351  void CalculateDataShapeDerivatives(
352  BoundedVector<double, TCoordsLocalSize>& rAreaDerivatives,
353  BoundedVector<double, TCoordsLocalSize>& rDetJDerivatives,
354  BoundedMatrix<double, TCoordsLocalSize, TDim>& rUnitNormalDerivatives,
355  const double Area,
356  const double DetJ,
357  const array_1d<double, TDim>& rUnitNormal) const;
358 
360 
361 private:
364 
365  double CalculateDetJ(const double Area) const;
366 
367  void CalculateDetJShapeDerivatives(
368  BoundedVector<double, TCoordsLocalSize>& rDetJDerivatives,
369  const BoundedVector<double, TCoordsLocalSize>& rAreaDerivatives) const;
370 
371  friend class Serializer;
372 
373  void save(Serializer& rSerializer) const override
374  {
376  }
377 
378  void load(Serializer& rSerializer) override
379  {
381  }
382 
384 
385 }; // Class AdjointMonolithicWallCondition
386 
390 
392 template <unsigned int TDim, unsigned int TNumNodes>
393 inline std::istream& operator>>(
394  std::istream& rIStream,
396 {
397  return rIStream;
398 }
399 
401 template <unsigned int TDim, unsigned int TNumNodes>
402 inline std::ostream& operator<<(
403  std::ostream& rOStream,
405 {
406  rThis.PrintInfo(rOStream);
407  rOStream << std::endl;
408  rThis.PrintData(rOStream);
409 
410  return rOStream;
411 }
412 
414 
416 
417 } // namespace Kratos.
418 
419 #endif // KRATOS_ADJOINT_MONOLITHIC_WALL_CONDITION_H_INCLUDED
Interface extensions for adjoint elements and conditions.
Definition: adjoint_extensions.h:37
Definition: adjoint_monolithic_wall_condition.h:41
AdjointMonolithicWallCondition(IndexType NewId=0)
Definition: adjoint_monolithic_wall_condition.h:115
AdjointMonolithicWallCondition & operator=(AdjointMonolithicWallCondition const &rOther)
Assignment operator.
Definition: adjoint_monolithic_wall_condition.h:154
Condition::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
It creates a new condition pointer.
Definition: adjoint_monolithic_wall_condition.h:164
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: adjoint_monolithic_wall_condition.h:314
AdjointMonolithicWallCondition(IndexType NewId, const NodesArrayType &ThisNodes)
Definition: adjoint_monolithic_wall_condition.h:120
void GetSecondDerivativesVector(Vector &Values, int Step=0) const override
std::string Info() const override
Turn back information as a string.
Definition: adjoint_monolithic_wall_condition.h:300
Condition::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override
It creates a new condition pointer.
Definition: adjoint_monolithic_wall_condition.h:173
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: adjoint_monolithic_wall_condition.h:308
void GetValuesVector(Vector &Values, int Step=0) const override
Condition::Pointer Clone(IndexType NewId, NodesArrayType const &rThisNodes) const override
It creates a new condition pointer and clones the previous condition data.
Definition: adjoint_monolithic_wall_condition.h:181
AdjointMonolithicWallCondition(AdjointMonolithicWallCondition const &rOther)
Definition: adjoint_monolithic_wall_condition.h:142
std::size_t IndexType
Definition: adjoint_monolithic_wall_condition.h:95
AdjointMonolithicWallCondition(IndexType NewId, GeometryType::Pointer pGeometry)
Definition: adjoint_monolithic_wall_condition.h:127
void Initialize(const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_monolithic_wall_condition.h:194
AdjointMonolithicWallCondition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Definition: adjoint_monolithic_wall_condition.h:134
~AdjointMonolithicWallCondition() override=default
void GetDofList(DofsVectorType &ConditionDofList, const ProcessInfo &CurrentProcessInfo) const override
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(AdjointMonolithicWallCondition)
Base class for all Conditions.
Definition: condition.h:59
std::size_t SizeType
Definition: condition.h:94
std::vector< std::size_t > EquationIdVectorType
Definition: condition.h:98
std::vector< DofType::Pointer > DofsVectorType
Definition: condition.h:100
Condition & operator=(Condition const &rOther)
Assignment operator.
Definition: condition.h:181
Geometry base class.
Definition: geometry.h:71
This object defines an indexed object.
Definition: indexed_object.h:54
Wrapper for a function which behaves like an arithmetic type.
Definition: indirect_scalar.h:45
This class defines the node.
Definition: node.h:65
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
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
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_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
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
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
Internals::Matrix< double, AMatrix::dynamic, 1 > Vector
Definition: amatrix_interface.h:472
Internals::Matrix< double, AMatrix::dynamic, AMatrix::dynamic > Matrix
Definition: amatrix_interface.h:470
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
def SetValue(entity, variable, value)
Definition: coupling_interface_data.py:256
def load(f)
Definition: ode_solve.py:307