13 #if !defined(KRATOS_RESIDUAL_BASED_ADJOINT_STEADY_SCHEME_H_INCLUDED)
14 #define KRATOS_RESIDUAL_BASED_ADJOINT_STEADY_SCHEME_H_INCLUDED
42 template <
class TSparseSpace,
class TDenseSpace>
89 const auto& r_const_elem_ref = rCurrentElement;
92 if (rRHS_Contribution.size() != rLHS_Contribution.size1())
93 rRHS_Contribution.resize(rLHS_Contribution.size1(),
false);
96 rCurrentElement, rLHS_Contribution, rRHS_Contribution, rCurrentProcessInfo);
98 noalias(rRHS_Contribution) = -rRHS_Contribution;
101 r_const_elem_ref.GetFirstDerivativesVector(this->
mAdjointValues[thread_id]);
104 r_const_elem_ref.EquationIdVector(rEquationId, rCurrentProcessInfo);
131 const auto& r_const_cond_ref = rCurrentCondition;
134 if (rRHS_Contribution.size() != rLHS_Contribution.size1())
135 rRHS_Contribution.resize(rLHS_Contribution.size1(),
false);
138 rCurrentCondition, rLHS_Contribution, rRHS_Contribution, rCurrentProcessInfo);
140 noalias(rRHS_Contribution) = -rRHS_Contribution;
143 r_const_cond_ref.GetFirstDerivativesVector(this->
mAdjointValues[thread_id]);
146 r_const_cond_ref.EquationIdVector(rEquationId, rCurrentProcessInfo);
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