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.
shared_points_mapper.h
Go to the documentation of this file.
1 //
2 // Project Name: Kratos
3 // Last Modified by: $Author: pooyan $
4 // Date: $Date: 2006-11-27 16:07:42 $
5 // Revision: $Revision: 1.1.1.1 $
6 //
7 //
8 
9 
10 #if !defined(KRATOS_SHARED_POINTS_MAPPER_H_INCLUDED )
11 #define KRATOS_SHARED_POINTS_MAPPER_H_INCLUDED
12 
13 
14 
15 // System includes
16 #include <set>
17 
18 
19 // External includes
20 #include "boost/smart_ptr.hpp"
21 
22 
23 // Project includes
24 #include "includes/define.h"
25 #include "includes/model_part.h"
27 #include "includes/node.h"
28 
29 
30 namespace Kratos
31 {
32 
35 
39 
43 
47 
51 
53 
56 {
57 public:
60 
63 
67 
69  //**************************************************************************************************
70  //**************************************************************************************************
72  const ModelPart::NodesContainerType& OriginNodes,
73  const ModelPart::NodesContainerType& DestinationNodes,
74  double tol = 1e-9)
75  {
77 
78  if( OriginNodes.size()!=0 && DestinationNodes.size()!=0)
79  {
80  //if( OriginNodes.size() != DestinationNodes.size() )
81  // KRATOS_THROW_ERROR(std::logic_error,"wrong number of nodes","");
82 
83  mOriginNodes.reserve( OriginNodes.size() );
84  mDestinationNodes.reserve( DestinationNodes.size() );
85 
86  for(ModelPart::NodesContainerType::const_iterator origin_it = OriginNodes.begin(); origin_it != OriginNodes.end(); origin_it++)
87  {
88  for(ModelPart::NodesContainerType::const_iterator destination_it = DestinationNodes.begin(); destination_it != DestinationNodes.end(); destination_it++)
89  {
90  if (
91  fabs(origin_it->X() - destination_it->X() ) < tol &&
92  fabs(origin_it->Y() - destination_it->Y() ) < tol &&
93  fabs(origin_it->Z() - destination_it->Z() ) < tol
94  )
95  {
96  mOriginNodes.push_back( *(origin_it.base() ) );
97  mDestinationNodes.push_back( *(destination_it.base() ) );
98  }
99  }
100 
101  }
102  }
103 
104  KRATOS_CATCH("")
105  }
106 
107 
109  virtual ~SharedPointsMapper() {}
110 
111 
115 
116 
120  //**************************************************************************************************
121  //**************************************************************************************************
122  void ScalarMap( const Variable<double>& rOriginVariable,
123  const Variable<double>& rDestinationVariable)
124  {
125  KRATOS_TRY
126 
127  PointerVector< Node >::iterator it_origin = mOriginNodes.begin();
128  PointerVector< Node >::iterator it_destination = mDestinationNodes.begin();
129 
130  for(unsigned int i = 0 ; i < mOriginNodes.size() ; ++i)
131  {
132  (it_destination++ )->FastGetSolutionStepValue(rDestinationVariable) =
133  (it_origin++ )->FastGetSolutionStepValue(rOriginVariable);
134 
135  }
136 
137  KRATOS_CATCH("")
138  }
139 
140  void InverseScalarMap( const Variable<double>& rOriginVariable,
141  const Variable<double>& rDestinationVariable)
142  {
143  KRATOS_TRY
144 
145  PointerVector< Node >::iterator it_origin = mOriginNodes.begin();
146  PointerVector< Node >::iterator it_destination = mDestinationNodes.begin();
147 
148  for(unsigned int i = 0 ; i < mOriginNodes.size() ; ++i)
149  {
150  (it_origin++ )->FastGetSolutionStepValue(rOriginVariable) =
151  (it_destination++ )->FastGetSolutionStepValue(rDestinationVariable);
152 
153  }
154 
155  KRATOS_CATCH("")
156  }
157 
158  //**************************************************************************************************
159  //**************************************************************************************************
160  void VectorMap( const Variable<array_1d<double,3> >& rOriginVariable,
161  const Variable<array_1d<double,3> >& rDestinationVariable)
162  {
163  KRATOS_TRY
164 
165  PointerVector< Node >::iterator it_origin = mOriginNodes.begin();
166  PointerVector< Node >::iterator it_destination = mDestinationNodes.begin();
167 
168  for(unsigned int i = 0 ; i < mOriginNodes.size() ; ++i)
169  {
170  (it_destination++ )->FastGetSolutionStepValue(rDestinationVariable) =
171  (it_origin++ )->FastGetSolutionStepValue(rOriginVariable);
172 
173  }
174 
175  KRATOS_CATCH("")
176  }
177 
178  //**************************************************************************************************
179  //**************************************************************************************************
180  void InverseVectorMap( const Variable<array_1d<double,3> >& rOriginVariable,
181  const Variable<array_1d<double,3> >& rDestinationVariable)
182  {
183  KRATOS_TRY
184 
185  PointerVector< Node >::iterator it_origin = mOriginNodes.begin();
186  PointerVector< Node >::iterator it_destination = mDestinationNodes.begin();
187 
188  for(unsigned int i = 0 ; i < mOriginNodes.size() ; ++i)
189  {
190  noalias((it_origin++ )->FastGetSolutionStepValue(rOriginVariable)) =
191  (it_destination++ )->FastGetSolutionStepValue(rDestinationVariable);
192 
193 
194  }
195 
196  KRATOS_CATCH("")
197  }
198 
199 
200 
204 
205 
209 
210 
214 
216  virtual void PrintInfo(std::ostream& OStream) const
217  {
218  }
219 
220 
222  virtual void PrintData(std::ostream& OStream) const
223  {
224  }
225 
226 
230 
231 
233 
234 protected:
237 
238 
243 
245 
246 
250 
251 
255 
256 
260 
261 
265 
266 
270 
271 
273 
274 private:
277 
278 
282 
283 
287 
288 
292 
293 
297 
298 
302 
303 
307 
309  SharedPointsMapper& operator=(const SharedPointsMapper& Other);
310 
312  //SharedPointsMapper(const SharedPointsMapper& Other);
313 
314 
316 
317 }; // Class SharedPointsMapper
318 
320 
323 
324 
328 
330 
331 
332 } // namespace Kratos.
333 
334 #endif // KRATOS_SHARED_POINTS_MAPPER_H_INCLUDED defined
335 
336 
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
boost::indirect_iterator< typename TContainerType::iterator > iterator
Definition: pointer_vector.h:89
Short class definition.
Definition: shared_points_mapper.h:56
void InverseScalarMap(const Variable< double > &rOriginVariable, const Variable< double > &rDestinationVariable)
Definition: shared_points_mapper.h:140
SharedPointsMapper(const ModelPart::NodesContainerType &OriginNodes, const ModelPart::NodesContainerType &DestinationNodes, double tol=1e-9)
Constructor with given array of Nodes.
Definition: shared_points_mapper.h:71
void ScalarMap(const Variable< double > &rOriginVariable, const Variable< double > &rDestinationVariable)
Definition: shared_points_mapper.h:122
virtual ~SharedPointsMapper()
Destructor.
Definition: shared_points_mapper.h:109
void InverseVectorMap(const Variable< array_1d< double, 3 > > &rOriginVariable, const Variable< array_1d< double, 3 > > &rDestinationVariable)
Definition: shared_points_mapper.h:180
virtual void PrintInfo(std::ostream &OStream) const
Print information about this object.
Definition: shared_points_mapper.h:216
void VectorMap(const Variable< array_1d< double, 3 > > &rOriginVariable, const Variable< array_1d< double, 3 > > &rDestinationVariable)
Definition: shared_points_mapper.h:160
virtual void PrintData(std::ostream &OStream) const
Print object's data.
Definition: shared_points_mapper.h:222
KRATOS_CLASS_POINTER_DEFINITION(SharedPointsMapper)
Counted pointer of SharedPointsMapper.
PointerVector< Node > mDestinationNodes
Definition: shared_points_mapper.h:244
PointerVector< Node > mOriginNodes
Definition: shared_points_mapper.h:242
#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
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
int tol
Definition: hinsberg_optimization.py:138
integer i
Definition: TensorModule.f:17
e
Definition: run_cpp_mpi_tests.py:31