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.
eliminate_isolated_nodes_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: Pooyan Dadvand
11 // Riccardo Rossi
12 //
13 
14 
15 
16 #if !defined(KRATOS_ELIMINATE_ISOLATED_NODES_PROCESS_INCLUDED )
17 #define KRATOS_ELIMINATE_ISOLATED_NODES_PROCESS_INCLUDED
18 
19 
20 
21 // System includes
22 #include <string>
23 #include <iostream>
24 #include <algorithm>
25 
26 // External includes
27 
28 
29 // Project includes
30 #include "includes/define.h"
31 #include "processes/process.h"
32 #include "includes/node.h"
33 #include "includes/element.h"
34 #include "includes/model_part.h"
35 
36 
37 namespace Kratos
38 {
39 
42 
46 
47 
51 
55 
59 
61 //erases the nodes marked as
66  : public Process
67 {
68 public:
71 
74 
78 
81  : mr_model_part(model_part)
82  {
84  KRATOS_CATCH("")
85  }
86 
89  {
90  }
91 
92 
96 
97  void operator()()
98  {
99  Execute();
100  }
101 
102 
106 
107  void Execute() override
108  {
109  KRATOS_TRY;
110 
111  ModelPart::NodesContainerType temp_nodes_container;
112  temp_nodes_container.reserve(mr_model_part.Nodes().size());
113 
114  temp_nodes_container.swap(mr_model_part.Nodes());
115 
116  for(ModelPart::NodesContainerType::iterator i_node = temp_nodes_container.begin() ; i_node != temp_nodes_container.end() ; i_node++)
117  {
118  if( i_node->GetValue(NEIGHBOUR_NODES).size() != 0 )
119  (mr_model_part.Nodes()).push_back(*(i_node.base()));
120  }
121 
122  KRATOS_CATCH("")
123  }
124 
125 
129 
130 
134 
135 
139 
141  std::string Info() const override
142  {
143  return "EliminateIsolatedNodesProcess";
144  }
145 
147  void PrintInfo(std::ostream& rOStream) const override
148  {
149  rOStream << "EliminateIsolatedNodesProcess";
150  }
151 
153  void PrintData(std::ostream& rOStream) const override
154  {
155  }
156 
157 
161 
162 
164 
165 protected:
168 
169 
173 
174 
178 
179 
183 
184 
188 
189 
193 
194 
198 
199 
201 
202 private:
205 
206 
210  ModelPart& mr_model_part;
211  PointerVector<Node > mTrashedNodes;
212 
213 
217 
221 
222 
226 
227 
231 
232 
236 
239 
241  //EliminateIsolatedNodesProcess(EliminateIsolatedNodesProcess const& rOther);
242 
243 
245 
246 }; // Class EliminateIsolatedNodesProcess
247 
249 
252 
253 
257 
258 
260 inline std::istream& operator >> (std::istream& rIStream,
262 
264 inline std::ostream& operator << (std::ostream& rOStream,
265  const EliminateIsolatedNodesProcess& rThis)
266 {
267  rThis.PrintInfo(rOStream);
268  rOStream << std::endl;
269  rThis.PrintData(rOStream);
270 
271  return rOStream;
272 }
274 
275 
276 } // namespace Kratos.
277 
278 #endif // KRATOS_ELIMINATE_ISOLATED_NODES_PROCESS_INCLUDED defined
279 
280 
Short class definition.
Definition: eliminate_isolated_nodes_process.h:67
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: eliminate_isolated_nodes_process.h:147
void operator()()
Definition: eliminate_isolated_nodes_process.h:97
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: eliminate_isolated_nodes_process.h:153
KRATOS_CLASS_POINTER_DEFINITION(EliminateIsolatedNodesProcess)
Pointer definition of EliminateIsolatedNodesProcess.
~EliminateIsolatedNodesProcess() override
Destructor.
Definition: eliminate_isolated_nodes_process.h:88
void Execute() override
Execute method is used to execute the Process algorithms.
Definition: eliminate_isolated_nodes_process.h:107
std::string Info() const override
Turn back information as a string.
Definition: eliminate_isolated_nodes_process.h:141
EliminateIsolatedNodesProcess(ModelPart &model_part)
Default constructor.
Definition: eliminate_isolated_nodes_process.h:80
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.h:72
boost::indirect_iterator< typename TContainerType::iterator > iterator
Definition: pointer_vector_set.h:95
void swap(PointerVectorSet &rOther)
Swaps the contents of this PointerVectorSet with another.
Definition: pointer_vector_set.h:532
iterator begin()
Returns an iterator pointing to the beginning of the container.
Definition: pointer_vector_set.h:278
void reserve(int reservedsize)
Reserves memory for a specified number of elements.
Definition: pointer_vector_set.h:733
iterator end()
Returns an iterator pointing to the end of the container.
Definition: pointer_vector_set.h:314
The base class for all processes in Kratos.
Definition: process.h:49
#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
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
model_part
Definition: face_heat.py:14