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.
residual_based_adjoint_steady_scheme.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:
11 //
12 
13 #if !defined(KRATOS_RESIDUAL_BASED_ADJOINT_STEADY_SCHEME_H_INCLUDED)
14 #define KRATOS_RESIDUAL_BASED_ADJOINT_STEADY_SCHEME_H_INCLUDED
15 
16 // System includes
17 #include <vector>
18 #include <string>
19 
20 // External includes
21 
22 // Project includes
23 #include "includes/define.h"
25 
26 namespace Kratos
27 {
30 
32 
42 template <class TSparseSpace, class TDenseSpace>
43 class ResidualBasedAdjointSteadyScheme : public ResidualBasedAdjointStaticScheme<TSparseSpace, TDenseSpace>
44 {
45 public:
48 
50 
52 
54 
56 
60 
62  explicit ResidualBasedAdjointSteadyScheme(AdjointResponseFunction::Pointer pResponseFunction)
63  : ResidualBasedAdjointStaticScheme<TSparseSpace, TDenseSpace>(pResponseFunction)
64  {
65  }
66 
69  {
70  }
71 
75 
79 
80  void CalculateSystemContributions(Element& rCurrentElement,
81  LocalSystemMatrixType& rLHS_Contribution,
82  LocalSystemVectorType& rRHS_Contribution,
83  Element::EquationIdVectorType& rEquationId,
84  const ProcessInfo& rCurrentProcessInfo) override
85  {
86  KRATOS_TRY;
87 
88  int thread_id = OpenMPUtils::ThisThread();
89  const auto& r_const_elem_ref = rCurrentElement;
90  rCurrentElement.CalculateFirstDerivativesLHS(rLHS_Contribution, rCurrentProcessInfo);
91 
92  if (rRHS_Contribution.size() != rLHS_Contribution.size1())
93  rRHS_Contribution.resize(rLHS_Contribution.size1(), false);
94 
95  this->mpResponseFunction->CalculateFirstDerivativesGradient(
96  rCurrentElement, rLHS_Contribution, rRHS_Contribution, rCurrentProcessInfo);
97 
98  noalias(rRHS_Contribution) = -rRHS_Contribution;
99 
100  // Calculate system contributions in residual form.
101  r_const_elem_ref.GetFirstDerivativesVector(this->mAdjointValues[thread_id]);
102  noalias(rRHS_Contribution) -= prod(rLHS_Contribution, this->mAdjointValues[thread_id]);
103 
104  r_const_elem_ref.EquationIdVector(rEquationId, rCurrentProcessInfo);
105 
106  KRATOS_CATCH("");
107  }
108 
109  void CalculateLHSContribution(Element& rCurrentElement,
110  LocalSystemMatrixType& rLHS_Contribution,
111  Element::EquationIdVectorType& rEquationId,
112  const ProcessInfo& rCurrentProcessInfo) override
113  {
114  KRATOS_TRY;
115 
116  rCurrentElement.CalculateFirstDerivativesLHS(rLHS_Contribution, rCurrentProcessInfo);
117  rCurrentElement.EquationIdVector(rEquationId, rCurrentProcessInfo);
118 
119  KRATOS_CATCH("");
120  }
121 
122  void CalculateSystemContributions(Condition& rCurrentCondition,
123  LocalSystemMatrixType& rLHS_Contribution,
124  LocalSystemVectorType& rRHS_Contribution,
125  Condition::EquationIdVectorType& rEquationId,
126  const ProcessInfo& rCurrentProcessInfo) override
127  {
128  KRATOS_TRY;
129 
130  int thread_id = OpenMPUtils::ThisThread();
131  const auto& r_const_cond_ref = rCurrentCondition;
132  rCurrentCondition.CalculateFirstDerivativesLHS(rLHS_Contribution, rCurrentProcessInfo);
133 
134  if (rRHS_Contribution.size() != rLHS_Contribution.size1())
135  rRHS_Contribution.resize(rLHS_Contribution.size1(), false);
136 
137  this->mpResponseFunction->CalculateFirstDerivativesGradient(
138  rCurrentCondition, rLHS_Contribution, rRHS_Contribution, rCurrentProcessInfo);
139 
140  noalias(rRHS_Contribution) = -rRHS_Contribution;
141 
142  // Calculate system contributions in residual form.
143  r_const_cond_ref.GetFirstDerivativesVector(this->mAdjointValues[thread_id]);
144  noalias(rRHS_Contribution) -= prod(rLHS_Contribution, this->mAdjointValues[thread_id]);
145 
146  r_const_cond_ref.EquationIdVector(rEquationId, rCurrentProcessInfo);
147 
148  KRATOS_CATCH("");
149  }
150 
151  void CalculateLHSContribution(Condition& rCurrentCondition,
152  LocalSystemMatrixType& rLHS_Contribution,
153  Condition::EquationIdVectorType& rEquationId,
154  const ProcessInfo& rCurrentProcessInfo) override
155  {
156  KRATOS_TRY;
157 
158  rCurrentCondition.CalculateFirstDerivativesLHS(rLHS_Contribution, rCurrentProcessInfo);
159  rCurrentCondition.EquationIdVector(rEquationId, rCurrentProcessInfo);
160 
161  KRATOS_CATCH("");
162  }
163 
167 
171 
175 
177 
178 protected:
181 
185 
189 
193 
197 
201 
205 
207 
208 private:
211 
215 
219 
223 
227 
231 
235 
237 
238 }; /* Class ResidualBasedAdjointSteadyScheme */
239 
241 
244 
246 
247 } /* namespace Kratos.*/
248 
249 #endif /* KRATOS_RESIDUAL_BASED_ADJOINT_STEADY_SCHEME_H_INCLUDED defined */
Base class for all Conditions.
Definition: condition.h:59
std::vector< std::size_t > EquationIdVectorType
Definition: condition.h:98
virtual void CalculateFirstDerivativesLHS(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:483
virtual void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const
Definition: condition.h:260
Base class for all Elements.
Definition: element.h:60
virtual void CalculateFirstDerivativesLHS(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:479
virtual void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const
Definition: element.h:258
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
static int ThisThread()
Wrapper for omp_get_thread_num().
Definition: openmp_utils.h:108
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
A scheme for static adjoint equations.
Definition: residual_based_adjoint_static_scheme.h:46
BaseType::LocalSystemMatrixType LocalSystemMatrixType
Definition: residual_based_adjoint_static_scheme.h:61
std::vector< LocalSystemVectorType > mAdjointValues
Definition: residual_based_adjoint_static_scheme.h:311
BaseType::LocalSystemVectorType LocalSystemVectorType
Definition: residual_based_adjoint_static_scheme.h:59
AdjointResponseFunction::Pointer mpResponseFunction
Definition: residual_based_adjoint_static_scheme.h:310
A scheme for steady adjoint equations.
Definition: residual_based_adjoint_steady_scheme.h:44
~ResidualBasedAdjointSteadyScheme() override
Destructor.
Definition: residual_based_adjoint_steady_scheme.h:68
ResidualBasedAdjointSteadyScheme(AdjointResponseFunction::Pointer pResponseFunction)
Constructor.
Definition: residual_based_adjoint_steady_scheme.h:62
void CalculateSystemContributions(Element &rCurrentElement, LocalSystemMatrixType &rLHS_Contribution, LocalSystemVectorType &rRHS_Contribution, Element::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo) override
This function is designed to be called in the builder and solver to introduce the selected time integ...
Definition: residual_based_adjoint_steady_scheme.h:80
void CalculateSystemContributions(Condition &rCurrentCondition, LocalSystemMatrixType &rLHS_Contribution, LocalSystemVectorType &rRHS_Contribution, Condition::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo) override
Functions totally analogous to the precedent but applied to the "condition" objects.
Definition: residual_based_adjoint_steady_scheme.h:122
void CalculateLHSContribution(Condition &rCurrentCondition, LocalSystemMatrixType &rLHS_Contribution, Condition::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo) override
Functions totally analogous to the precedent but applied to the "condition" objects.
Definition: residual_based_adjoint_steady_scheme.h:151
void CalculateLHSContribution(Element &rCurrentElement, LocalSystemMatrixType &rLHS_Contribution, Element::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo) override
This function is designed to calculate just the LHS contribution.
Definition: residual_based_adjoint_steady_scheme.h:109
BaseType::LocalSystemVectorType LocalSystemVectorType
Definition: residual_based_adjoint_steady_scheme.h:53
BaseType::LocalSystemMatrixType LocalSystemMatrixType
Definition: residual_based_adjoint_steady_scheme.h:55
KRATOS_CLASS_POINTER_DEFINITION(ResidualBasedAdjointSteadyScheme)
ResidualBasedAdjointStaticScheme< TSparseSpace, TDenseSpace > BaseType
Definition: residual_based_adjoint_steady_scheme.h:51
typename TDenseSpace::VectorType LocalSystemVectorType
Local system vector type definition.
Definition: scheme.h:80
typename TDenseSpace::MatrixType LocalSystemMatrixType
Local system matrix type definition.
Definition: scheme.h:77
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
AMatrix::MatrixProductExpression< TExpression1Type, TExpression2Type > prod(AMatrix::MatrixExpression< TExpression1Type, TCategory1 > const &First, AMatrix::MatrixExpression< TExpression2Type, TCategory2 > const &Second)
Definition: amatrix_interface.h:568
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484