meshToMesh0.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2020 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::meshToMesh0
29 
30 Description
31  Serial mesh to mesh interpolation class.
32 
33 Note
34  This class is due to be deprecated in favour of meshToMesh0New
35 
36 SourceFiles
37  meshToMesh0.C
38  calculateMeshToMesh0Addressing.C
39  calculateMeshToMesh0Weights.C
40  meshToMesh0Templates.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef meshToMesh0_H
45 #define meshToMesh0_H
46 
47 #include "fvMesh.H"
48 #include "HashTable.H"
49 #include "fvPatchMapper.H"
50 #include "scalarList.H"
51 #include "className.H"
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 namespace Foam
56 {
57 
58 // Forward Declarations
59 template<class Type> class indexedOctree;
60 class treeDataCell;
61 
62 /*---------------------------------------------------------------------------*\
63  Class meshToMesh0 Declaration
64 \*---------------------------------------------------------------------------*/
65 
66 class meshToMesh0
67 {
68  // Private Data
69 
70  // Mesh References
71 
72  const fvMesh& fromMesh_;
73  const fvMesh& toMesh_;
74 
75  //- fromMesh patch labels
76  HashTable<label> fromMeshPatches_;
77 
78  //- toMesh patch labels
79  HashTable<label> toMeshPatches_;
80 
81  //- Patch map
82  HashTable<word> patchMap_;
83 
84  //- toMesh patch labels which cut the from-mesh
85  HashTable<label> cuttingPatches_;
86 
87  //- Cell addressing
88  labelList cellAddressing_;
89 
90  //- Boundary addressing
91  labelListList boundaryAddressing_;
92 
93  //- Inverse-distance interpolation weights
94  mutable unique_ptr<scalarListList> inverseDistanceWeightsPtr_;
95 
96  //- Inverse-volume interpolation weights
97  mutable unique_ptr<scalarListList> inverseVolumeWeightsPtr_;
98 
99  //- Cell to cell overlap addressing
100  mutable unique_ptr<labelListList> cellToCellAddressingPtr_;
101 
102  //- Overlap volume
103  mutable scalar V_;
104 
105 
106  // Private Member Functions
107 
108  //- Calculates mesh to mesh addressing pattern.
109  // For each cell from one mesh find the closest cell centre
110  // in the other mesh
111  void calcAddressing();
112 
113  void cellAddresses
114  (
115  labelList& cells,
116  const pointField& points,
117  const fvMesh& fromMesh,
118  const List<bool>& boundaryCell,
120  ) const;
121 
122  void calculateInverseDistanceWeights() const;
123 
124  void calculateInverseVolumeWeights() const;
125 
126  void calculateCellToCellAddressing() const;
127 
128  const scalarListList& inverseDistanceWeights() const;
129 
130  const scalarListList& inverseVolumeWeights() const;
131 
132  const labelListList& cellToCellAddressing() const;
133 
134 
135  // Private Static Data Members
136 
137  //- Direct hit tolerance
138  static const scalar directHitTol;
139 
140 
141 public:
142 
143  // Declare name of the class and its debug switch
144  ClassName("meshToMesh0");
145 
146 
147  //- Enumeration specifying required accuracy
148  enum order
149  {
154  };
155 
156 
157  // Constructors
158 
159  //- Construct from the two meshes, the patch name map for the patches
160  // to be interpolated and the names of the toMesh-patches which
161  // cut the fromMesh
163  (
164  const fvMesh& fromMesh,
165  const fvMesh& toMesh,
166  const HashTable<word>& patchMap,
167  const wordList& cuttingPatchNames
168  );
169 
170  //- Construct from the two meshes assuming there is an exact mapping
171  // between the patches
173  (
174  const fvMesh& fromMesh,
175  const fvMesh& toMesh
176  );
177 
178 
179  //- Destructor
180  ~meshToMesh0() = default;
181 
182 
183  //- Patch-field interpolation class
185  :
186  public fvPatchFieldMapper
187  {
188  const labelList& directAddressing_;
189 
190 
191  public:
192 
193  // Constructors
194 
195  //- Construct given addressing
196  patchFieldInterpolator(const labelList& addr)
197  :
198  directAddressing_(addr)
199  {}
200 
201 
202  //- Destructor
203  virtual ~patchFieldInterpolator() = default;
204 
205 
206  // Member Functions
207 
208  label size() const
209  {
210  return directAddressing_.size();
211  }
212 
213  bool direct() const
214  {
215  return true;
216  }
217 
218  bool hasUnmapped() const
219  {
220  return false;
221  }
222 
223  const labelList& directAddressing() const
224  {
225  return directAddressing_;
226  }
227  };
228 
229 
230  // Member Functions
231 
232  // Access
233 
234  const fvMesh& fromMesh() const
235  {
236  return fromMesh_;
237  }
238 
239  const fvMesh& toMesh() const
240  {
241  return toMesh_;
242  }
243 
244  //- From toMesh cells to fromMesh cells
245  const labelList& cellAddressing() const
246  {
247  return cellAddressing_;
248  }
249 
250  //- Overlap volume
251  scalar V() const
252  {
253  return V_;
254  }
255 
256 
257  // Interpolation
258 
259  //- Map field
260  template<class Type, class CombineOp>
261  void mapField
262  (
263  Field<Type>&,
264  const Field<Type>&,
265  const labelList& adr,
266  const CombineOp& cop
267  ) const;
268 
269  //- Interpolate field using inverse-distance weights
270  template<class Type, class CombineOp>
271  void interpolateField
272  (
273  Field<Type>&,
275  const labelList& adr,
276  const scalarListList& weights,
277  const CombineOp& cop
278  ) const;
279 
280  //- Interpolate field using inverse-volume weights
281  template<class Type, class CombineOp>
282  void interpolateField
283  (
284  Field<Type>&,
286  const labelListList& adr,
287  const scalarListList& weights,
288  const CombineOp& cop
289  ) const;
290 
291 
292  //- Interpolate field using cell-point interpolation
293  template<class Type, class CombineOp>
294  void interpolateField
295  (
296  Field<Type>&,
298  const labelList& adr,
299  const vectorField& centres,
300  const CombineOp& cop
301  )const;
302 
303 
304  //- Interpolate internal volume field
305  template<class Type, class CombineOp>
307  (
308  Field<Type>&,
311  const CombineOp& cop = eqOp<Type>()
312  ) const;
313 
314  template<class Type, class CombineOp>
316  (
317  Field<Type>&,
320  const CombineOp& cop = eqOp<Type>()
321  ) const;
322 
323 
324  //- Interpolate volume field
325  template<class Type, class CombineOp>
326  void interpolate
327  (
331  const CombineOp& cop = eqOp<Type>()
332  ) const;
333 
334  template<class Type, class CombineOp>
335  void interpolate
336  (
340  const CombineOp& cop = eqOp<Type>()
341  ) const;
342 
343 
344  //- Interpolate volume field
345  template<class Type, class CombineOp>
347  (
350  const CombineOp& cop = eqOp<Type>()
351  ) const;
352 
353  template<class Type, class CombineOp>
355  (
358  const CombineOp& cop = eqOp<Type>()
359  ) const;
360 };
361 
362 
363 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
364 
365 } // End namespace Foam
366 
367 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
368 
369 #ifdef NoRepository
370  #include "meshToMesh0Templates.C"
371 #endif
372 
373 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
374 
375 #endif
376 
377 // ************************************************************************* //
Foam::meshToMesh0::patchFieldInterpolator
Patch-field interpolation class.
Definition: meshToMesh0.H:183
Foam::meshToMesh0::INTERPOLATE
Definition: meshToMesh0.H:150
Foam::meshToMesh0::interpolateField
void interpolateField(Field< Type > &, const GeometricField< Type, fvPatchField, volMesh > &, const labelList &adr, const scalarListList &weights, const CombineOp &cop) const
Interpolate field using inverse-distance weights.
Definition: meshToMesh0Templates.C:88
Foam::meshToMesh0::patchFieldInterpolator::patchFieldInterpolator
patchFieldInterpolator(const labelList &addr)
Construct given addressing.
Definition: meshToMesh0.H:195
HashTable.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::meshToMesh0::~meshToMesh0
~meshToMesh0()=default
Destructor.
Foam::meshToMesh0::patchFieldInterpolator::hasUnmapped
bool hasUnmapped() const
Are there unmapped values? I.e. do all size() elements get.
Definition: meshToMesh0.H:217
Foam::meshToMesh0::patchFieldInterpolator::direct
bool direct() const
Definition: meshToMesh0.H:212
scalarList.H
Foam::meshToMesh0::interpolate
void interpolate(GeometricField< Type, fvPatchField, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &, order=INTERPOLATE, const CombineOp &cop=eqOp< Type >()) const
Interpolate volume field.
Definition: meshToMesh0Templates.C:249
Foam::meshToMesh0::order
order
Enumeration specifying required accuracy.
Definition: meshToMesh0.H:147
Foam::meshToMesh0::fromMesh
const fvMesh & fromMesh() const
Definition: meshToMesh0.H:233
fvPatchMapper.H
Foam::meshToMesh0::meshToMesh0
meshToMesh0(const fvMesh &fromMesh, const fvMesh &toMesh, const HashTable< word > &patchMap, const wordList &cuttingPatchNames)
Construct from the two meshes, the patch name map for the patches.
Definition: meshToMesh0.C:45
Foam::meshToMesh0::mapField
void mapField(Field< Type > &, const Field< Type > &, const labelList &adr, const CombineOp &cop) const
Map field.
Definition: meshToMesh0Templates.C:38
Foam::Field< vector >
className.H
Macro definitions for declaring ClassName(), NamespaceName(), etc.
Foam::meshToMesh0::patchFieldInterpolator::size
label size() const
Definition: meshToMesh0.H:207
Foam::indexedOctree
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:50
Foam::meshToMesh0::CELL_POINT_INTERPOLATE
Definition: meshToMesh0.H:151
Foam::meshToMesh0::patchFieldInterpolator::~patchFieldInterpolator
virtual ~patchFieldInterpolator()=default
Destructor.
Foam::meshToMesh0::V
scalar V() const
Overlap volume.
Definition: meshToMesh0.H:250
Foam::eqOp
Definition: ops.H:71
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::meshToMesh0::interpolateInternalField
void interpolateInternalField(Field< Type > &, const GeometricField< Type, fvPatchField, volMesh > &, order=INTERPOLATE, const CombineOp &cop=eqOp< Type >()) const
Interpolate internal volume field.
Definition: meshToMesh0Templates.C:154
Foam::HashTable< label >
Foam::meshToMesh0::ClassName
ClassName("meshToMesh0")
Foam::meshToMesh0::MAP
Definition: meshToMesh0.H:149
meshToMesh0Templates.C
Foam::meshToMesh0::cellAddressing
const labelList & cellAddressing() const
From toMesh cells to fromMesh cells.
Definition: meshToMesh0.H:244
Foam::List< label >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::meshToMesh0::CELL_VOLUME_WEIGHT
Definition: meshToMesh0.H:152
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::meshToMesh0::patchFieldInterpolator::directAddressing
const labelList & directAddressing() const
Definition: meshToMesh0.H:222
Foam::meshToMesh0
Serial mesh to mesh interpolation class.
Definition: meshToMesh0.H:65
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::meshToMesh0::toMesh
const fvMesh & toMesh() const
Definition: meshToMesh0.H:238