13 #if !defined(KRATOS_STEADY_SCALAR_TRANSPORT_SCHEME)
14 #define KRATOS_STEADY_SCALAR_TRANSPORT_SCHEME
34 template <
class TSparseSpace,
class TDenseSpace>
36 :
public Scheme<TSparseSpace, TDenseSpace>
61 const double RelaxationFactor)
64 const int num_threads = OpenMPUtils::GetNumThreads();
85 mpDofUpdater->UpdateDofs(rDofSet, rDx);
92 this->mpDofUpdater->Clear();
107 rEquationIdVector, rCurrentProcessInfo,
k);
128 rEquationIdVector, rCurrentProcessInfo,
k);
148 rEquationIdVector, rCurrentProcessInfo,
k);
164 rEquationIdVector, rCurrentProcessInfo,
k);
175 std::string
Info()
const override
177 std::stringstream msg;
178 msg <<
"Using generic residual based steady scalar transport scheme "
197 template <
class TItemType>
202 typename TItemType::EquationIdVectorType& rEquationIdVector,
208 rItem.CalculateLocalSystem(rLHS_Contribution, rRHS_Contribution, rCurrentProcessInfo);
209 rItem.CalculateLocalVelocityContribution(
210 mDampingMatrix[ThreadId], rRHS_Contribution, rCurrentProcessInfo);
211 rItem.EquationIdVector(rEquationIdVector, rCurrentProcessInfo);
223 using DofUpdaterPointerType =
typename DofUpdaterType::UniquePointer;
225 DofUpdaterPointerType mpDofUpdater;
Base class for all Conditions.
Definition: condition.h:59
std::vector< std::size_t > EquationIdVectorType
Definition: condition.h:98
Base class for all Elements.
Definition: element.h:60
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
static int ThisThread()
Wrapper for omp_get_thread_num().
Definition: openmp_utils.h:108
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
Utility class to update the values of degree of freedom (Dof) variables after solving the system.
Definition: relaxed_dof_updater.h:46
This class provides the implementation of the basic tasks that are needed by the solution strategy.
Definition: scheme.h:56
typename TSparseSpace::MatrixType TSystemMatrixType
Matrix type definition.
Definition: scheme.h:71
typename TSparseSpace::VectorType TSystemVectorType
Vector type definition.
Definition: scheme.h:74
typename TDenseSpace::VectorType LocalSystemVectorType
Local system vector type definition.
Definition: scheme.h:80
ModelPart::DofsArrayType DofsArrayType
DoF array type definition.
Definition: scheme.h:86
typename TDenseSpace::MatrixType LocalSystemMatrixType
Local system matrix type definition.
Definition: scheme.h:77
Definition: steady_scalar_scheme.h:37
void Clear() override
Liberate internal storage.
Definition: steady_scalar_scheme.h:90
~SteadyScalarScheme() override=default
std::vector< LocalSystemMatrixType > mDampingMatrix
Definition: steady_scalar_scheme.h:190
void CalculateRHSContribution(Element &rElement, LocalSystemVectorType &rRHS_Contribution, Element::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo) override
This function is designed to calculate just the RHS contribution.
Definition: steady_scalar_scheme.h:137
std::string Info() const override
Turn back information as a string.
Definition: steady_scalar_scheme.h:175
void CalculateSystemContributions(Element &rElement, LocalSystemMatrixType &rLHS_Contribution, LocalSystemVectorType &rRHS_Contribution, Element::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo) override
This function is designed to be called in the builder and solver to introduce the selected time integ...
Definition: steady_scalar_scheme.h:95
SteadyScalarScheme(const double RelaxationFactor)
Definition: steady_scalar_scheme.h:60
KRATOS_CLASS_POINTER_DEFINITION(SteadyScalarScheme)
void Update(ModelPart &rModelPart, DofsArrayType &rDofSet, TSystemMatrixType &rA, TSystemVectorType &rDx, TSystemVectorType &rb) override
Performing the update of the solution.
Definition: steady_scalar_scheme.h:76
void CalculateDampingSystem(TItemType &rItem, LocalSystemMatrixType &rLHS_Contribution, LocalSystemVectorType &rRHS_Contribution, typename TItemType::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo, const int ThreadId)
Definition: steady_scalar_scheme.h:198
void CalculateRHSContribution(Condition &rCondition, LocalSystemVectorType &rRHS_Contribution, Element::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo) override
Functions totally analogous to the precedent but applied to the "condition" objects.
Definition: steady_scalar_scheme.h:153
double mRelaxationFactor
Definition: steady_scalar_scheme.h:191
void CalculateSystemContributions(Condition &rCondition, LocalSystemMatrixType &rLHS_Contribution, LocalSystemVectorType &rRHS_Contribution, Condition::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo) override
Functions totally analogous to the precedent but applied to the "condition" objects.
Definition: steady_scalar_scheme.h:116
#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
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
int k
Definition: quadrature.py:595
def num_threads
Definition: script.py:75