13 #ifndef ADJOINT_RESPONSE_FUNCTION_H
14 #define ADJOINT_RESPONSE_FUNCTION_H
89 const Matrix& rResidualGradient,
93 KRATOS_ERROR <<
"Calling base class response function method.\n";
105 const Matrix& rResidualGradient,
106 Vector& rResponseGradient,
109 KRATOS_ERROR <<
"Calling base class response function method.\n";
121 const Matrix& rResidualGradient,
122 Vector& rResponseGradient,
125 KRATOS_ERROR <<
"Calling base class response function method.\n";
137 const Matrix& rResidualGradient,
138 Vector& rResponseGradient,
141 KRATOS_ERROR <<
"Calling base class response function method.\n";
153 const Matrix& rResidualGradient,
154 Vector& rResponseGradient,
157 KRATOS_ERROR <<
"Calling base class response function method.\n";
169 const Matrix& rResidualGradient,
170 Vector& rResponseGradient,
173 KRATOS_ERROR <<
"Calling base class response function method.\n";
187 const Matrix& rSensitivityMatrix,
188 Vector& rSensitivityGradient,
191 KRATOS_ERROR <<
"Calling base class response function method.\n";
205 const Matrix& rSensitivityMatrix,
206 Vector& rSensitivityGradient,
209 KRATOS_ERROR <<
"Calling base class response function method.\n";
223 const Matrix& rSensitivityMatrix,
224 Vector& rSensitivityGradient,
227 KRATOS_ERROR <<
"Calling base class response function method.\n";
241 const Matrix& rSensitivityMatrix,
242 Vector& rSensitivityGradient,
245 KRATOS_ERROR <<
"Calling base class response function method.\n";
A base class for adjoint response functions.
Definition: adjoint_response_function.h:39
virtual void CalculateFirstDerivativesGradient(const Condition &rAdjointCondition, const Matrix &rResidualGradient, Vector &rResponseGradient, const ProcessInfo &rProcessInfo)
Calculate the local gradient w.r.t. first derivatives of primal solution.
Definition: adjoint_response_function.h:136
virtual void CalculatePartialSensitivity(Element &rAdjointElement, const Variable< double > &rVariable, const Matrix &rSensitivityMatrix, Vector &rSensitivityGradient, const ProcessInfo &rProcessInfo)
Calculate the partial sensitivity w.r.t. design variable.
Definition: adjoint_response_function.h:185
virtual void Initialize()
Definition: adjoint_response_function.h:68
KRATOS_CLASS_POINTER_DEFINITION(AdjointResponseFunction)
virtual void InitializeSolutionStep()
Definition: adjoint_response_function.h:72
virtual void CalculatePartialSensitivity(Condition &rAdjointCondition, const Variable< array_1d< double, 3 >> &rVariable, const Matrix &rSensitivityMatrix, Vector &rSensitivityGradient, const ProcessInfo &rProcessInfo)
Calculate the partial sensitivity w.r.t. design variable.
Definition: adjoint_response_function.h:239
virtual double CalculateValue(ModelPart &rModelPart)=0
Calculate the scalar valued response function.
virtual void CalculatePartialSensitivity(Condition &rAdjointCondition, const Variable< double > &rVariable, const Matrix &rSensitivityMatrix, Vector &rSensitivityGradient, const ProcessInfo &rProcessInfo)
Calculate the partial sensitivity w.r.t. design variable.
Definition: adjoint_response_function.h:203
virtual void CalculatePartialSensitivity(Element &rAdjointElement, const Variable< array_1d< double, 3 >> &rVariable, const Matrix &rSensitivityMatrix, Vector &rSensitivityGradient, const ProcessInfo &rProcessInfo)
Calculate the partial sensitivity w.r.t. design variable.
Definition: adjoint_response_function.h:221
virtual void CalculateGradient(const Condition &rAdjointCondition, const Matrix &rResidualGradient, Vector &rResponseGradient, const ProcessInfo &rProcessInfo)
Calculate the local gradient w.r.t. primal solution.
Definition: adjoint_response_function.h:104
virtual void CalculateFirstDerivativesGradient(const Element &rAdjointElement, const Matrix &rResidualGradient, Vector &rResponseGradient, const ProcessInfo &rProcessInfo)
Calculate the local gradient w.r.t. first derivatives of primal solution.
Definition: adjoint_response_function.h:120
AdjointResponseFunction()
Constructor.
Definition: adjoint_response_function.h:51
virtual void CalculateSecondDerivativesGradient(const Element &rAdjointElement, const Matrix &rResidualGradient, Vector &rResponseGradient, const ProcessInfo &rProcessInfo)
Calculate the local gradient w.r.t. second derivatives of primal solution.
Definition: adjoint_response_function.h:152
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
virtual void CalculateSecondDerivativesGradient(const Condition &rAdjointCondition, const Matrix &rResidualGradient, Vector &rResponseGradient, const ProcessInfo &rProcessInfo)
Calculate the local gradient w.r.t. second derivatives of primal solution.
Definition: adjoint_response_function.h:168
virtual void FinalizeSolutionStep()
Definition: adjoint_response_function.h:76
virtual ~AdjointResponseFunction()
Destructor.
Definition: adjoint_response_function.h:56
Base class for all Conditions.
Definition: condition.h:59
Base class for all Elements.
Definition: element.h:60
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
#define KRATOS_ERROR
Definition: exception.h:161
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21