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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "meshToMesh0.H"
29 #include "volFields.H"
30 #include "interpolationCellPoint.H"
31 #include "SubField.H"
32 #include "mixedFvPatchField.H"
33 
34 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35 
36 template<class Type, class CombineOp>
38 (
39  Field<Type>& toF,
40  const Field<Type>& fromVf,
41  const labelList& adr,
42  const CombineOp& cop
43 ) const
44 {
45  // Direct mapping of nearest-cell values
46 
47  forAll(toF, celli)
48  {
49  if (adr[celli] != -1)
50  {
51  cop(toF[celli], fromVf[adr[celli]]);
52  }
53  }
54 
55  //toF.map(fromVf, adr);
56 }
57 
58 
59 template<class Type, class CombineOp>
61 (
62  Field<Type>& toF,
64  const labelListList& adr,
65  const scalarListList& weights,
66  const CombineOp& cop
67 ) const
68 {
69  // Inverse volume weighted interpolation
70  forAll(toF, celli)
71  {
72  const labelList& overlapCells = adr[celli];
73  const scalarList& w = weights[celli];
74 
75  Type f = Zero;
76  forAll(overlapCells, i)
77  {
78  label fromCelli = overlapCells[i];
79  f += fromVf[fromCelli]*w[i];
80  cop(toF[celli], f);
81  }
82  }
83 }
84 
85 
86 template<class Type, class CombineOp>
88 (
89  Field<Type>& toF,
91  const labelList& adr,
92  const scalarListList& weights,
93  const CombineOp& cop
94 ) const
95 {
96  // Inverse distance weighted interpolation
97 
98  // get reference to cellCells
99  const labelListList& cc = fromMesh_.cellCells();
100 
101  forAll(toF, celli)
102  {
103  if (adr[celli] != -1)
104  {
105  const labelList& neighbours = cc[adr[celli]];
106  const scalarList& w = weights[celli];
107 
108  Type f = fromVf[adr[celli]]*w[0];
109 
110  for (label ni = 1; ni < w.size(); ni++)
111  {
112  f += fromVf[neighbours[ni - 1]]*w[ni];
113  }
114 
115  cop(toF[celli], f);
116  }
117  }
118 }
119 
120 
121 template<class Type, class CombineOp>
123 (
124  Field<Type>& toF,
126  const labelList& adr,
127  const vectorField& centres,
128  const CombineOp& cop
129 ) const
130 {
131  // Cell-Point interpolation
132  interpolationCellPoint<Type> interpolator(fromVf);
133 
134  forAll(toF, celli)
135  {
136  if (adr[celli] != -1)
137  {
138  cop
139  (
140  toF[celli],
141  interpolator.interpolate
142  (
143  centres[celli],
144  adr[celli]
145  )
146  );
147  }
148  }
149 }
150 
151 
152 template<class Type, class CombineOp>
154 (
155  Field<Type>& toF,
157  meshToMesh0::order ord,
158  const CombineOp& cop
159 ) const
160 {
161  if (fromVf.mesh() != fromMesh_)
162  {
164  << "the argument field does not correspond to the right mesh. "
165  << "Field size: " << fromVf.size()
166  << " mesh size: " << fromMesh_.nCells()
167  << exit(FatalError);
168  }
169 
170  if (toF.size() != toMesh_.nCells())
171  {
173  << "the argument field does not correspond to the right mesh. "
174  << "Field size: " << toF.size()
175  << " mesh size: " << toMesh_.nCells()
176  << exit(FatalError);
177  }
178 
179  switch(ord)
180  {
181  case MAP:
182  mapField(toF, fromVf, cellAddressing_, cop);
183  break;
184 
185  case INTERPOLATE:
186  {
187  interpolateField
188  (
189  toF,
190  fromVf,
191  cellAddressing_,
192  inverseDistanceWeights(),
193  cop
194  );
195  break;
196  }
197  case CELL_POINT_INTERPOLATE:
198  {
199  interpolateField
200  (
201  toF,
202  fromVf,
203  cellAddressing_,
204  toMesh_.cellCentres(),
205  cop
206  );
207 
208  break;
209  }
210  case CELL_VOLUME_WEIGHT:
211  {
212  const labelListList& cellToCell = cellToCellAddressing();
213  const scalarListList& invVolWeights = inverseVolumeWeights();
214 
215  interpolateField
216  (
217  toF,
218  fromVf,
219  cellToCell,
220  invVolWeights,
221  cop
222  );
223  break;
224  }
225  default:
227  << "unknown interpolation scheme " << ord
228  << exit(FatalError);
229  }
230 }
231 
232 
233 template<class Type, class CombineOp>
235 (
236  Field<Type>& toF,
238  meshToMesh0::order ord,
239  const CombineOp& cop
240 ) const
241 {
242  interpolateInternalField(toF, tfromVf(), ord, cop);
243  tfromVf.clear();
244 }
245 
246 
247 template<class Type, class CombineOp>
249 (
252  meshToMesh0::order ord,
253  const CombineOp& cop
254 ) const
255 {
256  interpolateInternalField(toVf, fromVf, ord, cop);
257 
259  Boundary& 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 (
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 (
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],
400  (
401  boundaryAddressing_[patchi]
402  )
403  )
404  );
405  }
406 
407 
408  // Create the complete field from the pieces
410  (
412  (
413  IOobject
414  (
415  "interpolated(" + fromVf.name() + ')',
416  toMesh_.time().timeName(),
417  toMesh_,
418  IOobject::NO_READ,
419  IOobject::NO_WRITE
420  ),
421  toMesh_,
422  fromVf.dimensions(),
423  internalField,
424  patchFields
425  )
426  );
427 
428  return ttoF;
429 }
430 
431 
432 template<class Type, class CombineOp>
435 (
437  meshToMesh0::order ord,
438  const CombineOp& cop
439 ) const
440 {
442  interpolate(tfromVf(), ord, cop);
443  tfromVf.clear();
444 
445  return tint;
446 }
447 
448 
449 // ************************************************************************* //
Foam::meshToMesh0::patchFieldInterpolator
Patch-field interpolation class.
Definition: meshToMesh0.H:183
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
volFields.H
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::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
SubField.H
meshToMesh0.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
mixedFvPatchField.H
Foam::fvPatch::name
virtual const word & name() const
Return name.
Definition: fvPatch.H:167
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
interpolationCellPoint.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::mapField
void mapField(Field< Type > &, const Field< Type > &, const labelList &adr, const CombineOp &cop) const
Map field.
Definition: meshToMesh0Templates.C:38
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::interpolationCellPoint
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
Definition: interpolationCellPoint.H:50
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::cellToCell
A topoSetCellSource to select all the cells from given cellSet(s).
Definition: cellToCell.H:159
Foam::fvPatch::Cf
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:119
Foam::FatalError
error FatalError
Foam::interpolationCellPoint::interpolate
Type interpolate(const cellPointWeight &cpw) const
Interpolate field for the given cellPointWeight.
Definition: interpolationCellPointI.H:32
Foam::mixedFvPatchField
This boundary condition provides a base class for 'mixed' type boundary conditions,...
Definition: mixedFvPatchField.H:123
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::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::isA
const TargetType * isA(const Type &t)
Check if dynamic_cast to TargetType is possible.
Definition: typeInfo.H:197
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
f
labelList f(nPoints)
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::List< label >
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54