Go to the documentation of this file.
34 template<
class Type,
class CombineOp>
49 cop(toF[celli], fromVf[adr[celli]]);
57 template<
class Type,
class CombineOp>
70 const labelList& overlapCells = adr[celli];
76 label fromCelli = overlapCells[i];
77 f += fromVf[fromCelli]*
w[i];
84 template<
class Type,
class CombineOp>
101 if (adr[celli] != -1)
103 const labelList& neighbours = cc[adr[celli]];
106 Type
f = fromVf[adr[celli]]*
w[0];
108 for (
label ni = 1; ni <
w.size(); ni++)
110 f += fromVf[neighbours[ni - 1]]*
w[ni];
119 template<
class Type,
class CombineOp>
134 if (adr[celli] != -1)
150 template<
class Type,
class CombineOp>
159 if (fromVf.mesh() != fromMesh_)
162 <<
"the argument field does not correspond to the right mesh. "
163 <<
"Field size: " << fromVf.size()
164 <<
" mesh size: " << fromMesh_.nCells()
168 if (toF.size() != toMesh_.nCells())
171 <<
"the argument field does not correspond to the right mesh. "
172 <<
"Field size: " << toF.size()
173 <<
" mesh size: " << toMesh_.nCells()
180 mapField(toF, fromVf, cellAddressing_, cop);
190 inverseDistanceWeights(),
195 case CELL_POINT_INTERPOLATE:
202 toMesh_.cellCentres(),
208 case CELL_VOLUME_WEIGHT:
225 <<
"unknown interpolation scheme " << ord
231 template<
class Type,
class CombineOp>
240 interpolateInternalField(toF, tfromVf(), ord, cop);
245 template<
class Type,
class CombineOp>
254 interpolateInternalField(toVf, fromVf, ord, cop);
260 if (cuttingPatches_.found(toPatch.
name()))
270 boundaryAddressing_[
patchi],
282 boundaryAddressing_[
patchi],
289 case CELL_POINT_INTERPOLATE:
295 boundaryAddressing_[
patchi],
301 case CELL_VOLUME_WEIGHT:
309 <<
"unknown interpolation scheme " << ord
315 refCast<mixedFvPatchField<Type> >
323 patchMap_.found(toPatch.
name())
324 && fromMeshPatches_.found(patchMap_.find(toPatch.
name())())
343 fromMeshPatches_.find(patchMap_.find(toPatch.
name())())()
345 boundaryAddressing_[
patchi],
353 template<
class Type,
class CombineOp>
367 template<
class Type,
class CombineOp>
382 if (fromMesh_.boundary().size() != toMesh_.boundary().size())
385 <<
"Incompatible meshes: different number of boundaries, "
386 "only internal field may be interpolated"
393 boundaryAddressing_.size()
396 forAll(boundaryAddressing_, patchI)
404 toMesh_.boundary()[patchI],
408 boundaryAddressing_[patchI]
422 "interpolated(" + fromVf.name() +
')',
423 toMesh_.time().timeName(),
439 template<
class Type,
class CombineOp>
Patch-field interpolation class.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
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.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
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))
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
tmp< surfaceScalarField > interpolate(const RhoType &rho)
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
void interpolate(GeometricField< Type, fvPatchField, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &, order=INTERPOLATE, const CombineOp &cop=eqOp< Type >()) const
Interpolate volume field.
order
Enumeration specifying required accuracy.
const word & name() const
Return name.
void mapField(Field< Type > &, const Field< Type > &, const labelList &adr, const CombineOp &cop) const
Map field.
bool set(const label) const
Is element set.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
A topoSetSource to select the cells from another cellSet.
const vectorField & Cf() const
Return face centres.
Type interpolate(const cellPointWeight &cpw) const
Interpolate field for the given cellPointWeight.
This boundary condition provides a base class for 'mixed' type boundary conditions,...
void interpolateInternalField(Field< Type > &, const GeometricField< Type, fvPatchField, volMesh > &, order=INTERPOLATE, const CombineOp &cop=eqOp< Type >()) const
Interpolate internal volume field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
conserve internalField()+
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Traits class for primitives.
Generic GeometricField class.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...