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.
find_nodal_neighbours_surface_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 author: Alex Jarauta
11 // Co-author : Elaf Mahrous
12 
13 
14 #if !defined(KRATOS_FIND_NODAL_NEIGHBOURS_SURFACE_PROCESS_H_INCLUDED )
15 #define KRATOS_FIND_NODAL_NEIGHBOURS_SURFACE_PROCESS_H_INCLUDED
16 
17 
18 
19 // System includes
20 #include <string>
21 #include <iostream>
22 
23 
24 // External includes
25 
26 
27 // Project includes
28 #include <pybind11/pybind11.h>
29 #include "includes/define.h"
30 #include "includes/define_python.h"
31 
32 #include "processes/process.h"
33 #include "includes/node.h"
34 // #include "includes/element.h"
35 #include "includes/condition.h"
36 #include "includes/model_part.h"
37 
38 
39 namespace Kratos
40 {
41 
44 
49 // typedef ModelPart::ElementsContainerType ElementsContainerType;
51 
52 
56 
60 
64 
66 
69  : public Process
70 {
71 public:
74 
77 
81 
86  FindNodalNeighboursSurfaceProcess(ModelPart& model_part, const int avg_conds = 8, const int avg_nodes = 8)
87  : mr_model_part(model_part)
88  {
89  mavg_conds = avg_conds;
90  mavg_nodes = avg_nodes;
91  }
92 
95  {
96  }
97 
98 
102 
103  void operator()()
104  {
105  Execute();
106  }
107 
108 
112 
113  void Execute() override
114  {
115  NodesContainerType& rNodes = mr_model_part.Nodes();
116  ConditionsContainerType& rConds = mr_model_part.Conditions();
117 
118  //first of all the neighbour nodes and conditions array are initialized to the guessed size
119  //and empties the old entries
120  for(NodesContainerType::iterator in = rNodes.begin(); in!=rNodes.end(); ++in)
121  {
122  (in->GetValue(NEIGHBOUR_NODES)).reserve(mavg_nodes);
123  GlobalPointersVector<Node >& rN = in->GetValue(NEIGHBOUR_NODES);
124  rN.erase(rN.begin(),rN.end() );
125 
126  (in->GetValue(NEIGHBOUR_CONDITIONS)).reserve(mavg_conds);
127  GlobalPointersVector<Condition >& rC = in->GetValue(NEIGHBOUR_CONDITIONS);
128  rC.erase(rC.begin(),rC.end() );
129  }
130 
131  //add the neighbour conditions to all the nodes in the mesh
132  for(ConditionsContainerType::iterator ic = rConds.begin(); ic!=rConds.end(); ++ic)
133  {
134  Condition::GeometryType& pGeom = ic->GetGeometry();
135  for(unsigned int i = 0; i < pGeom.size(); i++)
136  {
137  //KRATOS_WATCH( pGeom[i] );
138  (pGeom[i].GetValue(NEIGHBOUR_CONDITIONS)).push_back( Condition::WeakPointer( *(ic.base()) ) );
139  //KRATOS_WATCH( (pGeom[i].GetValue(NEIGHBOUR_ELEMENTS)).size() );
140  }
141  }
142 
143  //adding the neighbouring nodes
144  for(NodesContainerType::iterator in = rNodes.begin(); in!=rNodes.end(); ++in)
145  {
146  GlobalPointersVector<Condition >& rC = in->GetValue(NEIGHBOUR_CONDITIONS);
147 
148  for(unsigned int ic = 0; ic < rC.size(); ic++)
149  {
150  Condition::GeometryType& pGeom = rC[ic].GetGeometry();
151  for(unsigned int i = 0; i < pGeom.size(); i++)
152  {
153  if(pGeom[i].Id() != in->Id() )
154  {
155 // Element::NodeType::WeakPointer temp = pGeom(i);
156  //NOT SURE ABOUT THIS!!!!!!!!!
157  Condition::NodeType::WeakPointer temp = pGeom(i);
158  AddUniqueWeakPointer< Node >(in->GetValue(NEIGHBOUR_NODES), temp);
159  }
160  }
161  }
162  }
163  }
164 
166  {
167  NodesContainerType& rNodes = mr_model_part.Nodes();
168  for(NodesContainerType::iterator in = rNodes.begin(); in!=rNodes.end(); ++in)
169  {
170  GlobalPointersVector<Condition >& rC = in->GetValue(NEIGHBOUR_CONDITIONS);
171  rC.erase(rC.begin(),rC.end());
172 
173  GlobalPointersVector<Node >& rN = in->GetValue(NEIGHBOUR_NODES);
174  rN.erase(rN.begin(),rN.end() );
175  }
176  }
177 
181 
182 
186 
187 
191 
193  std::string Info() const override
194  {
195  return "FindNodalNeighboursSurfaceProcess";
196  }
197 
199  void PrintInfo(std::ostream& rOStream) const override
200  {
201  rOStream << "FindNodalNeighboursSurfaceProcess";
202  }
203 
205  void PrintData(std::ostream& rOStream) const override
206  {
207  }
208 
209 
213 
214 
216 
217 protected:
220 
221 
225 
226 
230 
231 
235 
236 
240 
241 
245 
246 
250 
251 
253 
254 private:
257 
258 
262  ModelPart& mr_model_part;
263  int mavg_conds;
264  int mavg_nodes;
265 
266 
270 
271  //******************************************************************************************
272  //******************************************************************************************
273  template< class TDataType > void AddUniqueWeakPointer
274  (GlobalPointersVector< TDataType >& v, const typename TDataType::WeakPointer candidate)
275  {
277  typename GlobalPointersVector< TDataType >::iterator endit = v.end();
278  while ( i != endit && (i)->Id() != (candidate.lock())->Id())
279  {
280  i++;
281  }
282  if( i == endit )
283  {
284  v.push_back(candidate);
285  }
286 
287  }
288 
292 
293 
297 
298 
302 
303 
307 
310 
312  //FindNodalNeighboursSurfaceProcess(FindNodalNeighboursSurfaceProcess const& rOther);
313 
314 
316 
317 }; // Class FindNodalNeighboursSurfaceProcess
318 
320 
323 
324 
328 
329 
331 inline std::istream& operator >> (std::istream& rIStream,
333 
335 inline std::ostream& operator << (std::ostream& rOStream,
337 {
338  rThis.PrintInfo(rOStream);
339  rOStream << std::endl;
340  rThis.PrintData(rOStream);
341 
342  return rOStream;
343 }
345 
346 
347 } // namespace Kratos.
348 
349 #endif // KRATOS_FIND_NODAL_NEIGHBOURS_SURFACE_PROCESS_H_INCLUDED defined
350 
351 
Short class definition.
Definition: find_nodal_neighbours_surface_process.h:70
FindNodalNeighboursSurfaceProcess(ModelPart &model_part, const int avg_conds=8, const int avg_nodes=8)
Definition: find_nodal_neighbours_surface_process.h:86
std::string Info() const override
Turn back information as a string.
Definition: find_nodal_neighbours_surface_process.h:193
~FindNodalNeighboursSurfaceProcess() override
Destructor.
Definition: find_nodal_neighbours_surface_process.h:94
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: find_nodal_neighbours_surface_process.h:199
KRATOS_CLASS_POINTER_DEFINITION(FindNodalNeighboursSurfaceProcess)
Pointer definition of FindNodalNeighboursSurfaceProcess.
void Execute() override
Execute method is used to execute the Process algorithms.
Definition: find_nodal_neighbours_surface_process.h:113
void operator()()
Definition: find_nodal_neighbours_surface_process.h:103
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: find_nodal_neighbours_surface_process.h:205
void ClearNeighbours()
Definition: find_nodal_neighbours_surface_process.h:165
Geometry base class.
Definition: geometry.h:71
SizeType size() const
Definition: geometry.h:518
IndexType const & Id() const
Id of this Geometry.
Definition: geometry.h:964
TVariableType::Type & GetValue(const TVariableType &rThisVariable)
Definition: geometry.h:627
This class is a vector which stores global pointers.
Definition: global_pointers_vector.h:61
boost::indirect_iterator< typename TContainerType::iterator > iterator
Definition: global_pointers_vector.h:79
iterator erase(iterator pos)
Definition: global_pointers_vector.h:351
iterator begin()
Definition: global_pointers_vector.h:221
size_type size() const
Definition: global_pointers_vector.h:307
iterator end()
Definition: global_pointers_vector.h:229
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
ConditionsContainerType & Conditions(IndexType ThisIndex=0)
Definition: model_part.h:1381
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
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
ModelPart::NodesContainerType NodesContainerType
Definition: find_conditions_neighbours_process.h:44
ModelPart::ConditionsContainerType ConditionsContainerType
Definition: find_conditions_neighbours_process.h:45
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
v
Definition: generate_convection_diffusion_explicit_element.py:114
float temp
Definition: rotating_cone.py:85
integer i
Definition: TensorModule.f:17