17 #ifndef KRATOS_FRACTIONAL_STEP_SETTINGS_FOR_CHIMERA_H
18 #define KRATOS_FRACTIONAL_STEP_SETTINGS_FOR_CHIMERA_H
43 #include "custom_utilities/solver_settings.h"
73 template<
class TSparseSpace,
101 const std::size_t ThisDomainSize,
102 const std::size_t ThisTimeOrder,
104 const bool MoveMeshFlag,
105 const bool reform_dof_set):
106 BaseType(r_model_part,ThisDomainSize,ThisTimeOrder,use_slip,MoveMeshFlag,reform_dof_set)
137 typename TLinearSolver::Pointer pLinearSolver,
138 const double Tolerance,
139 const unsigned int MaxIter)
override
148 bool calculate_reactions =
false;
149 bool calculate_norm_dx_flag =
true;
161 std::string fs_vel_mp_name = r_model_part.
Name()+
"fs_velocity_model_part";
166 p_build_and_solver->SetEchoLevel(strategy_echo_level);
168 SchemePointerType p_scheme;
173 SchemePointerType p_temp = Kratos::make_shared<ResidualBasedIncrementalUpdateStaticSchemeSlip< TSparseSpace, TDenseSpace>> (
domain_size,
domain_size);
174 p_scheme.swap(p_temp);
178 SchemePointerType p_temp = Kratos::make_shared<ResidualBasedIncrementalUpdateStaticScheme<TSparseSpace, TDenseSpace >>();
179 p_scheme.swap(p_temp);
183 BaseType::mStrategies[rStrategyLabel] = Kratos::make_shared< ResidualBasedLinearStrategy<TSparseSpace, TDenseSpace, TLinearSolver >>
184 (r_fs_velocity_model_part,
189 calculate_norm_dx_flag);
193 std::string fs_pressure_mp_name = r_model_part.
Name()+
"fs_pressure_model_part";
198 SchemePointerType p_scheme = Kratos::make_shared<ResidualBasedIncrementalUpdateStaticScheme< TSparseSpace, TDenseSpace >> ();
199 p_build_and_solver->SetEchoLevel(strategy_echo_level);
202 BaseType::mStrategies[rStrategyLabel] = Kratos::make_shared<ResidualBasedLinearStrategy<TSparseSpace, TDenseSpace, TLinearSolver >>(
203 r_fs_pressure_model_part,
208 calculate_norm_dx_flag);
212 KRATOS_ERROR <<
"Error in FractionalStepSettingsForChimera: Unknown strategy label."<<std::endl;
229 if(pThisStrategy !=
nullptr)
243 std::string
Info()
const override
245 std::stringstream buffer;
246 buffer <<
"FractionalStepSettingsForChimera" ;
251 void PrintInfo(std::ostream& rOStream)
const override {rOStream <<
"FractionalStepSettingsForChimera";}
254 void PrintData(std::ostream& rOStream)
const override {}
279 template<
class TDenseSpace,
class TSparseSpace,
class TLinearSolver >
287 template<
class TDenseSpace,
class TSparseSpace,
class TLinearSolver >
292 rOStream << std::endl;
The base class for assigning a value to scalar variables or array_1d components processes in Kratos.
Definition: fast_transfer_between_model_parts_process.h:56
void Execute() override
Execute method is used to execute the FastTransferBetweenModelPartsProcess algorithms.
Definition: fast_transfer_between_model_parts_process.cpp:60
Helper class to define solution strategies for FS_Strategy.
Definition: fractional_step_settings_for_chimera.h:78
void SetStrategy(StrategyLabel const &rStrategyLabel, typename TLinearSolver::Pointer pLinearSolver, const double Tolerance, const unsigned int MaxIter) override
Definition: fractional_step_settings_for_chimera.h:136
std::string Info() const override
Turn back information as a string.
Definition: fractional_step_settings_for_chimera.h:243
FractionalStepSettingsForChimera(ModelPart &r_model_part, const std::size_t ThisDomainSize, const std::size_t ThisTimeOrder, const bool use_slip, const bool MoveMeshFlag, const bool reform_dof_set)
Constructor.
Definition: fractional_step_settings_for_chimera.h:100
BaseType::StrategyType StrategyType
Definition: fractional_step_settings_for_chimera.h:87
BaseType::StrategyPointerType StrategyPointerType
Definition: fractional_step_settings_for_chimera.h:88
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: fractional_step_settings_for_chimera.h:254
bool FindStrategy(StrategyLabel const &rStrategyLabel, StrategyPointerType &pThisStrategy) override
Definition: fractional_step_settings_for_chimera.h:225
BaseType::ProcessPointerType ProcessPointerType
Definition: fractional_step_settings_for_chimera.h:89
ResidualBasedBlockBuilderAndSolver< TSparseSpace, TDenseSpace, TLinearSolver > ResidualBasedBlockBuilderAndSolverType
Definition: fractional_step_settings_for_chimera.h:91
FractionalStepSettingsForChimera()=delete
Default constructor.
ResidualBasedBlockBuilderAndSolverType::Pointer ResidualBasedBlockBuilderAndSolverPointerType
Definition: fractional_step_settings_for_chimera.h:92
FractionalStepSettingsForChimera & operator=(FractionalStepSettingsForChimera const &rOther)=delete
Assignment operator.
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: fractional_step_settings_for_chimera.h:251
FractionalStepSettingsForChimera(FractionalStepSettingsForChimera const &rOther)=delete
Copy constructor.
BaseType::StrategyLabel StrategyLabel
Definition: fractional_step_settings_for_chimera.h:90
~FractionalStepSettingsForChimera()=default
Destructor.
KRATOS_CLASS_POINTER_DEFINITION(FractionalStepSettingsForChimera)
Pointer definition of FractionalStepSettingsForChimera.
SolverSettings< TSparseSpace, TDenseSpace, TLinearSolver > BaseType
Definition: fractional_step_settings_for_chimera.h:86
Implicit solving strategy base class This is the base class from which we will derive all the implici...
Definition: implicit_solving_strategy.h:61
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
std::string & Name()
Definition: model_part.h:1811
ModelPart & CreateSubModelPart(std::string const &NewSubModelPartName)
Definition: model_part.cpp:2000
Current class provides an implementation for standard builder and solving operations.
Definition: residualbased_block_builder_and_solver.h:82
Current class provides an implementation for applying the chimera constraints that is enforcing conti...
Definition: residualbased_block_builder_and_solver_with_constraints_for_chimera.h:73
This class provides the implementation of the basic tasks that are needed by the solution strategy.
Definition: scheme.h:56
Helper class to define solution strategies for FS_Strategy.
Definition: solver_settings.h:53
StrategyLabel
Definition: solver_settings.h:66
@ Velocity
Definition: solver_settings.h:66
@ Pressure
Definition: solver_settings.h:66
virtual unsigned int GetDomainSize() const
Definition: solver_settings.h:146
ModelPart & GetModelPart()
Definition: solver_settings.h:285
bool GetReformDofSet()
Definition: solver_settings.h:226
StrategyType::Pointer StrategyPointerType
Definition: solver_settings.h:62
virtual bool UseSlipConditions() const
Definition: solver_settings.h:156
std::map< StrategyLabel, double > mTolerances
Definition: solver_settings.h:292
std::map< StrategyLabel, StrategyPointerType > mStrategies
Definition: solver_settings.h:290
std::map< StrategyLabel, unsigned int > mMaxIter
Definition: solver_settings.h:294
Process::Pointer ProcessPointerType
Definition: solver_settings.h:63
virtual unsigned int GetEchoLevel()
Definition: solver_settings.h:221
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_ERROR
Definition: exception.h:161
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
int domain_size
Definition: face_heat.py:4
echo_level
Definition: script.py:68