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_triangle_mesh_conditions.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: Pablo Agustin Becker
11 //
12 
13 #pragma once
14 
15 // NOTE: Before compute the remeshing it is necessary to compute the neighbours
16 
17 // System includes
18 
19 // External includes
20 
21 // Project includes
23 
24 namespace Kratos
25 {
40 
42  : public LocalRefineTriangleMeshGeneric<Condition, ModelPart::ConditionsContainerType>
43 {
44 public:
45 
51 
54  {
55 
56  }
57 
60  = default;
61 
65 
69 
80  ModelPart& rModelPart,
81  const compressed_matrix<int>& rCoord,
82  PointerVector< Element >& rNewElements,
83  bool InterpolateInternalVariables
84  ) override
85  {
86  }
87 
94  ModelPart& rModelPart,
95  const compressed_matrix<int>& rCoord
96  ) override
97  {
98  PointerVector< Condition > New_Conditions;
99 
101  rModelPart,
102  rModelPart.Conditions(),
103  rCoord,
104  New_Conditions,
105  false
106  );
107 
108  // Now update the elements in SubModelParts
109  if (New_Conditions.size() > 0)
110  {
111  KRATOS_WATCH(New_Conditions.size())
112  UpdateSubModelPartConditions(rModelPart.GetRootModelPart(), New_Conditions);
113  }
114  }
115 
116 
118  ModelPart& rModelPart,
119  compressed_matrix<int>& rCoord
120  ) override
121  {
122  KRATOS_TRY;
123 
124  ModelPart::ConditionsContainerType& rConditions = rModelPart.Conditions();
125  ModelPart::ConditionsContainerType::iterator it_begin = rConditions.ptr_begin();
126  ModelPart::ConditionsContainerType::iterator it_end = rConditions.ptr_end();
127 
128  this->SearchEdgeToBeRefinedGeneric<ModelPart::ConditionsContainerType::iterator>(it_begin,it_end,rCoord);
129 
130  KRATOS_CATCH("");
131  }
132 
133  void ResetFatherNodes(ModelPart &rModelPart) override
134  {
135  block_for_each(mModelPart.Nodes(), [&](Node& rNode)
136  {
137  if(!rNode.Has(FATHER_NODES)){
138  GlobalPointersVector<Node> empty_father_vector;
139  rNode.SetValue(FATHER_NODES, empty_father_vector);
140  }
141  });
142  }
143 
150  ModelPart& rModelPart,
151  PointerVector< Condition >& NewConditions
152  )
153  {
154  for (auto iSubModelPart = rModelPart.SubModelPartsBegin(); iSubModelPart != rModelPart.SubModelPartsEnd(); iSubModelPart++) {
155  unsigned int to_be_deleted = 0;
156  NewConditions.clear();
157 
158  // Create list of new elements in SubModelPart
159  // Count how many elements will be removed
160  for (auto iCond = iSubModelPart->ConditionsBegin(); iCond != iSubModelPart->ConditionsEnd(); iCond++) {
161  if( iCond->GetValue(SPLIT_ELEMENT) ) {
162  to_be_deleted++;
163  GlobalPointersVector< Condition >& rChildConditions = iCond->GetValue(NEIGHBOUR_CONDITIONS);
164 
165  for ( auto iChild = rChildConditions.ptr_begin(); iChild != rChildConditions.ptr_end(); iChild++ ) {
166  NewConditions.push_back((*iChild)->shared_from_this());
167  }
168  }
169  }
170 
171  // Add new elements to SubModelPart
172  iSubModelPart->Conditions().reserve( iSubModelPart->Conditions().size() + NewConditions.size() );
173  for (PointerVector< Condition >::iterator it_new = NewConditions.begin(); it_new != NewConditions.end(); it_new++) {
174  iSubModelPart->Conditions().push_back(*(it_new.base()));
175  }
176 
177  // Delete old elements
178  iSubModelPart->Conditions().Sort();
179  iSubModelPart->Conditions().erase(iSubModelPart->Conditions().end() - to_be_deleted, iSubModelPart->Conditions().end());
180 
181  //NEXT LEVEL
182  if (NewConditions.size() > 0) {
183  ModelPart &rSubModelPart = *iSubModelPart;
184  UpdateSubModelPartConditions(rSubModelPart,NewConditions);
185  }
186  }
187  }
188 
190  ModelPart& rModelPart,
191  compressed_matrix<int>& rCoord
192  ) override
193  {
194  KRATOS_TRY;
195 
196  NodesArrayType& pNodes = rModelPart.Nodes();
197  NodesArrayType::iterator it_begin = pNodes.ptr_begin();
198  NodesArrayType::iterator it_end = pNodes.ptr_end();
199 
200  rCoord.resize(pNodes.size(), pNodes.size(), false);
201 
202  for(auto i = it_begin; i!=it_end; i++) {
203  int index_i = mMapNodeIdToPos[i->Id()];
204  GlobalPointersVector< Node >& neighb_nodes = i->GetValue(NEIGHBOUR_CONDITION_NODES);
205 
206  std::vector<unsigned int> aux(neighb_nodes.size());
207  unsigned int active = 0;
208  for (auto inode = neighb_nodes.begin(); inode != neighb_nodes.end(); inode++) {
209  int index_j = mMapNodeIdToPos[inode->Id()];
210  if (index_j > index_i) {
211  aux[active] = index_j;
212  active++;
213  }
214  }
215 
216  std::sort(aux.begin(), aux.begin() + active);
217  for (unsigned int k = 0; k < active; k++) {
218  rCoord.push_back(index_i, aux[k], -1);
219  }
220  }
221 
222  KRATOS_CATCH("");
223  }
224 };
225 
226 } // namespace Kratos.
This class is a vector which stores global pointers.
Definition: global_pointers_vector.h:61
ptr_iterator ptr_begin()
Definition: global_pointers_vector.h:253
iterator begin()
Definition: global_pointers_vector.h:221
size_type size() const
Definition: global_pointers_vector.h:307
ptr_iterator ptr_end()
Definition: global_pointers_vector.h:261
iterator end()
Definition: global_pointers_vector.h:229
ModelPart::NodesContainerType NodesArrayType
Definition: local_refine_geometry_mesh.hpp:54
ModelPart & mModelPart
Definition: local_refine_geometry_mesh.hpp:238
Definition: local_refine_triangle_mesh_conditions.hpp:43
void CSRRowMatrix(ModelPart &rModelPart, compressed_matrix< int > &rCoord) override
Definition: local_refine_triangle_mesh_conditions.hpp:189
LocalRefineTriangleMeshConditions(ModelPart &model_part)
Default constructors.
Definition: local_refine_triangle_mesh_conditions.hpp:53
void EraseOldElementAndCreateNewElement(ModelPart &rModelPart, const compressed_matrix< int > &rCoord, PointerVector< Element > &rNewElements, bool InterpolateInternalVariables) override
Definition: local_refine_triangle_mesh_conditions.hpp:79
void SearchEdgeToBeRefined(ModelPart &rModelPart, compressed_matrix< int > &rCoord) override
Definition: local_refine_triangle_mesh_conditions.hpp:117
void ResetFatherNodes(ModelPart &rModelPart) override
Definition: local_refine_triangle_mesh_conditions.hpp:133
void UpdateSubModelPartConditions(ModelPart &rModelPart, PointerVector< Condition > &NewConditions)
Definition: local_refine_triangle_mesh_conditions.hpp:149
~LocalRefineTriangleMeshConditions() override=default
Destructor.
void EraseOldConditionsAndCreateNew(ModelPart &rModelPart, const compressed_matrix< int > &rCoord) override
Definition: local_refine_triangle_mesh_conditions.hpp:93
Definition: local_refine_triangle_mesh_generic.hpp:47
void EraseOldObjetcsAndCreateNewObjects(ModelPart &rThisModelPart, ModelPart::ConditionsContainerType &rObjects, const compressed_matrix< int > &rCoord, PointerVector< Condition > &rNewObjects, bool InterpolateInternalVariables)
Definition: local_refine_triangle_mesh_generic.hpp:81
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
SubModelPartIterator SubModelPartsEnd()
Definition: model_part.h:1708
SubModelPartIterator SubModelPartsBegin()
Definition: model_part.h:1698
ConditionsContainerType & Conditions(IndexType ThisIndex=0)
Definition: model_part.h:1381
ModelPart & GetRootModelPart()
Definition: model_part.cpp:510
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
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
void clear()
Definition: pointer_vector.h:309
size_type size() const
Definition: pointer_vector.h:255
iterator end()
Definition: pointer_vector.h:177
boost::indirect_iterator< typename TContainerType::iterator > iterator
Definition: pointer_vector.h:89
iterator begin()
Definition: pointer_vector.h:169
void push_back(const TPointerType &x)
Definition: pointer_vector.h:270
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_WATCH(variable)
Definition: define.h:806
#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
model_part
Definition: face_heat.py:14
active
Definition: generate_frictionless_components_mortar_condition.py:175
int k
Definition: quadrature.py:595
integer i
Definition: TensorModule.f:17