14 #if !defined(KRATOS_EXPLICIT_RUNGE_KUTTA_4_EULERIAN_CONVDIFF_STRATEGY)
15 #define KRATOS_EXPLICIT_RUNGE_KUTTA_4_EULERIAN_CONVDIFF_STRATEGY
55 template <
class TSparseSpace,
class TDenseSpace>
86 :
BaseType(rModelPart, ThisParameters)
98 typename ExplicitBuilderType::Pointer pExplicitBuilder,
100 int RebuildLevel = 0)
113 int RebuildLevel = 0)
146 auto& r_process_info = r_model_part.GetProcessInfo();
147 ConvectionDiffusionSettings::Pointer p_settings = r_process_info[CONVECTION_DIFFUSION_SETTINGS];
148 auto& r_settings = *p_settings;
153 if (r_process_info[OSS_SWITCH]) {
154 for (
auto& r_node : r_model_part.GetCommunicator().LocalMesh().Nodes()) {
155 r_node.SetValue(r_settings.GetProjectionVariable(), 0.0);
167 std::string
Info()
const override
169 return "ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion";
215 auto& r_process_info = r_model_part.GetProcessInfo();
218 if (r_process_info[OSS_SWITCH] == 1) {
219 CalculateOSSNodalProjections();
235 auto& r_process_info = r_model_part.GetProcessInfo();
238 if (r_process_info[OSS_SWITCH] == 1) {
239 CalculateOSSNodalProjections();
255 auto& r_process_info = r_model_part.GetProcessInfo();
256 if (r_process_info.GetValue(RUNGE_KUTTA_STEP) == 4) {
258 if (r_process_info[OSS_SWITCH] == 1) {
259 CalculateOSSNodalProjections();
306 virtual void CalculateOSSNodalProjections()
312 const auto& r_lumped_mass_vector = p_explicit_bs->GetLumpedMassMatrixVector();
316 auto& r_process_info = r_model_part.GetProcessInfo();
318 ConvectionDiffusionSettings::Pointer p_settings = r_process_info[CONVECTION_DIFFUSION_SETTINGS];
319 auto& r_settings = *p_settings;
323 rNode.GetValue(r_settings.GetProjectionVariable()) = 0.0;
329 rElement.Calculate(r_settings.GetProjectionVariable(), unknown_proj, r_process_info);
333 auto it_node = r_model_part.NodesBegin() + i_node;
334 const double mass = r_lumped_mass_vector(i_node);
335 it_node->FastGetSolutionStepValue(r_settings.GetProjectionVariable()) = it_node->GetValue(r_settings.GetProjectionVariable()) / mass;
Base class for all Elements.
Definition: element.h:60
ExplicitBuilderPointerType pGetExplicitBuilder()
Operations to get the pointer to the explicit builder and solver.
Definition: explicit_solving_strategy.h:286
bool SolveSolutionStep() override
Solves the current step. The function always return true as convergence is not checked in the explici...
Definition: explicit_solving_strategy.h:230
This strategy adds the orthogonal subgrid projections computation to the base explicit runge kutta 4 ...
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:57
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:173
ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion(ModelPart &rModelPart, Parameters ThisParameters)
Default constructor. (with parameters)
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:83
TDenseSpace::VectorType LocalSystemVectorType
The local vector definition.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:69
ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion(ModelPart &rModelPart, bool MoveMeshFlag=false, int RebuildLevel=0)
Default constructor.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:110
void Initialize() override
Initialization of variables.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:140
virtual bool SolveSolutionStep() override
Finalize the Runge Kutta explicit solver step and calculate the orthogonal subscale projections if re...
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:249
ExplicitSolvingStrategyRungeKutta4< TSparseSpace, TDenseSpace > BaseType
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:63
std::string Info() const override
Turn back information as a string.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:167
BaseType::ExplicitBuilderType ExplicitBuilderType
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:66
virtual void InitializeRungeKuttaIntermediateSubStep() override
Initialize the Runge-Kutta substep.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:208
ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion(ModelPart &rModelPart, typename ExplicitBuilderType::Pointer pExplicitBuilder, bool MoveMeshFlag=false, int RebuildLevel=0)
Default constructor.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:96
virtual void InitializeRungeKuttaLastSubStep() override
Initialize the Runge-Kutta substep.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:229
virtual ~ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion()=default
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:179
ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion(const ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion &Other)=delete
KRATOS_CLASS_POINTER_DEFINITION(ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion)
Family of explicit Runge-Kutta schemes.
Definition: explicit_solving_strategy_runge_kutta.h:66
virtual void InitializeRungeKuttaLastSubStep()
Initialize the Runge-Kutta last substep This method is intended to implement all the operations requi...
Definition: explicit_solving_strategy_runge_kutta.h:333
typename BaseType::ExplicitBuilderType ExplicitBuilderType
The explicit builder and solver definition.
Definition: explicit_solving_strategy_runge_kutta.h:81
void Initialize() override
Initialization of member variables and prior operations.
Definition: explicit_solving_strategy_runge_kutta.h:195
virtual void InitializeRungeKuttaIntermediateSubStep()
Initialize the Runge-Kutta intermediate substep This method is intended to implement all the operatio...
Definition: explicit_solving_strategy_runge_kutta.h:321
This class is useful for index iteration over containers.
Definition: parallel_utilities.h:451
void for_each(TUnaryFunction &&f)
Definition: parallel_utilities.h:514
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
This class defines the node.
Definition: node.h:65
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
Solving strategy base class This is the base class from which we will derive all the strategies (impl...
Definition: solving_strategy.h:64
ModelPart & GetModelPart()
Operations to get the pointer to the model.
Definition: solving_strategy.h:350
bool MoveMeshFlag()
This function returns the flag that says if the mesh is moved.
Definition: solving_strategy.h:290
#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
void block_for_each(TIterator itBegin, TIterator itEnd, TFunction &&rFunction)
Execute a functor on all items of a range in parallel.
Definition: parallel_utilities.h:299