Go to the documentation of this file.
61 mesh.pointsInstance(),
67 auto& oDelta = *oDelta_;
74 mesh.pointsInstance(),
80 auto& nDelta = *nDelta_;
95 mag(
n[facei] & (
C[owner[facei]] - Cf[facei]));
97 mag(
n[facei] & (
C[neighbour[facei]] - Cf[facei]));
104 const fvPatch& currPatch =
mesh.boundary()[patchi];
107 const vectorField nPatch(currPatch.Sf()/currPatch.magSf());
110 if (currPatch.coupled())
117 const label own = pOwner[facei];
120 oDelta.boundaryFieldRef()[patchi][facei] =
121 mag(nPatch[facei] & (pCf[facei] -
C[own]));
125 nDelta.boundaryFieldRef()[patchi] =
126 currPatch.weights()*oDelta.boundaryFieldRef()[patchi]
127 /(scalar(1) - currPatch.weights());
136 const label own = pOwner[facei];
139 oDelta.boundaryFieldRef()[patchi][facei] =
140 mag(nPatch[facei] & (pCf[facei] -
C[own]));
142 nDelta.boundaryFieldRef()[patchi][facei] =
143 mag(nPatch[facei] & (pCf[facei] -
C[own]));
165 "weightedFlux::interpolate(" + vf.name() +
')',
166 mesh.time().timeName(),
172 auto& result = tresult.ref();
179 const scalar sigmaDeltaO = sigma_[owner[facei]]/oDelta[facei];
180 const scalar sigmaDeltaN = sigma_[neighbour[facei]]/nDelta[facei];
183 (vf[owner[facei]]*sigmaDeltaO + vf[neighbour[facei]]*sigmaDeltaN)
184 /(sigmaDeltaO + sigmaDeltaN);
188 auto& bfld = result.boundaryFieldRef();
206 sigma_.boundaryField()[patchi].patchNeighbourField()
214 const label own = pOwner[facei];
216 const scalar sigmaDeltaO =
219 const scalar sigmaDeltaN =
223 (vfO[facei]*sigmaDeltaO + vfN[facei]*sigmaDeltaN)
224 /(sigmaDeltaO + sigmaDeltaN);
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Interpolate the cell values to faces.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
~weightedFlux()
Destructor.
void clearOut()
Clear all fields.
#define forAll(list, i)
Loop across all elements in list.
Field< vector > vectorField
Specialisation of Field<T> for vector.
PtrList< fvPatch > fvPatchList
container classes for fvPatch
void deleteDemandDrivenData(DataPtr &dataPtr)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Mesh data needed to do the Finite Volume discretisation.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Weighted flux interpolation scheme class.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const polyBoundaryMesh & patches
makeSurfaceInterpolationScheme(cellCoBlended)
UList< label > labelUList
A UList of labels.
Graphite solid properties.
virtual bool coupled() const
Return true if this patch field is coupled.
const Boundary & boundaryField() const
Return const-reference to the boundary field.