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 | Copyright (C) 2011-2015 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, class CombineOp>
36 (
37  Field<Type>& toF,
38  const Field<Type>& fromVf,
39  const labelList& adr,
40  const CombineOp& cop
41 ) const
42 {
43  // Direct mapping of nearest-cell values
44 
45  forAll(toF, celli)
46  {
47  if (adr[celli] != -1)
48  {
49  cop(toF[celli], fromVf[adr[celli]]);
50  }
51  }
52 
53  //toF.map(fromVf, adr);
54 }
55 
56 
57 template<class Type, class CombineOp>
59 (
60  Field<Type>& toF,
62  const labelListList& adr,
63  const scalarListList& weights,
64  const CombineOp& cop
65 ) const
66 {
67  // Inverse volume weighted interpolation
68  forAll(toF, celli)
69  {
70  const labelList& overlapCells = adr[celli];
71  const scalarList& w = weights[celli];
72 
73  Type f = pTraits<Type>::zero;
74  forAll(overlapCells, i)
75  {
76  label fromCelli = overlapCells[i];
77  f += fromVf[fromCelli]*w[i];
78  cop(toF[celli], f);
79  }
80  }
81 }
82 
83 
84 template<class Type, class CombineOp>
86 (
87  Field<Type>& toF,
89  const labelList& adr,
90  const scalarListList& weights,
91  const CombineOp& cop
92 ) const
93 {
94  // Inverse distance weighted interpolation
95 
96  // get reference to cellCells
97  const labelListList& cc = fromMesh_.cellCells();
98 
99  forAll(toF, celli)
100  {
101  if (adr[celli] != -1)
102  {
103  const labelList& neighbours = cc[adr[celli]];
104  const scalarList& w = weights[celli];
105 
106  Type f = fromVf[adr[celli]]*w[0];
107 
108  for (label ni = 1; ni < w.size(); ni++)
109  {
110  f += fromVf[neighbours[ni - 1]]*w[ni];
111  }
112 
113  cop(toF[celli], f);
114  }
115  }
116 }
117 
118 
119 template<class Type, class CombineOp>
121 (
122  Field<Type>& toF,
124  const labelList& adr,
125  const vectorField& centres,
126  const CombineOp& cop
127 ) const
128 {
129  // Cell-Point interpolation
130  interpolationCellPoint<Type> interpolator(fromVf);
131 
132  forAll(toF, celli)
133  {
134  if (adr[celli] != -1)
135  {
136  cop
137  (
138  toF[celli],
139  interpolator.interpolate
140  (
141  centres[celli],
142  adr[celli]
143  )
144  );
145  }
146  }
147 }
148 
149 
150 template<class Type, class CombineOp>
152 (
153  Field<Type>& toF,
155  meshToMesh0::order ord,
156  const CombineOp& cop
157 ) const
158 {
159  if (fromVf.mesh() != fromMesh_)
160  {
162  << "the argument field does not correspond to the right mesh. "
163  << "Field size: " << fromVf.size()
164  << " mesh size: " << fromMesh_.nCells()
165  << exit(FatalError);
166  }
167 
168  if (toF.size() != toMesh_.nCells())
169  {
171  << "the argument field does not correspond to the right mesh. "
172  << "Field size: " << toF.size()
173  << " mesh size: " << toMesh_.nCells()
174  << exit(FatalError);
175  }
176 
177  switch(ord)
178  {
179  case MAP:
180  mapField(toF, fromVf, cellAddressing_, cop);
181  break;
182 
183  case INTERPOLATE:
184  {
185  interpolateField
186  (
187  toF,
188  fromVf,
189  cellAddressing_,
190  inverseDistanceWeights(),
191  cop
192  );
193  break;
194  }
195  case CELL_POINT_INTERPOLATE:
196  {
197  interpolateField
198  (
199  toF,
200  fromVf,
201  cellAddressing_,
202  toMesh_.cellCentres(),
203  cop
204  );
205 
206  break;
207  }
208  case CELL_VOLUME_WEIGHT:
209  {
210  const labelListList& cellToCell = cellToCellAddressing();
211  const scalarListList& invVolWeights = inverseVolumeWeights();
212 
213  interpolateField
214  (
215  toF,
216  fromVf,
217  cellToCell,
218  invVolWeights,
219  cop
220  );
221  break;
222  }
223  default:
225  << "unknown interpolation scheme " << ord
226  << exit(FatalError);
227  }
228 }
229 
230 
231 template<class Type, class CombineOp>
233 (
234  Field<Type>& toF,
236  meshToMesh0::order ord,
237  const CombineOp& cop
238 ) const
239 {
240  interpolateInternalField(toF, tfromVf(), ord, cop);
241  tfromVf.clear();
242 }
243 
244 
245 template<class Type, class CombineOp>
247 (
250  meshToMesh0::order ord,
251  const CombineOp& cop
252 ) const
253 {
254  interpolateInternalField(toVf, fromVf, ord, cop);
255 
256  forAll(toMesh_.boundaryMesh(), patchi)
257  {
258  const fvPatch& toPatch = toMesh_.boundary()[patchi];
259 
260  if (cuttingPatches_.found(toPatch.name()))
261  {
262  switch(ord)
263  {
264  case MAP:
265  {
266  mapField
267  (
268  toVf.boundaryField()[patchi],
269  fromVf,
270  boundaryAddressing_[patchi],
271  cop
272  );
273  break;
274  }
275 
276  case INTERPOLATE:
277  {
278  interpolateField
279  (
280  toVf.boundaryField()[patchi],
281  fromVf,
282  boundaryAddressing_[patchi],
283  toPatch.Cf(),
284  cop
285  );
286  break;
287  }
288 
289  case CELL_POINT_INTERPOLATE:
290  {
291  interpolateField
292  (
293  toVf.boundaryField()[patchi],
294  fromVf,
295  boundaryAddressing_[patchi],
296  toPatch.Cf(),
297  cop
298  );
299  break;
300  }
301  case CELL_VOLUME_WEIGHT:
302  {
303  // Do nothing
304  break;
305  }
306 
307  default:
309  << "unknown interpolation scheme " << ord
310  << exit(FatalError);
311  }
312 
314  {
315  refCast<mixedFvPatchField<Type> >
316  (
317  toVf.boundaryField()[patchi]
318  ).refValue() = toVf.boundaryField()[patchi];
319  }
320  }
321  else if
322  (
323  patchMap_.found(toPatch.name())
324  && fromMeshPatches_.found(patchMap_.find(toPatch.name())())
325  )
326  {
327  /*
328  toVf.boundaryField()[patchi].map
329  (
330  fromVf.boundaryField()
331  [
332  fromMeshPatches_.find(patchMap_.find(toPatch.name())())()
333  ],
334  boundaryAddressing_[patchi]
335  );
336  */
337 
338  mapField
339  (
340  toVf.boundaryField()[patchi],
341  fromVf.boundaryField()
342  [
343  fromMeshPatches_.find(patchMap_.find(toPatch.name())())()
344  ],
345  boundaryAddressing_[patchi],
346  cop
347  );
348  }
349  }
350 }
351 
352 
353 template<class Type, class CombineOp>
355 (
358  meshToMesh0::order ord,
359  const CombineOp& cop
360 ) const
361 {
362  interpolate(toVf, tfromVf(), ord, cop);
363  tfromVf.clear();
364 }
365 
366 
367 template<class Type, class CombineOp>
370 (
372  meshToMesh0::order ord,
373  const CombineOp& cop
374 ) const
375 {
376  // Create and map the internal-field values
377  Field<Type> internalField(toMesh_.nCells());
378  interpolateInternalField(internalField, fromVf, ord, cop);
379 
380  // check whether both meshes have got the same number
381  // of boundary patches
382  if (fromMesh_.boundary().size() != toMesh_.boundary().size())
383  {
385  << "Incompatible meshes: different number of boundaries, "
386  "only internal field may be interpolated"
387  << exit(FatalError);
388  }
389 
390  // Create and map the patch field values
391  PtrList<fvPatchField<Type> > patchFields
392  (
393  boundaryAddressing_.size()
394  );
395 
396  forAll(boundaryAddressing_, patchI)
397  {
398  patchFields.set
399  (
400  patchI,
402  (
403  fromVf.boundaryField()[patchI],
404  toMesh_.boundary()[patchI],
407  (
408  boundaryAddressing_[patchI]
409  )
410  )
411  );
412  }
413 
414 
415  // Create the complete field from the pieces
417  (
419  (
420  IOobject
421  (
422  "interpolated(" + fromVf.name() + ')',
423  toMesh_.time().timeName(),
424  toMesh_,
425  IOobject::NO_READ,
426  IOobject::NO_WRITE
427  ),
428  toMesh_,
429  fromVf.dimensions(),
431  patchFields
432  )
433  );
434 
435  return ttoF;
436 }
437 
438 
439 template<class Type, class CombineOp>
442 (
444  meshToMesh0::order ord,
445  const CombineOp& cop
446 ) const
447 {
449  interpolate(tfromVf(), ord, cop);
450  tfromVf.clear();
451 
452  return tint;
453 }
454 
455 
456 // ************************************************************************* //
Foam::meshToMesh0::patchFieldInterpolator
Patch-field interpolation class.
Definition: meshToMesh0.H:178
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
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:86
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
w
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
SubField.H
meshToMesh0.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
mixedFvPatchField.H
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:247
Foam::meshToMesh0::order
order
Enumeration specifying required accuracy.
Definition: meshToMesh0.H:142
Foam::fvPatch::name
const word & name() const
Return name.
Definition: fvPatch.H:149
Foam::meshToMesh0::mapField
void mapField(Field< Type > &, const Field< Type > &, const labelList &adr, const CombineOp &cop) const
Map field.
Definition: meshToMesh0Templates.C:36
Foam::PtrList::set
bool set(const label) const
Is element set.
Foam::label
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
Foam::Field< Type >
Foam::interpolationCellPoint
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
Definition: interpolationCellPoint.H:48
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::cellToCell
A topoSetSource to select the cells from another cellSet.
Definition: cellToCell.H:48
Foam::fvPatch::Cf
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:99
Foam::FatalError
error FatalError
Foam::interpolationCellPoint::interpolate
Type interpolate(const cellPointWeight &cpw) const
Interpolate field for the given cellPointWeight.
Definition: interpolationCellPointI.H:30
Foam::mixedFvPatchField
This boundary condition provides a base class for 'mixed' type boundary conditions,...
Definition: mixedFvPatchField.H:126
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:152
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
internalField
conserve internalField()+
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
f
labelList f(nPoints)
Foam::isA
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51