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-2023 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 #include "forwardFieldMapper.H"
32 
33 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  Field<Type>& toF,
39  const Field<Type>& fromVf,
40  const labelList& adr
41 ) const
42 {
43  // Direct mapping of nearest-cell values
44 
45  forAll(toF, celli)
46  {
47  if (adr[celli] != -1)
48  {
49  toF[celli] = fromVf[adr[celli]];
50  }
51  }
52 
53  // toF.map(fromVf, adr);
54 }
55 
56 
57 template<class Type>
59 (
60  Field<Type>& toF,
61  const VolField<Type>& fromVf,
62  const labelListList& adr,
63  const scalarListList& weights
64 ) const
65 {
66  // Inverse volume weighted interpolation
67  forAll(toF, celli)
68  {
69  const labelList& overlapCells = adr[celli];
70  const scalarList& w = weights[celli];
71 
72  Type f = Zero;
73  forAll(overlapCells, i)
74  {
75  label fromCelli = overlapCells[i];
76  f += fromVf[fromCelli]*w[i];
77  toF[celli] = f;
78  }
79  }
80 }
81 
82 
83 template<class Type>
85 (
86  Field<Type>& toF,
87  const VolField<Type>& fromVf,
88  const labelList& adr,
89  const scalarListList& weights
90 ) const
91 {
92  // Inverse distance weighted interpolation
93 
94  // get reference to cellCells
95  const labelListList& cc = fromMesh_.cellCells();
96 
97  forAll(toF, celli)
98  {
99  if (adr[celli] != -1)
100  {
101  const labelList& neighbours = cc[adr[celli]];
102  const scalarList& w = weights[celli];
103 
104  Type f = fromVf[adr[celli]]*w[0];
105 
106  for (label ni = 1; ni < w.size(); ni++)
107  {
108  f += fromVf[neighbours[ni - 1]]*w[ni];
109  }
110 
111  toF[celli] = f;
112  }
113  }
114 }
115 
116 
117 template<class Type>
119 (
120  Field<Type>& toF,
121  const VolField<Type>& fromVf,
122  const labelList& adr,
123  const vectorField& centres
124 ) const
125 {
126  // Cell-Point interpolation
127  interpolationCellPoint<Type> interpolator(fromVf);
128 
129  forAll(toF, celli)
130  {
131  if (adr[celli] != -1)
132  {
133  toF[celli] = interpolator.interpolate(centres[celli], adr[celli]);
134  }
135  }
136 }
137 
138 
139 template<class Type>
141 (
142  Field<Type>& toF,
143  const VolField<Type>& fromVf,
145 ) const
146 {
147  if (fromVf.mesh() != fromMesh_)
148  {
150  << "the argument field does not correspond to the right mesh. "
151  << "Field size: " << fromVf.size()
152  << " mesh size: " << fromMesh_.nCells()
153  << exit(FatalError);
154  }
155 
156  if (toF.size() != toMesh_.nCells())
157  {
159  << "the argument field does not correspond to the right mesh. "
160  << "Field size: " << toF.size()
161  << " mesh size: " << toMesh_.nCells()
162  << exit(FatalError);
163  }
164 
165  switch(ord)
166  {
167  case MAP:
168  mapField(toF, fromVf, cellAddressing_);
169  break;
170 
171  case INTERPOLATE:
172  {
173  interpolateField
174  (
175  toF,
176  fromVf,
177  cellAddressing_,
178  inverseDistanceWeights()
179  );
180  break;
181  }
182  case CELL_POINT_INTERPOLATE:
183  {
184  interpolateField
185  (
186  toF,
187  fromVf,
188  cellAddressing_,
189  toMesh_.cellCentres()
190  );
191 
192  break;
193  }
194  case CELL_VOLUME_WEIGHT:
195  {
196  const labelListList& cellToCell = cellToCellAddressing();
197  const scalarListList& invVolWeights = inverseVolumeWeights();
198 
199  interpolateField
200  (
201  toF,
202  fromVf,
203  cellToCell,
204  invVolWeights
205  );
206  break;
207  }
208  default:
210  << "unknown interpolation scheme " << ord
211  << exit(FatalError);
212  }
213 }
214 
215 
216 template<class Type>
218 (
219  Field<Type>& toF,
220  const tmp<VolField<Type>>& tfromVf,
222 ) const
223 {
224  interpolateInternalField(toF, tfromVf(), ord);
225  tfromVf.clear();
226 }
227 
228 
229 template<class Type>
231 (
232  VolField<Type>& toVf,
233  const VolField<Type>& fromVf,
235 ) const
236 {
237  interpolateInternalField(toVf, fromVf, ord);
238 
239  typename VolField<Type>::
240  Boundary& toVfBf = toVf.boundaryFieldRef();
241 
242  forAll(toMesh_.boundaryMesh(), patchi)
243  {
244  const fvPatch& toPatch = toMesh_.boundary()[patchi];
245 
246  if (cuttingPatches_.found(toPatch.name()))
247  {
248  switch(ord)
249  {
250  case MAP:
251  {
252  mapField
253  (
254  toVfBf[patchi],
255  fromVf,
256  boundaryAddressing_[patchi]
257  );
258  break;
259  }
260 
261  case INTERPOLATE:
262  {
263  interpolateField
264  (
265  toVfBf[patchi],
266  fromVf,
267  boundaryAddressing_[patchi],
268  toPatch.Cf()
269  );
270  break;
271  }
272 
273  case CELL_POINT_INTERPOLATE:
274  {
275  interpolateField
276  (
277  toVfBf[patchi],
278  fromVf,
279  boundaryAddressing_[patchi],
280  toPatch.Cf()
281  );
282  break;
283  }
284  case CELL_VOLUME_WEIGHT:
285  {
286  break;
287  }
288 
289  default:
291  << "unknown interpolation scheme " << ord
292  << exit(FatalError);
293  }
294 
295  if (isA<mixedFvPatchField<Type>>(toVfBf[patchi]))
296  {
297  refCast<mixedFvPatchField<Type>>
298  (
299  toVfBf[patchi]
300  ).refValue() = toVfBf[patchi];
301  }
302  }
303  else if
304  (
305  patchMap_.found(toPatch.name())
306  && fromMeshPatches_.found(patchMap_.find(toPatch.name())())
307  )
308  {
309  mapField
310  (
311  toVfBf[patchi],
312  fromVf.boundaryField()
313  [
314  fromMeshPatches_.find(patchMap_.find(toPatch.name())())()
315  ],
316  boundaryAddressing_[patchi]
317  );
318  }
319  }
320 }
321 
322 
323 template<class Type>
325 (
326  VolField<Type>& toVf,
327  const tmp<VolField<Type>>& tfromVf,
329 ) const
330 {
331  interpolate(toVf, tfromVf(), ord);
332  tfromVf.clear();
333 }
334 
335 
336 template<class Type>
339 (
340  const VolField<Type>& fromVf,
342 ) const
343 {
344  // Create and map the internal-field values
345  Field<Type> internalField(toMesh_.nCells());
346  interpolateInternalField(internalField, fromVf, ord);
347 
348  // check whether both meshes have got the same number
349  // of boundary patches
350  if (fromMesh_.boundary().size() != toMesh_.boundary().size())
351  {
353  << "Incompatible meshes: different number of boundaries, "
354  "only internal field may be interpolated"
355  << exit(FatalError);
356  }
357 
358  // Create and map the patch field values
359  PtrList<fvPatchField<Type>> patchFields
360  (
361  boundaryAddressing_.size()
362  );
363 
364  forAll(boundaryAddressing_, patchi)
365  {
366  patchFields.set
367  (
368  patchi,
370  (
371  fromVf.boundaryField()[patchi],
372  toMesh_.boundary()[patchi],
374  forwardFieldMapper(boundaryAddressing_[patchi])
375  )
376  );
377  }
378 
379 
380  // Create the complete field from the pieces
381  tmp<VolField<Type>> ttoF
382  (
383  new VolField<Type>
384  (
385  IOobject
386  (
387  "interpolated(" + fromVf.name() + ')',
388  toMesh_.time().name(),
389  toMesh_,
392  false
393  ),
394  toMesh_,
395  fromVf.dimensions(),
396  internalField,
397  patchFields
398  )
399  );
400 
401  return ttoF;
402 }
403 
404 
405 template<class Type>
408 (
409  const tmp<VolField<Type>>& tfromVf,
411 ) const
412 {
413  tmp<VolField<Type>> tint =
414  interpolate(tfromVf(), ord);
415  tfromVf.clear();
416 
417  return tint;
418 }
419 
420 
421 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
static tmp< fvPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
void interpolateField(Field< Type > &, const VolField< Type > &, const labelList &adr, const scalarListList &weights) const
Interpolate field using inverse-distance weights.
void interpolateInternalField(Field< Type > &, const VolField< Type > &, order=INTERPOLATE) const
Interpolate internal volume field.
order
Enumeration specifying required accuracy.
Definition: meshToMesh0.H:143
void mapField(Field< Type > &, const Field< Type > &, const labelList &adr) const
Map field.
void interpolate(VolField< Type > &, const VolField< Type > &, order=INTERPOLATE) const
Interpolate volume field.
A class for managing temporary objects.
Definition: tmp.H:55
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
List< scalarList > scalarListList
Definition: scalarList.H:51
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:166
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
Field< vector > vectorField
Specialisation of Field<T> for vector.
error FatalError
labelList f(nPoints)