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.
laplacian_meshmoving_strategy.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license:
9 // kratos/license.txt
10 //
11 // Main authors: Andreas Winterstein (a.winterstein@tum.de)
12 //
13 
14 #if !defined(KRATOS_NEW_LAPLACIAN_MESHMOVING_STRATEGY)
15 #define KRATOS_NEW_LAPLACIAN_MESHMOVING_STRATEGY
16 
17 /* System includes */
18 
19 /* External includes */
20 
21 /* Project includes */
30 
31 
32 namespace Kratos {
33 
55 
57 template <class TSparseSpace, class TDenseSpace, class TLinearSolver>
59  : public ImplicitSolvingStrategy<TSparseSpace, TDenseSpace, TLinearSolver> {
60 
61 public:
67 
71 
78  typename TLinearSolver::Pointer pNewLinearSolver,
79  int TimeOrder = 1,
80  bool ReformDofSetAtEachStep = false,
81  bool ComputeReactions = false,
82  bool CalculateMeshVelocities = true,
83  int EchoLevel = 0,
84  const bool ReInitializeModelPartEachStep = false)
85  : ImplicitSolvingStrategy<TSparseSpace, TDenseSpace, TLinearSolver>(ModelPart) {
86 
87  KRATOS_TRY;
88 
89  mreform_dof_set_at_each_step = ReformDofSetAtEachStep || ReInitializeModelPartEachStep;
90  mecho_level = EchoLevel;
91  mcompute_reactions = ComputeReactions;
92  mcalculate_mesh_velocities = CalculateMeshVelocities;
93  mtime_order = TimeOrder;
94  bool calculate_norm_dx_flag = false;
95  mreinitialize_model_part_at_each_step = ReInitializeModelPartEachStep;
96 
97  typename SchemeType::Pointer pscheme = typename SchemeType::Pointer(
99  TDenseSpace>());
100 
101  const std::string element_type = "LaplacianMeshMovingElement";
102  mpmesh_model_part = MoveMeshUtilities::GenerateMeshPart(
104 
105  typedef Variable<double> VarComponent;
106 
107  mpbuilder_and_solver_x = typename TBuilderAndSolverType::Pointer(
109  TSparseSpace, TDenseSpace, TLinearSolver, VarComponent>(
110  pNewLinearSolver, MESH_DISPLACEMENT_X));
111 
112  mpbuilder_and_solver_y = typename TBuilderAndSolverType::Pointer(
114  TSparseSpace, TDenseSpace, TLinearSolver, VarComponent>(
115  pNewLinearSolver, MESH_DISPLACEMENT_Y));
116 
117  mpbuilder_and_solver_z = typename TBuilderAndSolverType::Pointer(
119  TSparseSpace, TDenseSpace, TLinearSolver, VarComponent>(
120  pNewLinearSolver, MESH_DISPLACEMENT_Z));
121 
122  mstrategy_x = typename BaseType::Pointer(new ResidualBasedLinearStrategy<TSparseSpace, TDenseSpace,TLinearSolver>(
123  *mpmesh_model_part,
124  pscheme,
125  mpbuilder_and_solver_x,
126  ComputeReactions,
127  mreform_dof_set_at_each_step,
128  calculate_norm_dx_flag));
129 
130  mstrategy_x->SetEchoLevel(mecho_level);
131 
132  mstrategy_y = typename BaseType::Pointer(new ResidualBasedLinearStrategy<TSparseSpace, TDenseSpace,TLinearSolver>(
133  *mpmesh_model_part,
134  pscheme,
135  mpbuilder_and_solver_y,
136  ComputeReactions,
137  mreform_dof_set_at_each_step,
138  calculate_norm_dx_flag));
139 
140  mstrategy_y->SetEchoLevel(mecho_level);
141 
142  mstrategy_z = typename BaseType::Pointer(new ResidualBasedLinearStrategy<TSparseSpace, TDenseSpace,TLinearSolver>(
143  *mpmesh_model_part,
144  pscheme,
145  mpbuilder_and_solver_z,
146  ComputeReactions,
147  mreform_dof_set_at_each_step,
148  calculate_norm_dx_flag));
149 
150  mstrategy_z->SetEchoLevel(mecho_level);
151 
152  KRATOS_CATCH("")
153  }
154 
156 
157  };
158 
159  void Initialize() override {
160  FindNodalNeighboursProcess find_nodal_neighbours_process(*mpmesh_model_part);
161  find_nodal_neighbours_process.Execute();
162  }
163 
164  double Solve() override {
165  KRATOS_TRY;
166 
167  if (mreinitialize_model_part_at_each_step) {
169  *mpmesh_model_part, BaseType::GetModelPart(),
170  mpmesh_model_part->pGetProperties(0), "LaplacianMeshMovingElement");
171  }
172 
173  ProcessInfo &rCurrentProcessInfo = (mpmesh_model_part)->GetProcessInfo();
174 
176  mpmesh_model_part->GetCommunicator().LocalMesh().Nodes());
177 
178  unsigned int dimension =
179  BaseType::GetModelPart().GetProcessInfo()[DOMAIN_SIZE];
180 
181  if (dimension == 2) {
182  // X DIRECTION
183  rCurrentProcessInfo[LAPLACIAN_DIRECTION] = 1;
184  mstrategy_x->Solve();
185  // Y DIRECTION
186  rCurrentProcessInfo[LAPLACIAN_DIRECTION] = 2;
187  mstrategy_y->Solve();
188  } else {
189  // X DIRECTION
190  rCurrentProcessInfo[LAPLACIAN_DIRECTION] = 1;
191  mstrategy_x->Solve();
192  // Y DIRECTION
193  rCurrentProcessInfo[LAPLACIAN_DIRECTION] = 2;
194  mstrategy_y->Solve();
195  // Z DIRECTION
196  rCurrentProcessInfo[LAPLACIAN_DIRECTION] = 3;
197  mstrategy_z->Solve();
198  }
199 
200  if (mreform_dof_set_at_each_step == true)
201  {
202  mstrategy_x->Clear();
203  mstrategy_y->Clear();
204  mstrategy_z->Clear();
205  }
206 
207  return 0.0;
208 
209  KRATOS_CATCH("");
210  }
211 
235 protected:
265 private:
272  ModelPart* mpmesh_model_part;
273 
274  typename BaseType::Pointer mstrategy_x;
275  typename BaseType::Pointer mstrategy_y;
276  typename BaseType::Pointer mstrategy_z;
277  typename TBuilderAndSolverType::Pointer mpbuilder_and_solver_x;
278  typename TBuilderAndSolverType::Pointer mpbuilder_and_solver_y;
279  typename TBuilderAndSolverType::Pointer mpbuilder_and_solver_z;
280 
281  bool mreform_dof_set_at_each_step;
282  bool mcompute_reactions;
283  int mtime_order;
284  int mecho_level;
285  bool mcalculate_mesh_velocities;
286  bool mreinitialize_model_part_at_each_step;
287 
310 
313 }; /* Class LaplacianMeshMovingStrategy */
314 
321 }
322 /* namespace Kratos.*/
323 
324 #endif /* KRATOS_NEW_LAPLACIAN_MESHMOVING_STRATEGY defined */
MeshType & LocalMesh()
Returns the reference to the mesh storing all local entities.
Definition: communicator.cpp:245
This method allows to look for neighbours in a triangular or tetrahedral mesh.
Definition: find_nodal_neighbours_process.h:59
void Execute() override
This method esxecutes the neighbour search.
Definition: find_nodal_neighbours_process.cpp:49
Implicit solving strategy base class This is the base class from which we will derive all the implici...
Definition: implicit_solving_strategy.h:61
BuilderAndSolver< TSparseSpace, TDenseSpace, TLinearSolver > TBuilderAndSolverType
Definition: implicit_solving_strategy.h:84
Short class definition.
Definition: laplacian_meshmoving_strategy.h:59
KRATOS_CLASS_POINTER_DEFINITION(LaplacianMeshMovingStrategy)
ImplicitSolvingStrategy< TSparseSpace, TDenseSpace, TLinearSolver > BaseType
Definition: laplacian_meshmoving_strategy.h:68
double Solve() override
The problem of interest is solved.
Definition: laplacian_meshmoving_strategy.h:164
Scheme< TSparseSpace, TDenseSpace > SchemeType
Definition: laplacian_meshmoving_strategy.h:70
LaplacianMeshMovingStrategy(ModelPart &ModelPart, typename TLinearSolver::Pointer pNewLinearSolver, int TimeOrder=1, bool ReformDofSetAtEachStep=false, bool ComputeReactions=false, bool CalculateMeshVelocities=true, int EchoLevel=0, const bool ReInitializeModelPartEachStep=false)
Definition: laplacian_meshmoving_strategy.h:77
void Initialize() override
Initialization of member variables and prior operations.
Definition: laplacian_meshmoving_strategy.h:159
BaseType::TBuilderAndSolverType TBuilderAndSolverType
Definition: laplacian_meshmoving_strategy.h:69
virtual ~LaplacianMeshMovingStrategy()
Definition: laplacian_meshmoving_strategy.h:155
NodesContainerType & Nodes()
Definition: mesh.h:346
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
PropertiesType::Pointer pGetProperties(IndexType PropertiesId, IndexType MeshIndex=0)
Returns the Properties::Pointer corresponding to it's identifier.
Definition: model_part.cpp:664
Communicator & GetCommunicator()
Definition: model_part.h:1821
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
Definition: residualbased_elimination_builder_and_solver_componentwise.h:95
This class provides the implementation of a static scheme.
Definition: residualbased_incrementalupdate_static_scheme.h:57
This is a very simple strategy to solve linearly the problem.
Definition: residualbased_linear_strategy.h:64
This class provides the implementation of the basic tasks that are needed by the solution strategy.
Definition: scheme.h:56
ModelPart & GetModelPart()
Operations to get the pointer to the model.
Definition: solving_strategy.h:350
This class implements a set of auxiliar, already parallelized, methods to perform some common tasks r...
Definition: variable_utils.h:63
void UpdateCurrentToInitialConfiguration(const ModelPart::NodesContainerType &rNodes)
This method updates the current nodal coordinates back to the initial coordinates.
Definition: variable_utils.cpp:245
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
static int EchoLevel
Definition: co_sim_EMPIRE_API.h:42
void CalculateMeshVelocities(ModelPart &rModelPart, const TimeDiscretization::BDF1 &rBDF)
Definition: mesh_velocity_calculation.cpp:56
void InitializeMeshPartWithElements(ModelPart &rDestinationModelPart, ModelPart &rOriginModelPart, Properties::Pointer pProperties, const std::string &rElementName)
Definition: move_mesh_utilities.cpp:176
ModelPart * GenerateMeshPart(ModelPart &rModelPart, const std::string &rElementName)
Definition: move_mesh_utilities.cpp:143
ProcessInfo & GetProcessInfo(ModelPart &rModelPart)
Definition: add_model_part_to_python.cpp:54
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
int dimension
Definition: isotropic_damage_automatic_differentiation.py:123
element_type
Choosing element type for lagrangian_model_part.
Definition: lagrangian_droplet_test.py:56