35 oDelta_.reset(
nullptr);
36 nDelta_.reset(
nullptr);
52void Foam::weightedFlux<Type>::makeDeltas()
const
56 oDelta_ = std::make_unique<surfaceScalarField>
61 mesh.pointsInstance(),
67 auto& oDelta = *oDelta_;
69 nDelta_ = std::make_unique<surfaceScalarField>
74 mesh.pointsInstance(),
80 auto& nDelta = *nDelta_;
82 const labelUList& owner =
mesh.owner();
83 const labelUList& neighbour =
mesh.neighbour();
85 const surfaceVectorField
n(
mesh.Sf()/
mesh.magSf());
87 const vectorField&
C =
mesh.cellCentres();
88 const vectorField& Cf =
mesh.faceCentres();
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]));
154 const GeometricField<Type, fvPatchField, volMesh>& vf
157 const fvMesh&
mesh = vf.mesh();
161 auto tresult = tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>::New
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();
192 fvsPatchField<Type>& pfld = bfld[patchi];
197 pfld = vf.boundaryField()[patchi];
206 sigma_.boundaryField()[patchi].patchNeighbourField()
209 Field<Type> vfO(vf.boundaryField()[patchi].patchInternalField());
210 Field<Type> vfN(vf.boundaryField()[patchi].patchNeighbourField());
214 const label own = pOwner[facei];
216 const scalar sigmaDeltaO =
217 sigma_[own]/oDelta.boundaryField()[patchi][facei];
219 const scalar sigmaDeltaN =
220 sigmaN[facei]/nDelta.boundaryField()[patchi][facei];
223 (vfO[facei]*sigmaDeltaO + vfN[facei]*sigmaDeltaN)
224 /(sigmaDeltaO + sigmaDeltaN);
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
Graphite solid properties.
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Generic GeometricField class.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & name() const noexcept
Return the object name.
Mesh data needed to do the Finite Volume discretisation.
virtual bool coupled() const
True if the patch field is coupled.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
const fvMesh & mesh() const
Return mesh reference.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
A class for managing temporary objects.
Weighted flux interpolation scheme class.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Interpolate the cell values to faces.
virtual ~weightedFlux()
Destructor.
const surfaceScalarField & oDelta() const
Return the distance between face and owner cell.
const surfaceScalarField & nDelta() const
Return the distance between face and neighbour cell.
void clearOut()
Clear all fields.
const polyBoundaryMesh & patches
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
UList< label > labelUList
A UList of labels.
#define forAll(list, i)
Loop across all elements in list.
#define makeSurfaceInterpolationScheme(SS)