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.
local_refine_geometry_mesh.hpp
Go to the documentation of this file.
1 // KRATOS __ __ _____ ____ _ _ ___ _ _ ____
2 // | \/ | ____/ ___|| | | |_ _| \ | |/ ___|
3 // | |\/| | _| \___ \| |_| || || \| | | _
4 // | | | | |___ ___) | _ || || |\ | |_| |
5 // |_| |_|_____|____/|_| |_|___|_| \_|\____| APPLICATION
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Nelson Lafontaine
11 // Jordi Cotela Dalmau
12 // Riccardo Rossi
13 // Co-authors: Vicente Mataix Ferrandiz
14 //
15 
16 #pragma once
17 
18 
19 // NOTE: Before compute the remeshing it is necessary to compute the neighbours
20 
21 // System includes
22 
23 // Extrenal includes
24 
25 // Project includes
26 #include "includes/variables.h"
28 #include "includes/model_part.h"
30 #include "geometries/geometry.h"
31 
32 namespace Kratos
33 {
48 class KRATOS_API(MESHING_APPLICATION) LocalRefineGeometryMesh
49 {
50 public:
51 
57  typedef std::vector<Matrix> Matrix_Order_Tensor;
58  typedef std::vector<Vector> Vector_Order_Tensor;
59  typedef std::vector<Vector_Order_Tensor> Node_Vector_Order_Tensor;
60  typedef Node PointType;
61  typedef Node ::Pointer PointPointerType;
62  typedef std::vector<PointType::Pointer> PointVector;
63  typedef PointVector::iterator PointIterator;
64 
68 
71  {
72 
73  }
74 
77  = default;
78 
82 
86 
93  virtual void LocalRefineMesh(
94  bool refine_on_reference,
95  bool interpolate_internal_variables
96  );
97 
104  virtual void CSRRowMatrix(
105  ModelPart& this_model_part,
106  compressed_matrix<int>& Coord
107  );
108 
115  virtual void SearchEdgeToBeRefined(
116  ModelPart& this_model_part,
117  compressed_matrix<int>& Coord
118  );
119 
128  virtual void CreateListOfNewNodes(
129  ModelPart& this_model_part,
130  compressed_matrix<int>& Coord,
131  std::vector<int> &List_New_Nodes,
132  std::vector<array_1d<int, 2 > >& Position_Node
133  );
134 
143  virtual void CalculateCoordinateAndInsertNewNodes(
144  ModelPart& this_model_part,
145  const std::vector<array_1d<int, 2 > >& Position_Node,
146  const std::vector<int> &List_New_Nodes
147  );
148 
157  virtual void EraseOldElementAndCreateNewElement(
158  ModelPart& this_model_part,
159  const compressed_matrix<int>& Coord,
160  PointerVector< Element >& New_Elements,
161  bool interpolate_internal_variables
162  );
163 
170  virtual void EraseOldConditionsAndCreateNew(
171  ModelPart& this_model_part,
172  const compressed_matrix<int>& Coord
173  );
174 
184  virtual void CalculateEdges(
185  Element::GeometryType& geom,
186  const compressed_matrix<int>& Coord,
187  int* edge_ids,
188  std::vector<int> & aux
189  );
190 
198  inline void CreatePartition(
199  unsigned int number_of_threads,
200  const int number_of_rows,
201  vector<unsigned int>& partitions
202  );
203 
212  template<typename TGeometricalObjectPointerType>
214  const int& number_elem,
215  const TGeometricalObjectPointerType father_elem,
216  TGeometricalObjectPointerType child_elem,
217  const ProcessInfo& rCurrentProcessInfo
218  )
219  {
220  // NOTE: Right now there is not an interpolation at all, it just copying the values
221  std::vector<Vector> values;
222  father_elem->CalculateOnIntegrationPoints(INTERNAL_VARIABLES, values, rCurrentProcessInfo);
223  child_elem->SetValuesOnIntegrationPoints(INTERNAL_VARIABLES, values, rCurrentProcessInfo);
224  }
225 
226  virtual void UpdateSubModelPartNodes(ModelPart &rModelPart);
227 
228  virtual void ResetFatherNodes(ModelPart &rModelPart);
229 
230 protected:
233 
237 
240  std::unordered_map<std::size_t, unsigned int> mMapNodeIdToPos;
241  std::vector<std::size_t> mMapPosToNodeId;
245 
246  template<typename TIteratorType>
248  TIteratorType GeometricalObjectsBegin,
249  TIteratorType GeometricalObjectsEnd,
250  compressed_matrix<int>& rCoord
251  )
252  {
253  KRATOS_TRY;
254 
255  for (TIteratorType it = GeometricalObjectsBegin; it != GeometricalObjectsEnd; ++it) {
256  if (it->GetValue(SPLIT_ELEMENT)) {
257  auto& r_geom = it->GetGeometry(); // Nodes of the element
258  for (unsigned int i = 0; i < r_geom.size(); i++) {
259  int index_i = mMapNodeIdToPos[r_geom[i].Id()];
260  for (unsigned int j = 0; j < r_geom.size(); j++) {
261  int index_j = mMapNodeIdToPos[r_geom[j].Id()];
262  if (index_j > index_i)
263  {
264  rCoord(index_i, index_j) = -2;
265  }
266  }
267  }
268  }
269  }
270 
271  KRATOS_CATCH("");
272  }
273 
277 
281 
285 
290 
291 private:
294 
298 
302 
306 
310 
314 
318 
320 };
321 
322 } // namespace Kratos.
Geometry base class.
Definition: geometry.h:71
Definition: local_refine_geometry_mesh.hpp:49
LocalRefineGeometryMesh(ModelPart &model_part)
Default constructors.
Definition: local_refine_geometry_mesh.hpp:70
std::vector< std::size_t > mMapPosToNodeId
Definition: local_refine_geometry_mesh.hpp:241
ModelPart::NodesContainerType NodesArrayType
Definition: local_refine_geometry_mesh.hpp:54
ModelPart & mModelPart
Definition: local_refine_geometry_mesh.hpp:238
std::unordered_map< std::size_t, unsigned int > mMapNodeIdToPos
The current refinement level.
Definition: local_refine_geometry_mesh.hpp:240
void InterpolateInteralVariables(const int &number_elem, const TGeometricalObjectPointerType father_elem, TGeometricalObjectPointerType child_elem, const ProcessInfo &rCurrentProcessInfo)
Definition: local_refine_geometry_mesh.hpp:213
ModelPart::ElementsContainerType ElementsArrayType
Definition: local_refine_geometry_mesh.hpp:55
std::vector< Matrix > Matrix_Order_Tensor
Definition: local_refine_geometry_mesh.hpp:57
virtual ~LocalRefineGeometryMesh()=default
Destructor.
Node ::Pointer PointPointerType
Definition: local_refine_geometry_mesh.hpp:61
std::vector< PointType::Pointer > PointVector
Definition: local_refine_geometry_mesh.hpp:62
std::vector< Vector_Order_Tensor > Node_Vector_Order_Tensor
Definition: local_refine_geometry_mesh.hpp:59
ModelPart::ConditionsContainerType ConditionsArrayType
Definition: local_refine_geometry_mesh.hpp:56
Node PointType
Definition: local_refine_geometry_mesh.hpp:60
std::vector< Vector > Vector_Order_Tensor
Definition: local_refine_geometry_mesh.hpp:58
PointVector::iterator PointIterator
Definition: local_refine_geometry_mesh.hpp:63
void SearchEdgeToBeRefinedGeneric(TIteratorType GeometricalObjectsBegin, TIteratorType GeometricalObjectsEnd, compressed_matrix< int > &rCoord)
Definition: local_refine_geometry_mesh.hpp:247
int mCurrentRefinementLevel
The model part to be refined.
Definition: local_refine_geometry_mesh.hpp:239
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
MeshType::ConditionsContainerType ConditionsContainerType
Condintions container. A vector set of Conditions with their Id's as key.
Definition: model_part.h:183
MeshType::ElementsContainerType ElementsContainerType
Element container. A vector set of Elements with their Id's as key.
Definition: model_part.h:168
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
This class defines the node.
Definition: node.h:65
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
#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
list values
Definition: bombardelli_test.py:42
model_part
Definition: face_heat.py:14
int j
Definition: quadrature.py:648
integer i
Definition: TensorModule.f:17