meshToMesh0Templates.C
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-2016 OpenFOAM Foundation
9  Copyright (C) 2022 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 \*---------------------------------------------------------------------------*/
28 
29 #include "meshToMesh0.H"
30 #include "volFields.H"
31 #include "interpolationCellPoint.H"
32 #include "SubField.H"
33 #include "mixedFvPatchField.H"
34 
35 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
36 
37 template<class Type, class CombineOp>
39 (
40  Field<Type>& toF,
41  const Field<Type>& fromVf,
42  const labelList& adr,
43  const CombineOp& cop
44 ) const
45 {
46  // Direct mapping of nearest-cell values
47 
48  forAll(toF, celli)
49  {
50  if (adr[celli] != -1)
51  {
52  cop(toF[celli], fromVf[adr[celli]]);
53  }
54  }
55 
56  //toF.map(fromVf, adr);
57 }
58 
59 
60 template<class Type, class CombineOp>
62 (
63  Field<Type>& toF,
64  const VolumeField<Type>& fromVf,
65  const labelListList& adr,
66  const scalarListList& weights,
67  const CombineOp& cop
68 ) const
69 {
70  // Inverse volume weighted interpolation
71  forAll(toF, celli)
72  {
73  const labelList& overlapCells = adr[celli];
74  const scalarList& w = weights[celli];
75 
76  Type f = Zero;
77  forAll(overlapCells, i)
78  {
79  label fromCelli = overlapCells[i];
80  f += fromVf[fromCelli]*w[i];
81  cop(toF[celli], f);
82  }
83  }
84 }
85 
86 
87 template<class Type, class CombineOp>
89 (
90  Field<Type>& toF,
91  const VolumeField<Type>& fromVf,
92  const labelList& adr,
93  const scalarListList& weights,
94  const CombineOp& cop
95 ) const
96 {
97  // Inverse distance weighted interpolation
98 
99  // get reference to cellCells
100  const labelListList& cc = fromMesh_.cellCells();
101 
102  forAll(toF, celli)
103  {
104  if (adr[celli] != -1)
105  {
106  const labelList& neighbours = cc[adr[celli]];
107  const scalarList& w = weights[celli];
108 
109  Type f = fromVf[adr[celli]]*w[0];
110 
111  for (label ni = 1; ni < w.size(); ni++)
112  {
113  f += fromVf[neighbours[ni - 1]]*w[ni];
114  }
115 
116  cop(toF[celli], f);
117  }
118  }
119 }
120 
121 
122 template<class Type, class CombineOp>
124 (
125  Field<Type>& toF,
126  const VolumeField<Type>& fromVf,
127  const labelList& adr,
128  const vectorField& centres,
129  const CombineOp& cop
130 ) const
131 {
132  // Cell-Point interpolation
133  interpolationCellPoint<Type> interpolator(fromVf);
134 
135  forAll(toF, celli)
136  {
137  if (adr[celli] != -1)
138  {
139  cop
140  (
141  toF[celli],
142  interpolator.interpolate
143  (
144  centres[celli],
145  adr[celli]
146  )
147  );
148  }
149  }
150 }
151 
152 
153 template<class Type, class CombineOp>
155 (
156  Field<Type>& toF,
157  const VolumeField<Type>& fromVf,
158  meshToMesh0::order ord,
159  const CombineOp& cop
160 ) const
161 {
162  if (fromVf.mesh() != fromMesh_)
163  {
165  << "the argument field does not correspond to the right mesh. "
166  << "Field size: " << fromVf.size()
167  << " mesh size: " << fromMesh_.nCells()
168  << exit(FatalError);
169  }
170 
171  if (toF.size() != toMesh_.nCells())
172  {
174  << "the argument field does not correspond to the right mesh. "
175  << "Field size: " << toF.size()
176  << " mesh size: " << toMesh_.nCells()
177  << exit(FatalError);
178  }
179 
180  switch(ord)
181  {
182  case MAP:
183  mapField(toF, fromVf, cellAddressing_, cop);
184  break;
185 
186  case INTERPOLATE:
187  {
188  interpolateField
189  (
190  toF,
191  fromVf,
192  cellAddressing_,
193  inverseDistanceWeights(),
194  cop
195  );
196  break;
197  }
198  case CELL_POINT_INTERPOLATE:
199  {
200  interpolateField
201  (
202  toF,
203  fromVf,
204  cellAddressing_,
205  toMesh_.cellCentres(),
206  cop
207  );
208 
209  break;
210  }
211  case CELL_VOLUME_WEIGHT:
212  {
213  const labelListList& cellToCell = cellToCellAddressing();
214  const scalarListList& invVolWeights = inverseVolumeWeights();
215 
216  interpolateField
217  (
218  toF,
219  fromVf,
220  cellToCell,
221  invVolWeights,
222  cop
223  );
224  break;
225  }
226  default:
228  << "unknown interpolation scheme " << ord
230  }
231 }
232 
233 
234 template<class Type, class CombineOp>
236 (
237  Field<Type>& toF,
238  const tmp<VolumeField<Type>>& tfromVf,
239  meshToMesh0::order ord,
240  const CombineOp& cop
241 ) const
242 {
243  interpolateInternalField(toF, tfromVf(), ord, cop);
244  tfromVf.clear();
245 }
246 
247 
248 template<class Type, class CombineOp>
250 (
251  VolumeField<Type>& toVf,
252  const VolumeField<Type>& fromVf,
253  meshToMesh0::order ord,
254  const CombineOp& cop
255 ) const
256 {
257  interpolateInternalField(toVf, fromVf, ord, cop);
258 
259  auto& toVfBf = toVf.boundaryFieldRef();
260 
261  forAll(toMesh_.boundaryMesh(), patchi)
262  {
263  const fvPatch& toPatch = toMesh_.boundary()[patchi];
264 
265  if (cuttingPatches_.found(toPatch.name()))
266  {
267  switch(ord)
268  {
269  case MAP:
270  {
271  mapField
272  (
273  toVfBf[patchi],
274  fromVf,
275  boundaryAddressing_[patchi],
276  cop
277  );
278  break;
279  }
280 
281  case INTERPOLATE:
282  {
283  interpolateField
284  (
285  toVfBf[patchi],
286  fromVf,
287  boundaryAddressing_[patchi],
288  toPatch.Cf(),
289  cop
290  );
291  break;
292  }
293 
294  case CELL_POINT_INTERPOLATE:
295  {
296  interpolateField
297  (
298  toVfBf[patchi],
299  fromVf,
300  boundaryAddressing_[patchi],
301  toPatch.Cf(),
302  cop
303  );
304  break;
305  }
306  case CELL_VOLUME_WEIGHT:
307  {
308  break;
309  }
310 
311  default:
313  << "unknown interpolation scheme " << ord
314  << exit(FatalError);
315  }
316 
317  if (isA<mixedFvPatchField<Type>>(toVfBf[patchi]))
318  {
319  refCast<mixedFvPatchField<Type>>
320  (
321  toVfBf[patchi]
322  ).refValue() = toVfBf[patchi];
323  }
324  }
325  else if
326  (
327  patchMap_.found(toPatch.name())
328  && fromMeshPatches_.found(patchMap_.find(toPatch.name())())
329  )
330  {
331  mapField
332  (
333  toVfBf[patchi],
334  fromVf.boundaryField()
335  [
336  fromMeshPatches_.find(patchMap_.find(toPatch.name())())()
337  ],
338  boundaryAddressing_[patchi],
339  cop
340  );
341  }
342  }
343 }
344 
345 
346 template<class Type, class CombineOp>
348 (
349  VolumeField<Type>& toVf,
350  const tmp<VolumeField<Type>>& tfromVf,
351  meshToMesh0::order ord,
352  const CombineOp& cop
353 ) const
354 {
355  interpolate(toVf, tfromVf(), ord, cop);
356  tfromVf.clear();
357 }
358 
359 
360 template<class Type, class CombineOp>
363 (
364  const VolumeField<Type>& fromVf,
365  meshToMesh0::order ord,
366  const CombineOp& cop
367 ) const
368 {
369  // Create and map the internal-field values
370  Field<Type> internalField(toMesh_.nCells());
371  interpolateInternalField(internalField, fromVf, ord, cop);
372 
373  // check whether both meshes have got the same number
374  // of boundary patches
375  if (fromMesh_.boundary().size() != toMesh_.boundary().size())
376  {
378  << "Incompatible meshes: different number of boundaries, "
379  "only internal field may be interpolated"
380  << exit(FatalError);
381  }
382 
383  // Create and map the patch field values
384  PtrList<fvPatchField<Type>> patchFields
385  (
386  boundaryAddressing_.size()
387  );
388 
389  forAll(boundaryAddressing_, patchi)
390  {
391  patchFields.set
392  (
393  patchi,
395  (
396  fromVf.boundaryField()[patchi],
397  toMesh_.boundary()[patchi],
399  patchFieldInterpolator
400  (
401  boundaryAddressing_[patchi]
402  )
403  )
404  );
405  }
406 
407 
408  // Create the complete field from the pieces
409  return
410  tmp<VolumeField<Type>>::New
411  (
412  IOobject
413  (
414  "interpolated(" + fromVf.name() + ')',
415  toMesh_.time().timeName(),
416  toMesh_,
419  ),
420  toMesh_,
421  fromVf.dimensions(),
422  internalField,
423  patchFields
424  );
425 }
426 
427 
428 template<class Type, class CombineOp>
431 (
432  const tmp<VolumeField<Type>>& tfromVf,
433  meshToMesh0::order ord,
434  const CombineOp& cop
435 ) const
436 {
437  tmp<VolumeField<Type>> tint = interpolate(tfromVf(), ord, cop);
438  tfromVf.clear();
439 
440  return tint;
441 }
442 
443 
444 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
void interpolate(VolumeField< Type > &, const VolumeField< Type > &, order=INTERPOLATE, const CombineOp &cop=eqOp< Type >()) const
Interpolate volume field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
void mapField(Field< Type > &, const Field< Type > &, const labelList &adr, const CombineOp &cop) const
Map field.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
order
Enumeration specifying required accuracy.
Definition: meshToMesh0.H:170
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
Ignore writing from objectRegistry::writeObject()
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static const DimensionedField< Type, GeoMesh > & null()
Return a NullObjectRef DimensionedField.
Generic templated field type.
Definition: Field.H:62
List< scalarList > scalarListList
List of scalarList.
Definition: scalarList.H:35
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:113
labelList f(nPoints)
void interpolateField(Field< Type > &, const VolumeField< Type > &, const labelList &adr, const scalarListList &weights, const CombineOp &cop) const
Interpolate field using inverse-distance weights.
Nothing to be read.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Field< vector > vectorField
Specialisation of Field<T> for vector.
virtual const word & name() const
Return name.
Definition: fvPatch.H:210
const Type * isA(const U &obj)
Check if dynamic_cast to Type is possible.
Definition: typeInfo.H:88
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void interpolateInternalField(Field< Type > &, const VolumeField< Type > &, order=INTERPOLATE, const CombineOp &cop=eqOp< Type >()) const
Interpolate internal volume field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
static tmp< fvPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
const dimensionSet & dimensions() const noexcept
Return dimensions.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127