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_nodal_reaction_response_function.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: Martin Fusseder, https://github.com/MFusseder
10 //
11 
12 #pragma once
13 
14 
15 // System includes
16 
17 // External includes
18 
19 // Project includes
21 
22 
23 namespace Kratos
24 {
25 
28 
32 
36 
40 
44 
50 class KRATOS_API(STRUCTURAL_MECHANICS_APPLICATION) AdjointNodalReactionResponseFunction : public AdjointStructuralResponseFunction
51 {
52 public:
55 
57  typedef Node::Pointer PointTypePointer;
58  typedef matrix_column< Matrix > MatrixColumnType;
59  typedef matrix_row< Matrix > MatrixRowType;
60 
65 
69 
71  AdjointNodalReactionResponseFunction(ModelPart& rModelPart, Parameters ResponseSettings);
72 
75 
79 
83 
84  void InitializeSolutionStep() override;
85 
86  void FinalizeSolutionStep() override;
87 
89 
90  void CalculateGradient(const Element& rAdjointElement,
91  const Matrix& rResidualGradient,
92  Vector& rResponseGradient,
93  const ProcessInfo& rProcessInfo) override;
94 
95  void CalculateFirstDerivativesGradient(const Element& rAdjointElement,
96  const Matrix& rResidualGradient,
97  Vector& rResponseGradient,
98  const ProcessInfo& rProcessInfo) override;
99 
100  void CalculateFirstDerivativesGradient(const Condition& rAdjointCondition,
101  const Matrix& rResidualGradient,
102  Vector& rResponseGradient,
103  const ProcessInfo& rProcessInfo) override;
104 
105  void CalculateSecondDerivativesGradient(const Element& rAdjointElement,
106  const Matrix& rResidualGradient,
107  Vector& rResponseGradient,
108  const ProcessInfo& rProcessInfo) override;
109 
110  void CalculateSecondDerivativesGradient(const Condition& rAdjointCondition,
111  const Matrix& rResidualGradient,
112  Vector& rResponseGradient,
113  const ProcessInfo& rProcessInfo) override;
114 
115  void CalculatePartialSensitivity(Element& rAdjointElement,
116  const Variable<double>& rVariable,
117  const Matrix& rSensitivityMatrix,
118  Vector& rSensitivityGradient,
119  const ProcessInfo& rProcessInfo) override;
120 
121  void CalculatePartialSensitivity(Condition& rAdjointCondition,
122  const Variable<double>& rVariable,
123  const Matrix& rSensitivityMatrix,
124  Vector& rSensitivityGradient,
125  const ProcessInfo& rProcessInfo) override;
126 
127  void CalculatePartialSensitivity(Element& rAdjointElement,
128  const Variable<array_1d<double, 3>>& rVariable,
129  const Matrix& rSensitivityMatrix,
130  Vector& rSensitivityGradient,
131  const ProcessInfo& rProcessInfo) override;
132 
133  void CalculatePartialSensitivity(Condition& rAdjointCondition,
134  const Variable<array_1d<double, 3>>& rVariable,
135  const Matrix& rSensitivityMatrix,
136  Vector& rSensitivityGradient,
137  const ProcessInfo& rProcessInfo) override;
138 
139 
140  double CalculateValue(ModelPart& rModelPart) override;
141 
145 
149 
153 
157 
159 
160 protected:
163 
167 
171 
175 
179 
183 
187 
189 
190 private:
191 
194 
198 
199  std::string mTracedDisplacementLabel;
200  std::string mTracedReactionLabel;
201  PointTypePointer mpTracedNode;
202  GlobalPointersVector<Element> mpNeighborElements;
203  GlobalPointersVector<Condition> mpNeighborConditions;
204  bool mAdjustAdjointDisplacement = false;
205 
209 
213 
214 
215  void CalculateContributionToPartialSensitivity(Element& rAdjointElement,
216  const Matrix& rSensitivityMatrix,
217  Vector& rSensitivityGradient,
218  const ProcessInfo& rProcessInfo);
219 
220  void CalculateContributionToPartialSensitivity(Condition& rAdjointCondition,
221  const Matrix& rSensitivityMatrix,
222  Vector& rSensitivityGradient,
223  const ProcessInfo& rProcessInfo);
224 
225 
226  template <typename TObjectType>
227  size_t GetDofIndex(TObjectType& rAdjointObject, const ProcessInfo& rProcessInfo)
228  {
229  KRATOS_TRY;
230 
231  const Variable<double>& r_corresponding_adjoint_dof =
232  KratosComponents<Variable<double>>::Get(std::string("ADJOINT_") + mTracedDisplacementLabel);
233 
234  DofsVectorType dof_list;
235  rAdjointObject.GetDofList(dof_list, rProcessInfo);
236 
237  size_t dof_index = 0;
238  for(IndexType dof_it = 0; dof_it < dof_list.size(); ++dof_it)
239  {
240  if (dof_list[dof_it]->Id() == mpTracedNode->Id() &&
241  dof_list[dof_it]->GetVariable() == r_corresponding_adjoint_dof)
242  {
243  dof_index = dof_it;
244  break;
245  }
246  }
247  return dof_index;
248 
249  KRATOS_CATCH("");
250  }
251 
252  Vector GetColumnCopy(const Matrix& rMatrix, size_t ColumnIndex);
253 
254  Vector GetRowCopy(const Matrix& rMatrix, size_t RowIndex);
255 
256  std::string GetCorrespondingDisplacementLabel(std::string& rReactionLabel) const;
257 
258  void PerformResponseVariablesCheck();
259 
263 
267 
271 
273 
274 }; // Class AdjointNodalReactionResponseFunction
275 
277 
280 
284 
286 
287 } // namespace Kratos.
AdjointNodalReactionResponseFunction.
Definition: adjoint_nodal_reaction_response_function.h:51
matrix_column< Matrix > MatrixColumnType
Definition: adjoint_nodal_reaction_response_function.h:58
matrix_row< Matrix > MatrixRowType
Definition: adjoint_nodal_reaction_response_function.h:59
Node::Pointer PointTypePointer
Definition: adjoint_nodal_reaction_response_function.h:57
KRATOS_CLASS_POINTER_DEFINITION(AdjointNodalReactionResponseFunction)
Element::DofsVectorType DofsVectorType
Definition: adjoint_nodal_reaction_response_function.h:56
AdjointStructuralResponseFunction.
Definition: adjoint_structural_response_function.h:39
std::size_t IndexType
Definition: adjoint_structural_response_function.h:46
virtual void CalculateGradient(const Element &rAdjointElement, const Matrix &rResidualGradient, Vector &rResponseGradient, const ProcessInfo &rProcessInfo)
Calculate the local gradient w.r.t. primal solution.
Definition: adjoint_response_function.h:88
Base class for all Conditions.
Definition: condition.h:59
Base class for all Elements.
Definition: element.h:60
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
KratosComponents class encapsulates a lookup table for a family of classes in a generic way.
Definition: kratos_components.h:49
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
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
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