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.
sensitivity_builder.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: Martin Fusseder, https://github.com/MFusseder
11 // Michael Andre, https://github.com/msandre
12 // Suneth Warnakulasuriya, https://github.com/sunethwarna
13 //
14 
15 #pragma once
16 
17 // System includes
18 #include <string>
19 #include <vector>
20 #include <unordered_map>
21 
22 // External includes
23 
24 // Project includes
25 #include "includes/define.h"
26 #include "includes/model_part.h"
29 
30 namespace Kratos
31 {
34 
35 class KRATOS_API(KRATOS_CORE) SensitivityBuilder
36 {
37 public:
40 
42 
61  template <class TDataType>
63  {
64  const Variable<TDataType>* pDesignVariable = nullptr;
65  const Variable<TDataType>* pOutputVariable = nullptr;
66 
67  explicit SensitivityVariables(const std::string& rName)
68  {
69  KRATOS_TRY;
70  const std::string output_suffix = "_SENSITIVITY";
71  pDesignVariable = &KratosComponents<Variable<TDataType>>::Get(rName);
72  if (rName.size() > output_suffix.size() &&
73  std::equal(output_suffix.rbegin(), output_suffix.rend(), rName.rbegin()))
74  {
75  pOutputVariable = pDesignVariable;
76  }
77  else
78  {
79  pOutputVariable =
80  &KratosComponents<Variable<TDataType>>::Get(rName + output_suffix);
81  }
82  KRATOS_CATCH("");
83  }
84  };
85 
86 
98  template <class TDataType>
99  using THomogeneousSensitivityVariables = std::vector<SensitivityVariables<TDataType>>;
100 
105  using TSensitivityVariables = std::tuple<
108  >;
109 
113 
116  ModelPart& rModelPart,
117  AdjointResponseFunction::Pointer pResponseFunction);
118 
121  ModelPart& rModelPart,
122  AdjointResponseFunction::Pointer pResponseFunction,
123  SensitivityBuilderScheme::Pointer pSensitivityBuilderScheme);
124 
128 
140  void SetResponseFunction(AdjointResponseFunction::Pointer pResponseFunction);
141 
142  void Initialize();
143 
144  void InitializeSolutionStep();
145 
146  void UpdateSensitivities();
147 
148  void FinalizeSolutionStep();
149 
150  void Finalize();
151 
152  void Clear();
153 
155  void ClearFlags();
156 
158  void ClearSensitivities();
159 
163 
164  static void CalculateNodalSolutionStepSensitivities(
165  const std::vector<std::string>& rVariables,
166  ModelPart& rModelPart,
167  AdjointResponseFunction& rResponseFunction,
168  double ScalingFactor);
169 
170  static void CalculateNodalSolutionStepSensitivities(
171  const TSensitivityVariables& rVariables,
172  ModelPart& rModelPart,
173  AdjointResponseFunction& rResponseFunction,
174  SensitivityBuilderScheme& rSensitivityBuilderScheme,
175  double ScalingFactor);
176 
177  static void CalculateNonHistoricalSensitivities(
178  const std::vector<std::string>& rVariables,
180  AdjointResponseFunction& rResponseFunction,
181  const ProcessInfo& rProcessInfo,
182  double ScalingFactor);
183 
184  static void CalculateNonHistoricalSensitivities(
185  const TSensitivityVariables& rVariables,
187  AdjointResponseFunction& rResponseFunction,
188  SensitivityBuilderScheme& rSensitivityBuilderScheme,
189  const ProcessInfo& rProcessInfo,
190  double ScalingFactor);
191 
192  static void CalculateNonHistoricalSensitivities(
193  const std::vector<std::string>& rVariables,
195  AdjointResponseFunction& rResponseFunction,
196  const ProcessInfo& rProcessInfo,
197  double ScalingFactor);
198 
199  static void CalculateNonHistoricalSensitivities(
200  const TSensitivityVariables& rVariables,
202  AdjointResponseFunction& rResponseFunction,
203  SensitivityBuilderScheme& rSensitivityBuilderScheme,
204  const ProcessInfo& rProcessInfo,
205  double ScalingFactor);
206 
207 private:
210 
211  ModelPart* mpModelPart = nullptr;
212  ModelPart* mpSensitivityModelPart = nullptr;
213  AdjointResponseFunction::Pointer mpResponseFunction;
214  SensitivityBuilderScheme::Pointer mpSensitivityBuilderScheme;
215 
216  TSensitivityVariables mNodalSolutionStepSensitivityVariablesList;
217  TSensitivityVariables mElementDataValueSensitivityVariablesList;
218  TSensitivityVariables mConditionDataValueSensitivityVariablesList;
219 
220  std::string mBuildMode = "static";
221  bool mNodalSolutionStepSensitivityCalculationIsThreadSafe = false;
222 
226 
227  static TSensitivityVariables GetVariableLists(const std::vector<std::string>& rVariableNames);
228 
230 };
231 
233 
234 } /* namespace Kratos.*/
A base class for adjoint response functions.
Definition: adjoint_response_function.h:39
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
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
Definition: sensitivity_builder.h:36
std::tuple< THomogeneousSensitivityVariables< double >, THomogeneousSensitivityVariables< array_1d< double, 3 > > > TSensitivityVariables
This type holds list of homogeneous variable lists.
Definition: sensitivity_builder.h:108
KRATOS_CLASS_POINTER_DEFINITION(SensitivityBuilder)
std::vector< SensitivityVariables< TDataType > > THomogeneousSensitivityVariables
This type holds list of variables for the sensitivity analysis.
Definition: sensitivity_builder.h:99
Scheme used in the Sensitivity Builder.
Definition: sensitivity_builder_scheme.h:53
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
#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
Contains the sensitivity design and output variables.
Definition: sensitivity_builder.h:63
SensitivityVariables(const std::string &rName)
Definition: sensitivity_builder.h:67