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.
steady_scalar_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: Suneth Warnakulasuriya
11 //
12 
13 #if !defined(KRATOS_STEADY_SCALAR_TRANSPORT_SCHEME)
14 #define KRATOS_STEADY_SCALAR_TRANSPORT_SCHEME
15 
16 // System includes
17 #include <iomanip>
18 #include <sstream>
19 
20 // Project includes
21 #include "includes/define.h"
22 #include "includes/model_part.h"
24 #include "utilities/openmp_utils.h"
25 
26 // Application includes
28 
29 namespace Kratos
30 {
33 
34 template <class TSparseSpace, class TDenseSpace>
36  : public Scheme<TSparseSpace, TDenseSpace>
37 {
38 public:
41 
43 
45 
47 
49 
51 
53 
55 
59 
61  const double RelaxationFactor)
62  : mRelaxationFactor(RelaxationFactor)
63  {
64  const int num_threads = OpenMPUtils::GetNumThreads();
66 
67  mpDofUpdater = Kratos::make_unique<DofUpdaterType>(mRelaxationFactor);
68  }
69 
70  ~SteadyScalarScheme() override = default;
71 
75 
76  void Update(
77  ModelPart& rModelPart,
78  DofsArrayType& rDofSet,
80  TSystemVectorType& rDx,
81  TSystemVectorType& rb) override
82  {
83  KRATOS_TRY;
84 
85  mpDofUpdater->UpdateDofs(rDofSet, rDx);
86 
87  KRATOS_CATCH("");
88  }
89 
90  void Clear() override
91  {
92  this->mpDofUpdater->Clear();
93  }
94 
96  Element& rElement,
97  LocalSystemMatrixType& rLHS_Contribution,
98  LocalSystemVectorType& rRHS_Contribution,
99  Element::EquationIdVectorType& rEquationIdVector,
100  const ProcessInfo& rCurrentProcessInfo) override
101  {
102  KRATOS_TRY;
103 
104  const int k = OpenMPUtils::ThisThread();
105 
106  this->CalculateDampingSystem(rElement, rLHS_Contribution, rRHS_Contribution,
107  rEquationIdVector, rCurrentProcessInfo, k);
108 
109  if (mDampingMatrix[k].size1() != 0) {
110  noalias(rLHS_Contribution) += mDampingMatrix[k];
111  }
112 
113  KRATOS_CATCH("");
114  }
115 
117  Condition& rCondition,
118  LocalSystemMatrixType& rLHS_Contribution,
119  LocalSystemVectorType& rRHS_Contribution,
120  Condition::EquationIdVectorType& rEquationIdVector,
121  const ProcessInfo& rCurrentProcessInfo) override
122  {
123  KRATOS_TRY;
124 
125  const int k = OpenMPUtils::ThisThread();
126 
127  this->CalculateDampingSystem(rCondition, rLHS_Contribution, rRHS_Contribution,
128  rEquationIdVector, rCurrentProcessInfo, k);
129 
130  if (mDampingMatrix[k].size1() != 0) {
131  noalias(rLHS_Contribution) += mDampingMatrix[k];
132  }
133 
134  KRATOS_CATCH("");
135  }
136 
138  Element& rElement,
139  LocalSystemVectorType& rRHS_Contribution,
140  Element::EquationIdVectorType& rEquationIdVector,
141  const ProcessInfo& rCurrentProcessInfo) override
142  {
143  KRATOS_TRY;
144 
145  const int k = OpenMPUtils::ThisThread();
146 
147  this->CalculateDampingSystem(rElement, mDampingMatrix[k], rRHS_Contribution,
148  rEquationIdVector, rCurrentProcessInfo, k);
149 
150  KRATOS_CATCH("");
151  }
152 
154  Condition& rCondition,
155  LocalSystemVectorType& rRHS_Contribution,
156  Element::EquationIdVectorType& rEquationIdVector,
157  const ProcessInfo& rCurrentProcessInfo) override
158  {
159  KRATOS_TRY;
160 
161  const int k = OpenMPUtils::ThisThread();
162 
163  this->CalculateDampingSystem(rCondition, mDampingMatrix[k], rRHS_Contribution,
164  rEquationIdVector, rCurrentProcessInfo, k);
165 
166  KRATOS_CATCH("");
167  }
168 
172 
174 
175  std::string Info() const override
176  {
177  std::stringstream msg;
178  msg << "Using generic residual based steady scalar transport scheme "
179  "with\n"
180  << " Relaxation factor : " << this->mRelaxationFactor;
181  return msg.str();
182  }
183 
185 
186 protected:
189 
190  std::vector<LocalSystemMatrixType> mDampingMatrix;
192 
196 
197  template <class TItemType>
199  TItemType& rItem,
200  LocalSystemMatrixType& rLHS_Contribution,
201  LocalSystemVectorType& rRHS_Contribution,
202  typename TItemType::EquationIdVectorType& rEquationIdVector,
203  const ProcessInfo& rCurrentProcessInfo,
204  const int ThreadId)
205  {
206  KRATOS_TRY;
207 
208  rItem.CalculateLocalSystem(rLHS_Contribution, rRHS_Contribution, rCurrentProcessInfo);
209  rItem.CalculateLocalVelocityContribution(
210  mDampingMatrix[ThreadId], rRHS_Contribution, rCurrentProcessInfo);
211  rItem.EquationIdVector(rEquationIdVector, rCurrentProcessInfo);
212 
213  KRATOS_CATCH("");
214  }
215 
217 
218 private:
221 
222  using DofUpdaterType = RelaxedDofUpdater<TSparseSpace>;
223  using DofUpdaterPointerType = typename DofUpdaterType::UniquePointer;
224 
225  DofUpdaterPointerType mpDofUpdater;
226 
228 };
229 
231 
232 } // namespace Kratos
233 
234 #endif /* KRATOS_STEADY_SCALAR_TRANSPORT_SCHEME defined */
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