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_semi_analytic_base_condition.h
Go to the documentation of this file.
1 // KRATOS ___| | | |
2 // \___ \ __| __| | | __| __| | | __| _` | |
3 // | | | | | ( | | | | ( | |
4 // _____/ \__|_| \__,_|\___|\__|\__,_|_| \__,_|_| MECHANICS
5 //
6 // License: BSD License
7 // license: StructuralMechanicsApplication/license.txt
8 //
9 // Main authors: Armin Geiser, https://github.com/armingeiser
10 //
11 
12 // System includes
13 #pragma once
14 
15 // System includes
16 
17 // External includes
18 
19 // Project includes
20 #include "includes/condition.h"
21 
22 namespace Kratos
23 {
24 
27 
31 
35 
39 
43 
51 template <typename TPrimalCondition>
53  : public Condition
54 {
55 public:
60 
64 
66  : Condition(NewId),
67  mpPrimalCondition(Kratos::make_intrusive<TPrimalCondition>(NewId, pGetGeometry()))
68  {
69  }
70 
71  AdjointSemiAnalyticBaseCondition(IndexType NewId, GeometryType::Pointer pGeometry)
72  : Condition(NewId, pGeometry),
73  mpPrimalCondition(Kratos::make_intrusive<TPrimalCondition>(NewId, pGeometry))
74  {
75  }
76 
78  GeometryType::Pointer pGeometry,
79  PropertiesType::Pointer pProperties)
80  : Condition(NewId, pGeometry, pProperties),
81  mpPrimalCondition(Kratos::make_intrusive<TPrimalCondition>(NewId, pGeometry, pProperties))
82  {
83  }
84 
88 
92 
96 
97  Condition::Pointer Create(IndexType NewId,
98  NodesArrayType const& ThisNodes,
99  PropertiesType::Pointer pProperties) const override
100  {
101  return Kratos::make_intrusive<AdjointSemiAnalyticBaseCondition<TPrimalCondition>>(
102  NewId, GetGeometry().Create(ThisNodes), pProperties);
103  }
104 
105  Condition::Pointer Create(IndexType NewId,
106  GeometryType::Pointer pGeometry,
107  PropertiesType::Pointer pProperties) const override
108  {
109  return Kratos::make_intrusive<AdjointSemiAnalyticBaseCondition<TPrimalCondition>>(
110  NewId, pGeometry, pProperties);
111  }
112 
113  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo ) const override;
114 
115  void GetDofList(DofsVectorType& ElementalDofList, const ProcessInfo& rCurrentProcessInfo ) const override;
116 
118  {
119  return mpPrimalCondition->GetIntegrationMethod();
120  }
121 
122  void GetValuesVector(Vector& rValues, int Step = 0 ) const override;
123 
124  void Initialize(const ProcessInfo& rCurrentProcessInfo) override
125  {
126  mpPrimalCondition->Initialize(rCurrentProcessInfo);
127  }
128 
129  void ResetConstitutiveLaw() override
130  {
131  mpPrimalCondition->ResetConstitutiveLaw();
132  }
133 
134  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override
135  {
136  mpPrimalCondition->InitializeSolutionStep(rCurrentProcessInfo);
137  }
138 
139  void InitializeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override
140  {
141  mpPrimalCondition->InitializeNonLinearIteration(rCurrentProcessInfo);
142  }
143 
144  void FinalizeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override
145  {
146  mpPrimalCondition->FinalizeNonLinearIteration(rCurrentProcessInfo);
147  }
148 
149  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override
150  {
151  mpPrimalCondition->FinalizeSolutionStep(rCurrentProcessInfo);
152  }
153 
154  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,
155  VectorType& rRightHandSideVector,
156  const ProcessInfo& rCurrentProcessInfo) override
157  {
158  mpPrimalCondition->CalculateLocalSystem(rLeftHandSideMatrix,
159  rRightHandSideVector,
160  rCurrentProcessInfo);
161  }
162 
163  void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix,
164  const ProcessInfo& rCurrentProcessInfo) override
165  {
166  mpPrimalCondition->CalculateLeftHandSide(rLeftHandSideMatrix,
167  rCurrentProcessInfo);
168  }
169 
170  void CalculateRightHandSide(VectorType& rRightHandSideVector,
171  const ProcessInfo& rCurrentProcessInfo) override
172  {
173  mpPrimalCondition->CalculateRightHandSide(rRightHandSideVector,
174  rCurrentProcessInfo);
175  }
176 
178  VectorType& rRightHandSideVector,
179  const ProcessInfo& rCurrentProcessInfo) override
180  {
181  mpPrimalCondition->CalculateFirstDerivativesContributions(rLeftHandSideMatrix,
182  rRightHandSideVector,
183  rCurrentProcessInfo);
184  }
185 
186  void CalculateFirstDerivativesLHS(MatrixType& rLeftHandSideMatrix,
187  const ProcessInfo& rCurrentProcessInfo) override
188  {
189  mpPrimalCondition->CalculateFirstDerivativesLHS(rLeftHandSideMatrix,
190  rCurrentProcessInfo);
191  }
192 
193  void CalculateFirstDerivativesRHS(VectorType& rRightHandSideVector,
194  const ProcessInfo& rCurrentProcessInfo) override
195  {
196  mpPrimalCondition->CalculateFirstDerivativesRHS(rRightHandSideVector,
197  rCurrentProcessInfo);
198  }
199 
201  VectorType& rRightHandSideVector,
202  const ProcessInfo& rCurrentProcessInfo) override
203  {
204  mpPrimalCondition->CalculateSecondDerivativesContributions(rLeftHandSideMatrix,
205  rRightHandSideVector,
206  rCurrentProcessInfo);
207  }
208 
209  void CalculateSecondDerivativesLHS(MatrixType& rLeftHandSideMatrix,
210  const ProcessInfo& rCurrentProcessInfo) override
211  {
212  mpPrimalCondition->CalculateSecondDerivativesLHS(rLeftHandSideMatrix,
213  rCurrentProcessInfo);
214  }
215 
216  void CalculateSecondDerivativesRHS(VectorType& rRightHandSideVector,
217  const ProcessInfo& rCurrentProcessInfo) override
218  {
219  mpPrimalCondition->CalculateSecondDerivativesRHS(rRightHandSideVector,
220  rCurrentProcessInfo);
221  }
222 
223  void CalculateMassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo) override
224  {
225  mpPrimalCondition->CalculateMassMatrix(rMassMatrix, rCurrentProcessInfo);
226  }
227 
228  void CalculateDampingMatrix(MatrixType& rDampingMatrix, const ProcessInfo& rCurrentProcessInfo) override
229  {
230  mpPrimalCondition->CalculateDampingMatrix(rDampingMatrix, rCurrentProcessInfo);
231  }
232 
233  void Calculate(const Variable<double >& rVariable,
234  double& Output,
235  const ProcessInfo& rCurrentProcessInfo) override
236  {
237  KRATOS_ERROR << "Calculate of the adjoint base condition is called!" << std::endl;
238  }
239 
240  void Calculate(const Variable< array_1d<double,3> >& rVariable,
241  array_1d<double,3>& Output,
242  const ProcessInfo& rCurrentProcessInfo) override
243  {
244  KRATOS_ERROR << "Calculate of the adjoint base condition is called!" << std::endl;
245  }
246 
247  void Calculate(const Variable<Vector >& rVariable,
248  Vector& Output,
249  const ProcessInfo& rCurrentProcessInfo) override
250  {
251  KRATOS_ERROR << "Calculate of the adjoint base condition is called!" << std::endl;
252  }
253 
254  void Calculate(const Variable<Matrix >& rVariable,
255  Matrix& Output,
256  const ProcessInfo& rCurrentProcessInfo) override
257  {
258  KRATOS_ERROR << "Calculate of the adjoint base condition is called!" << std::endl;
259  }
260 
261 
262  void CalculateOnIntegrationPoints(const Variable<double>& rVariable,
263  std::vector<double>& rOutput,
264  const ProcessInfo& rCurrentProcessInfo) override;
265 
267  std::vector< array_1d<double, 3 > >& Output,
268  const ProcessInfo& rCurrentProcessInfo) override;
269 
271  std::vector< Vector >& Output,
272  const ProcessInfo& rCurrentProcessInfo) override
273  {
274  KRATOS_ERROR << "CalculateOnIntegrationPoints of the adjoint base condition is called!" << std::endl;
275  }
276 
278  std::vector< Matrix >& Output,
279  const ProcessInfo& rCurrentProcessInfo) override
280  {
281  KRATOS_ERROR << "CalculateOnIntegrationPoints of the adjoint base condition is called!" << std::endl;
282  }
283 
284  int Check( const ProcessInfo& rCurrentProcessInfo ) const override;
285 
289  void CalculateSensitivityMatrix(const Variable<double>& rDesignVariable,
290  Matrix& rOutput,
291  const ProcessInfo& rCurrentProcessInfo) override;
292 
296  void CalculateSensitivityMatrix(const Variable<array_1d<double,3> >& rDesignVariable,
297  Matrix& rOutput,
298  const ProcessInfo& rCurrentProcessInfo) override;
299 
300  Condition::Pointer pGetPrimalCondition()
301  {
302  return mpPrimalCondition;
303  }
304 
305  const Condition::Pointer pGetPrimalCondition() const
306  {
307  return mpPrimalCondition;
308  }
309 
313 
314 
318 
319 
323 
324 
328 
329 
331 
332 protected:
335 
336 
340 
341  Condition::Pointer mpPrimalCondition;
342 
346 
347 
351 
355 
356 
360 
361 
365 
367 
368 private:
371 
372 
376 
377 
378 
382 
386 
390  double GetPerturbationSize(const Variable<double>& rDesignVariable, const ProcessInfo& rCurrentProcessInfo) const;
391 
395  double GetPerturbationSize(const Variable<array_1d<double,3>>& rDesignVariable, const ProcessInfo& rCurrentProcessInfo) const;
396 
403  virtual double GetPerturbationSizeModificationFactor(const Variable<double>& rVariable) const;
404 
411  virtual double GetPerturbationSizeModificationFactor(const Variable<array_1d<double,3>>& rDesignVariable) const;
412 
413 
417 
418 
422 
426 
427  friend class Serializer;
428 
429  void save( Serializer& rSerializer ) const override
430  {
432  rSerializer.save("mpPrimalCondition", mpPrimalCondition);
433  }
434 
435  void load( Serializer& rSerializer ) override
436  {
438  rSerializer.load("mpPrimalCondition", mpPrimalCondition);
439  }
440 
444 
446 
447 }; // Class AdjointSemiAnalyticBaseCondition
448 
452 
453 
457 
459 
460 } // namespace Kratos.
461 
462 
AdjointSemiAnalyticBaseCondition.
Definition: adjoint_semi_analytic_base_condition.h:54
void ResetConstitutiveLaw() override
Definition: adjoint_semi_analytic_base_condition.h:129
AdjointSemiAnalyticBaseCondition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Definition: adjoint_semi_analytic_base_condition.h:77
const Condition::Pointer pGetPrimalCondition() const
Definition: adjoint_semi_analytic_base_condition.h:305
AdjointSemiAnalyticBaseCondition(IndexType NewId=0)
Definition: adjoint_semi_analytic_base_condition.h:65
void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:163
void Calculate(const Variable< array_1d< double, 3 > > &rVariable, array_1d< double, 3 > &Output, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:240
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:170
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:154
void CalculateSecondDerivativesRHS(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:216
Condition::Pointer mpPrimalCondition
Definition: adjoint_semi_analytic_base_condition.h:341
void CalculateDampingMatrix(MatrixType &rDampingMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:228
void CalculateFirstDerivativesContributions(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:177
void CalculateSecondDerivativesLHS(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:209
void Initialize(const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:124
void CalculateOnIntegrationPoints(const Variable< Matrix > &rVariable, std::vector< Matrix > &Output, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:277
Condition::Pointer pGetPrimalCondition()
Definition: adjoint_semi_analytic_base_condition.h:300
void CalculateFirstDerivativesLHS(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:186
Condition::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
It creates a new condition pointer.
Definition: adjoint_semi_analytic_base_condition.h:97
void CalculateOnIntegrationPoints(const Variable< double > &rVariable, std::vector< double > &rOutput, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.cpp:146
void Calculate(const Variable< double > &rVariable, double &Output, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:233
void GetValuesVector(Vector &rValues, int Step=0) const override
Definition: adjoint_semi_analytic_base_condition.cpp:96
int Check(const ProcessInfo &rCurrentProcessInfo) const override
Definition: adjoint_semi_analytic_base_condition.cpp:295
AdjointSemiAnalyticBaseCondition(IndexType NewId, GeometryType::Pointer pGeometry)
Definition: adjoint_semi_analytic_base_condition.h:71
void CalculateSecondDerivativesContributions(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:200
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
Definition: adjoint_semi_analytic_base_condition.cpp:33
void GetDofList(DofsVectorType &ElementalDofList, const ProcessInfo &rCurrentProcessInfo) const override
Definition: adjoint_semi_analytic_base_condition.cpp:64
void CalculateSensitivityMatrix(const Variable< double > &rDesignVariable, Matrix &rOutput, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.cpp:176
void CalculateMassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:223
void InitializeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:139
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(AdjointSemiAnalyticBaseCondition)
void InitializeSolutionStep(const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:134
void CalculateOnIntegrationPoints(const Variable< Vector > &rVariable, std::vector< Vector > &Output, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:270
void Calculate(const Variable< Matrix > &rVariable, Matrix &Output, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:254
IntegrationMethod GetIntegrationMethod() const override
Definition: adjoint_semi_analytic_base_condition.h:117
void FinalizeSolutionStep(const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:149
void FinalizeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:144
void CalculateFirstDerivativesRHS(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:193
void Calculate(const Variable< Vector > &rVariable, Vector &Output, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_semi_analytic_base_condition.h:247
Condition::Pointer Create(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties) const override
It creates a new condition pointer.
Definition: adjoint_semi_analytic_base_condition.h:105
Base class for all Conditions.
Definition: condition.h:59
std::vector< std::size_t > EquationIdVectorType
Definition: condition.h:98
std::vector< DofType::Pointer > DofsVectorType
Definition: condition.h:100
std::size_t IndexType
Definition: flags.h:74
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
IntegrationMethod
Definition: geometry_data.h:76
virtual Pointer Create(PointsArrayType const &rThisPoints) const
Creates a new geometry pointer.
Definition: geometry.h:813
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
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
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
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
#define KRATOS_ERROR
Definition: exception.h:161
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
intrusive_ptr< C > make_intrusive(Args &&...args)
Definition: smart_pointers.h:36
def load(f)
Definition: ode_solve.py:307