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.
parallel_distance_calculation_process.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: Riccardo Rossi
11 //
12 //
13 
14 #if !defined(KRATOS_PARALLEL_DISTANCE_CALCULATOR_H_INCLUDED )
15 #define KRATOS_PARALLEL_DISTANCE_CALCULATOR_H_INCLUDED
16 
17 // System includes
18 #include <string>
19 #include <iostream>
20 
21 // External includes
22 
23 // Project includes
24 #include "includes/define.h"
25 #include "includes/model_part.h"
26 #include "process.h"
27 
28 namespace Kratos
29 {
30 
33 
34 
38 
39 
43 
44 
48 
49 
53 
54 
56 
58 template <unsigned int TDim>
59 class KRATOS_API(KRATOS_CORE) ParallelDistanceCalculationProcess : public Process
60 {
61 public:
64 
66  enum class NodalDatabase {
67  NodeHistorical,
68  NodeNonHistorical
69  };
70 
72  using NodeType = typename ModelPart::NodeType;
73 
75  using NodeScalarGetFunctionType = std::function<double&(NodeType& rNode)>;
76 
79 
83 
85  ModelPart &rModelPart,
86  Parameters Settings);
87 
89  Model &rModel,
90  Parameters Settings);
91 
94 
97 
98  const Parameters GetDefaultParameters() const override;
99 
100  int Check() override;
101 
102  void Execute() override;
103 
104  //TODO: This method must be removed as it does not match the base process API
105  double FindMaximumEdgeSize();
106 
110 
111 
115 
116 
120 
121 
125 
126 
130 
132  std::string Info() const override
133  {
134  std::stringstream buffer;
135  buffer << "ParallelDistanceCalculationProcess" << TDim << "D";
136  return buffer.str();
137  };
138 
140  void PrintInfo(std::ostream& rOStream) const override
141  {
142  rOStream << "ParallelDistanceCalculationProcess" << TDim << "D";
143  };
144 
146  void PrintData(std::ostream& rOStream) const override {};
147 
151 
152 
154 protected:
157 
158 
162 
166  unsigned int mMaxLevels;
167  double mMaxDistance;
169 
170 
174 
175 
179 
180 
184 
185 
189 
190 
194 
195 
197 private:
200 
201 
205 
206  const Variable<double>* mpAuxDistanceVar;
207 
208  NodalDatabase mDistanceDatabase;
209 
210  NodeScalarGetFunctionType mDistanceGetter;
211  NodeScalarGetFunctionType mAuxDistanceGetter;
212  NodeScalarGetFunctionType mNodalAreaGetter;
213 
217 
218 
222 
223  void ResetVariables();
224 
225  void CalculateExactDistancesOnDividedElements();
226 
227  void ExtendDistancesByLayer();
228 
229  void AssignDistanceSign();
230 
231  bool IsDivided(const array_1d<double,TDim+1>& rDistance);
232 
233  bool IsActive(const array_1d<double,TDim+1>& rVisited);
234 
235  void AddDistanceToNodes(
236  Geometry<NodeType>& rGeometry,
237  const BoundedMatrix<double,TDim+1,TDim>& rDN_DX,
238  const double& Volume);
239 
243 
244 
248 
249 
253 
254 
256 }; // Class ParallelDistanceCalculationProcess
257 
261 
263 template<unsigned int TDim>
264 inline std::istream& operator >> (
265  std::istream& rIStream,
267 {
268  return rIStream;
269 }
270 
272 template<unsigned int TDim>
273 inline std::ostream& operator << (
274  std::ostream& rOStream,
276 {
277  rThis.PrintInfo(rOStream);
278  rOStream << std::endl;
279  rThis.PrintData(rOStream);
280 
281  return rOStream;
282 }
284 
285 } // namespace Kratos.
286 
287 #endif // KRATOS_PARALLEL_DISTANCE_CALCULATOR_H_INCLUDED defined
Geometry base class.
Definition: geometry.h:71
Definition: amatrix_interface.h:41
This class aims to manage different model parts across multi-physics simulations.
Definition: model.h:60
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
Node NodeType
Definition: model_part.h:117
Short class definition.
Definition: parallel_distance_calculation_process.h:60
double mMaxDistance
Definition: parallel_distance_calculation_process.h:167
unsigned int mMaxLevels
Definition: parallel_distance_calculation_process.h:166
KRATOS_CLASS_POINTER_DEFINITION(ParallelDistanceCalculationProcess)
Pointer definition of ParallelDistanceCalculationProcess.
bool mCalculateExactDistancesToPlane
Definition: parallel_distance_calculation_process.h:168
ParallelDistanceCalculationProcess(ParallelDistanceCalculationProcess< TDim > const &rOther)=delete
Copy constructor.
std::function< double &(NodeType &rNode)> NodeScalarGetFunctionType
Lambda type to retrieve a scalar variable from a node.
Definition: parallel_distance_calculation_process.h:75
const Variable< double > * mpAreaVar
Definition: parallel_distance_calculation_process.h:165
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: parallel_distance_calculation_process.h:140
ModelPart & mrModelPart
Definition: parallel_distance_calculation_process.h:163
virtual ~ParallelDistanceCalculationProcess()
Destructor.
Definition: parallel_distance_calculation_process.h:96
const Variable< double > * mpDistanceVar
Definition: parallel_distance_calculation_process.h:164
std::string Info() const override
Turn back information as a string.
Definition: parallel_distance_calculation_process.h:132
NodalDatabase
Nodal databases auxiliary enum.
Definition: parallel_distance_calculation_process.h:66
typename ModelPart::NodeType NodeType
Node type from model part.
Definition: parallel_distance_calculation_process.h:72
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: parallel_distance_calculation_process.h:146
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
The base class for all processes in Kratos.
Definition: process.h:49
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