meshToMesh0Templates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "meshToMesh0.H"
27 #include "volFields.H"
28 #include "interpolationCellPoint.H"
29 #include "SubField.H"
30 #include "mixedFvPatchField.H"
31 
32 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
33 
34 template<class Type>
36 (
37  Field<Type>& toF,
38  const Field<Type>& fromVf,
39  const labelList& adr
40 ) const
41 {
42  // Direct mapping of nearest-cell values
43 
44  forAll(toF, celli)
45  {
46  if (adr[celli] != -1)
47  {
48  toF[celli] = fromVf[adr[celli]];
49  }
50  }
51 
52  // toF.map(fromVf, adr);
53 }
54 
55 
56 template<class Type>
58 (
59  Field<Type>& toF,
61  const labelListList& adr,
62  const scalarListList& weights
63 ) const
64 {
65  // Inverse volume weighted interpolation
66  forAll(toF, celli)
67  {
68  const labelList& overlapCells = adr[celli];
69  const scalarList& w = weights[celli];
70 
71  Type f = Zero;
72  forAll(overlapCells, i)
73  {
74  label fromCelli = overlapCells[i];
75  f += fromVf[fromCelli]*w[i];
76  toF[celli] = f;
77  }
78  }
79 }
80 
81 
82 template<class Type>
84 (
85  Field<Type>& toF,
87  const labelList& adr,
88  const scalarListList& weights
89 ) const
90 {
91  // Inverse distance weighted interpolation
92 
93  // get reference to cellCells
94  const labelListList& cc = fromMesh_.cellCells();
95 
96  forAll(toF, celli)
97  {
98  if (adr[celli] != -1)
99  {
100  const labelList& neighbours = cc[adr[celli]];
101  const scalarList& w = weights[celli];
102 
103  Type f = fromVf[adr[celli]]*w[0];
104 
105  for (label ni = 1; ni < w.size(); ni++)
106  {
107  f += fromVf[neighbours[ni - 1]]*w[ni];
108  }
109 
110  toF[celli] = f;
111  }
112  }
113 }
114 
115 
116 template<class Type>
118 (
119  Field<Type>& toF,
121  const labelList& adr,
122  const vectorField& centres
123 ) const
124 {
125  // Cell-Point interpolation
126  interpolationCellPoint<Type> interpolator(fromVf);
127 
128  forAll(toF, celli)
129  {
130  if (adr[celli] != -1)
131  {
132  toF[celli] = interpolator.interpolate(centres[celli], adr[celli]);
133  }
134  }
135 }
136 
137 
138 template<class Type>
140 (
141  Field<Type>& toF,
144 ) const
145 {
146  if (fromVf.mesh() != fromMesh_)
147  {
149  << "the argument field does not correspond to the right mesh. "
150  << "Field size: " << fromVf.size()
151  << " mesh size: " << fromMesh_.nCells()
152  << exit(FatalError);
153  }
154 
155  if (toF.size() != toMesh_.nCells())
156  {
158  << "the argument field does not correspond to the right mesh. "
159  << "Field size: " << toF.size()
160  << " mesh size: " << toMesh_.nCells()
161  << exit(FatalError);
162  }
163 
164  switch(ord)
165  {
166  case MAP:
167  mapField(toF, fromVf, cellAddressing_);
168  break;
169 
170  case INTERPOLATE:
171  {
172  interpolateField
173  (
174  toF,
175  fromVf,
176  cellAddressing_,
177  inverseDistanceWeights()
178  );
179  break;
180  }
181  case CELL_POINT_INTERPOLATE:
182  {
183  interpolateField
184  (
185  toF,
186  fromVf,
187  cellAddressing_,
188  toMesh_.cellCentres()
189  );
190 
191  break;
192  }
193  case CELL_VOLUME_WEIGHT:
194  {
195  const labelListList& cellToCell = cellToCellAddressing();
196  const scalarListList& invVolWeights = inverseVolumeWeights();
197 
198  interpolateField
199  (
200  toF,
201  fromVf,
202  cellToCell,
203  invVolWeights
204  );
205  break;
206  }
207  default:
209  << "unknown interpolation scheme " << ord
210  << exit(FatalError);
211  }
212 }
213 
214 
215 template<class Type>
217 (
218  Field<Type>& toF,
221 ) const
222 {
223  interpolateInternalField(toF, tfromVf(), ord);
224  tfromVf.clear();
225 }
226 
227 
228 template<class Type>
230 (
234 ) const
235 {
236  interpolateInternalField(toVf, fromVf, ord);
237 
239  Boundary& toVfBf = toVf.boundaryFieldRef();
240 
241  forAll(toMesh_.boundaryMesh(), patchi)
242  {
243  const fvPatch& toPatch = toMesh_.boundary()[patchi];
244 
245  if (cuttingPatches_.found(toPatch.name()))
246  {
247  switch(ord)
248  {
249  case MAP:
250  {
251  mapField
252  (
253  toVfBf[patchi],
254  fromVf,
255  boundaryAddressing_[patchi]
256  );
257  break;
258  }
259 
260  case INTERPOLATE:
261  {
262  interpolateField
263  (
264  toVfBf[patchi],
265  fromVf,
266  boundaryAddressing_[patchi],
267  toPatch.Cf()
268  );
269  break;
270  }
271 
272  case CELL_POINT_INTERPOLATE:
273  {
274  interpolateField
275  (
276  toVfBf[patchi],
277  fromVf,
278  boundaryAddressing_[patchi],
279  toPatch.Cf()
280  );
281  break;
282  }
283  case CELL_VOLUME_WEIGHT:
284  {
285  break;
286  }
287 
288  default:
290  << "unknown interpolation scheme " << ord
291  << exit(FatalError);
292  }
293 
294  if (isA<mixedFvPatchField<Type>>(toVfBf[patchi]))
295  {
296  refCast<mixedFvPatchField<Type>>
297  (
298  toVfBf[patchi]
299  ).refValue() = toVfBf[patchi];
300  }
301  }
302  else if
303  (
304  patchMap_.found(toPatch.name())
305  && fromMeshPatches_.found(patchMap_.find(toPatch.name())())
306  )
307  {
308  mapField
309  (
310  toVfBf[patchi],
311  fromVf.boundaryField()
312  [
313  fromMeshPatches_.find(patchMap_.find(toPatch.name())())()
314  ],
315  boundaryAddressing_[patchi]
316  );
317  }
318  }
319 }
320 
321 
322 template<class Type>
324 (
328 ) const
329 {
330  interpolate(toVf, tfromVf(), ord);
331  tfromVf.clear();
332 }
333 
334 
335 template<class Type>
338 (
341 ) const
342 {
343  // Create and map the internal-field values
344  Field<Type> internalField(toMesh_.nCells());
345  interpolateInternalField(internalField, fromVf, ord);
346 
347  // check whether both meshes have got the same number
348  // of boundary patches
349  if (fromMesh_.boundary().size() != toMesh_.boundary().size())
350  {
352  << "Incompatible meshes: different number of boundaries, "
353  "only internal field may be interpolated"
354  << exit(FatalError);
355  }
356 
357  // Create and map the patch field values
358  PtrList<fvPatchField<Type>> patchFields
359  (
360  boundaryAddressing_.size()
361  );
362 
363  forAll(boundaryAddressing_, patchi)
364  {
365  patchFields.set
366  (
367  patchi,
369  (
370  fromVf.boundaryField()[patchi],
371  toMesh_.boundary()[patchi],
374  (
375  boundaryAddressing_[patchi]
376  )
377  )
378  );
379  }
380 
381 
382  // Create the complete field from the pieces
384  (
386  (
387  IOobject
388  (
389  "interpolated(" + fromVf.name() + ')',
390  toMesh_.time().timeName(),
391  toMesh_,
392  IOobject::NO_READ,
393  IOobject::NO_WRITE,
394  false
395  ),
396  toMesh_,
397  fromVf.dimensions(),
398  internalField,
399  patchFields
400  )
401  );
402 
403  return ttoF;
404 }
405 
406 
407 template<class Type>
410 (
413 ) const
414 {
416  interpolate(tfromVf(), ord);
417  tfromVf.clear();
418 
419  return tint;
420 }
421 
422 
423 // ************************************************************************* //
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const word & name() const
Return name.
Definition: IOobject.H:315
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
This boundary condition provides a base class for &#39;mixed&#39; type boundary conditions, i.e. conditions that mix fixed value and patch-normal gradient conditions.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
order
Enumeration specifying required accuracy.
Definition: meshToMesh0.H:143
void interpolateInternalField(Field< Type > &, const GeometricField< Type, fvPatchField, volMesh > &, order=INTERPOLATE) const
Interpolate internal volume field.
Generic GeometricField class.
A topoSetSource to select the cells from another cellSet.
Definition: cellToCell.H:48
const dimensionSet & dimensions() const
Return dimensions.
Pre-declare SubField and related Field type.
Definition: Field.H:56
void mapField(Field< Type > &, const Field< Type > &, const labelList &adr) const
Map field.
Type interpolate(const cellPointWeight &cpw) const
Interpolate field for the given cellPointWeight.
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:105
static const zero Zero
Definition: zero.H:97
const Mesh & mesh() const
Return mesh.
labelList f(nPoints)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
label patchi
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
Patch-field interpolation class.
Definition: meshToMesh0.H:179
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void interpolateField(Field< Type > &, const GeometricField< Type, fvPatchField, volMesh > &, const labelList &adr, const scalarListList &weights) const
Interpolate field using inverse-distance weights.
virtual const word & name() const
Return name.
Definition: fvPatch.H:145
void interpolate(GeometricField< Type, fvPatchField, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &, order=INTERPOLATE) const
Interpolate volume field.
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98