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.
search_utilities.h
Go to the documentation of this file.
1 #ifndef DEM_SEARCH_UTILITIES_H
2 #define DEM_SEARCH_UTILITIES_H
3 
4 /* System includes */
5 #include <iostream>
6 
7 /* External includes */
8 #ifdef _OPENMP
9 #include <omp.h>
10 #endif
11 
12 /* Project includes */
13 #include "includes/define.h"
14 #include "includes/model_part.h"
15 #include "includes/variables.h"
16 #include "utilities/openmp_utils.h"
17 #include "utilities/timer.h"
18 
20 
21 namespace Kratos
22 {
23 
25  {
26  public:
27 
28  typedef SpatialSearch::Pointer SpatialSearchPtrType;
29 
32 
34 
38 
39  typedef std::vector<double> RadiusArrayType;
40 
42 
44 
46  {
47  mSpatialSearch = pSpatialSearch;
48  }
49 
51 
52  virtual ~DemSearchUtilities(){}
53 
54  //************************************************************************
55  // Elemental Distance Calcualtion
56  //************************************************************************
57 
58 
59 
60  //************************************************************************
61  // Nodal Distance Calcualtion
62  //************************************************************************
63 
71  template<class TVariableType>
72  void SearchNodeNeigboursDistances(ModelPart& rSearchModelPart, ModelPart& rBinsModelPart, const double& rSearchRadius, const TVariableType& rDistanceVar)
73  {
75 
76  NodesArrayType& rSearchNodes = rSearchModelPart.GetCommunicator().LocalMesh().Nodes();
77  NodesArrayType& rBinsNodes = rBinsModelPart.GetCommunicator().LocalMesh().Nodes();
78 
79  SearchNodeNeigboursDistances(rSearchNodes,rBinsNodes,rSearchRadius,rDistanceVar);
80 
81  KRATOS_CATCH("")
82  }
83 
91  template<class TVariableType>
92  void SearchNodeNeigboursDistances(ModelPart& rSearchModelPart, NodesArrayType& rBinsNodes, const double& rSearchRadius, const TVariableType& rDistanceVar)
93  {
95 
96  NodesArrayType& rSearchNodes = rSearchModelPart.GetCommunicator().LocalMesh().Nodes();
97 
98  SearchNodeNeigboursDistances(rSearchNodes,rBinsNodes,rSearchRadius,rDistanceVar);
99 
100  KRATOS_CATCH("")
101  }
102 
110  template<class TVariableType>
111  void SearchNodeNeigboursDistances(NodesArrayType& rSearchNodes, ModelPart& rBinsModelPart, const double& rSearchRadius, const TVariableType& rDistanceVar)
112  {
113  KRATOS_TRY
114 
115  NodesArrayType& rBinsNodes = rBinsModelPart.GetCommunicator().LocalMesh().Nodes();
116 
117  SearchNodeNeigboursDistances(rSearchNodes,rBinsNodes,rSearchRadius,rDistanceVar);
118 
119  KRATOS_CATCH("")
120  }
121 
130  template<class TVariableType>
131  void SearchNodeNeigboursDistances(NodesArrayType& rSearchNodes, NodesArrayType& rBinsNodes, const double& rSearchRadius, const TVariableType& rDistanceVar)
132  {
133  KRATOS_TRY
134  KRATOS_ERROR << "This function uses FastGetSolutionStepValue(RADIUS) instead of the list of radii!" << std::endl;
135 
136  if(rSearchNodes.size() && rBinsNodes.size())
137  {
138  std::size_t node_size = rSearchNodes.size();
139 
140  mResultsDistances.resize(node_size);
141  mSearchRadii.resize(node_size);
142  mNodesResults.resize(node_size);
143 
144  mResultsDistances.clear();
145  mSearchRadii.clear();
146  mNodesResults.clear();
147 
148  for(NodesArrayType::iterator it = rSearchNodes.begin(); it != rSearchNodes.end(); ++it)
149  {
150  mSearchRadii[it-rSearchNodes.begin()] = rSearchRadius;
151  }
152 
153  mSpatialSearch->SearchNodesInRadiusExclusive(rBinsNodes,rSearchNodes,mSearchRadii,mNodesResults,mResultsDistances);
154 
155  for(NodesArrayType::iterator it = rSearchNodes.begin(); it != rSearchNodes.end(); ++it)
156  {
157  int i = it - rSearchNodes.begin();
158  double minDist = 1;
159 
160  if(mResultsDistances[i].size())
161  {
162  minDist = mResultsDistances[i][0] - mNodesResults[i][0]->FastGetSolutionStepValue(RADIUS);
163 
164  for(std::size_t j = 0; j < mResultsDistances[i].size() ; j++)
165  {
166  mResultsDistances[i][j] = mResultsDistances[i][j] - mNodesResults[i][j]->FastGetSolutionStepValue(RADIUS);
167  minDist = minDist < mResultsDistances[i][j] ? minDist : mResultsDistances[i][j];
168  }
169  }
170 
171  it->FastGetSolutionStepValue(rDistanceVar) = minDist;
172  }
173  }
174 
175  KRATOS_CATCH("")
176  }
177 
178  //************************************************************************
179  // Condition Distance Calculation
180  //************************************************************************
181 
182  //TO BE IMPLEMENTED
183 
184  virtual std::string Info() const
185  {
186  return "";
187  }
188 
190 
191  virtual void PrintInfo(std::ostream& rOStream) const
192  {
193  }
194 
196 
197  virtual void PrintData(std::ostream& rOStream) const
198  {
199  }
200 
201 
202  protected:
203 
207 
210 
212 
213  std::vector<unsigned int> mPartition;
214 
215  private:
216 
218  DemSearchUtilities & operator=(DemSearchUtilities const& rOther);
219 
220 
221  }; // Class DemSearchUtilities
222 
223 
224 } // namespace Kratos.
225 
226 #endif // DEM_SEARCH_UTILITIES_H
MeshType & LocalMesh()
Returns the reference to the mesh storing all local entities.
Definition: communicator.cpp:245
Definition: search_utilities.h:25
SpatialSearchPtrType mSpatialSearch
Definition: search_utilities.h:211
virtual ~DemSearchUtilities()
Destructor.
Definition: search_utilities.h:52
virtual void PrintData(std::ostream &rOStream) const
Print object's data.
Definition: search_utilities.h:197
virtual std::string Info() const
Definition: search_utilities.h:184
virtual void PrintInfo(std::ostream &rOStream) const
Print information about this object.
Definition: search_utilities.h:191
SpatialSearch::NodesContainerType NodesArrayType
Definition: search_utilities.h:31
SpatialSearch::ElementsContainerType ElementsArrayType
Definition: search_utilities.h:30
SpatialSearch::Pointer SpatialSearchPtrType
Definition: search_utilities.h:28
VectorResultNodesContainerType mNodesResults
Definition: search_utilities.h:208
DemSearchUtilities(SpatialSearchPtrType pSpatialSearch)
Default constructor.
Definition: search_utilities.h:45
SpatialSearch::VectorResultNodesContainerType VectorResultNodesContainerType
Definition: search_utilities.h:37
SpatialSearch::VectorResultElementsContainerType VectorResultElementsContainerType
Definition: search_utilities.h:36
void SearchNodeNeigboursDistances(NodesArrayType &rSearchNodes, ModelPart &rBinsModelPart, const double &rSearchRadius, const TVariableType &rDistanceVar)
Definition: search_utilities.h:111
VectorDistanceType mResultsDistances
Definition: search_utilities.h:209
std::vector< unsigned int > mPartition
Definition: search_utilities.h:213
void SearchNodeNeigboursDistances(ModelPart &rSearchModelPart, ModelPart &rBinsModelPart, const double &rSearchRadius, const TVariableType &rDistanceVar)
Definition: search_utilities.h:72
std::vector< double > RadiusArrayType
Definition: search_utilities.h:39
RadiusArrayType mSearchRadii
Definition: search_utilities.h:206
KRATOS_CLASS_POINTER_DEFINITION(DemSearchUtilities)
void SearchNodeNeigboursDistances(NodesArrayType &rSearchNodes, NodesArrayType &rBinsNodes, const double &rSearchRadius, const TVariableType &rDistanceVar)
Definition: search_utilities.h:131
SpatialSearch::VectorDistanceType VectorDistanceType
Definition: search_utilities.h:35
SpatialSearch::NodesContainerType::ContainerType NodesContainerType
Definition: search_utilities.h:33
void SearchNodeNeigboursDistances(ModelPart &rSearchModelPart, NodesArrayType &rBinsNodes, const double &rSearchRadius, const TVariableType &rDistanceVar)
Definition: search_utilities.h:92
NodesContainerType & Nodes()
Definition: mesh.h:346
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
Communicator & GetCommunicator()
Definition: model_part.h:1821
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
size_type size() const
Returns the number of elements in the container.
Definition: pointer_vector_set.h:502
TContainerType ContainerType
Definition: pointer_vector_set.h:90
iterator begin()
Returns an iterator pointing to the beginning of the container.
Definition: pointer_vector_set.h:278
iterator end()
Returns an iterator pointing to the end of the container.
Definition: pointer_vector_set.h:314
std::vector< DistanceType > VectorDistanceType
Definition: spatial_search.h:97
std::vector< ResultNodesContainerType > VectorResultNodesContainerType
Definition: spatial_search.h:82
std::vector< ResultElementsContainerType > VectorResultElementsContainerType
Definition: spatial_search.h:87
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_ERROR
Definition: exception.h:161
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
int j
Definition: quadrature.py:648
integer i
Definition: TensorModule.f:17