48 const volScalarField& vsf,
52 const fvMesh&
mesh = vsf.mesh();
54 tmp<volVectorField> tGrad = basicGradScheme_().calcGrad(vsf,
name);
61 volVectorField&
g = tGrad.ref();
63 const labelUList& owner =
mesh.owner();
64 const labelUList& neighbour =
mesh.neighbour();
66 const volVectorField&
C =
mesh.C();
67 const surfaceVectorField& Cf =
mesh.Cf();
69 const scalar rk = (1.0/k_ - 1.0);
73 const label own = owner[facei];
74 const label nei = neighbour[facei];
76 const scalar vsfOwn = vsf[own];
77 const scalar vsfNei = vsf[nei];
79 scalar maxFace = max(vsfOwn, vsfNei);
80 scalar minFace = min(vsfOwn, vsfNei);
84 const scalar maxMinFace = rk*(maxFace - minFace);
85 maxFace += maxMinFace;
86 minFace -= maxMinFace;
90 cellMDLimitedGrad<scalar>::limitFace
99 cellMDLimitedGrad<scalar>::limitFace
108 const volScalarField::Boundary& bsf = vsf.boundaryField();
112 const fvPatchScalarField& psf = bsf[patchi];
114 const labelUList& pOwner =
mesh.boundary()[patchi].faceCells();
115 const vectorField& pCf = Cf.boundaryField()[patchi];
119 const scalarField psfNei(psf.patchNeighbourField());
123 const label own = pOwner[pFacei];
125 scalar vsfOwn = vsf[own];
126 scalar vsfNei = psfNei[pFacei];
128 scalar maxFace = max(vsfOwn, vsfNei);
129 scalar minFace = min(vsfOwn, vsfNei);
133 const scalar maxMinFace = rk*(maxFace - minFace);
134 maxFace += maxMinFace;
135 minFace -= maxMinFace;
138 cellMDLimitedGrad<scalar>::limitFace
147 else if (psf.fixesValue())
151 const label own = pOwner[pFacei];
153 const scalar vsfOwn = vsf[own];
154 const scalar vsfNei = psf[pFacei];
156 scalar maxFace = max(vsfOwn, vsfNei);
157 scalar minFace = min(vsfOwn, vsfNei);
161 const scalar maxMinFace = rk*(maxFace - minFace);
162 maxFace += maxMinFace;
163 minFace -= maxMinFace;
166 cellMDLimitedGrad<scalar>::limitFace
177 g.correctBoundaryConditions();
178 gaussGrad<scalar>::correctBoundaryConditions(vsf,
g);
209 const scalar rk = (1.0/k_ - 1.0);
213 const label own = owner[facei];
214 const label nei = neighbour[facei];
216 const vector& vvfOwn = vvf[own];
217 const vector& vvfNei = vvf[nei];
224 const vector maxMinFace(rk*(maxFace - minFace));
225 maxFace += maxMinFace;
226 minFace -= maxMinFace;
230 cellMDLimitedGrad<vector>::limitFace
240 cellMDLimitedGrad<vector>::limitFace
261 const vectorField psfNei(psf.patchNeighbourField());
265 const label own = pOwner[pFacei];
267 const vector& vvfOwn = vvf[own];
268 const vector& vvfNei = psfNei[pFacei];
275 const vector maxMinFace(rk*(maxFace - minFace));
276 maxFace += maxMinFace;
277 minFace -= maxMinFace;
280 cellMDLimitedGrad<vector>::limitFace
283 maxFace - vvfOwn, minFace - vvfOwn,
288 else if (psf.fixesValue())
292 const label own = pOwner[pFacei];
294 const vector& vvfOwn = vvf[own];
295 const vector& vvfNei = psf[pFacei];
302 const vector maxMinFace(rk*(maxFace - minFace));
303 maxFace += maxMinFace;
304 minFace -= maxMinFace;
307 cellMDLimitedGrad<vector>::limitFace
318 g.correctBoundaryConditions();
319 gaussGrad<vector>::correctBoundaryConditions(vvf,
g);
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
const uniformDimensionedVectorField & g
Graphite solid properties.
const Mesh & mesh() const noexcept
Return const reference to mesh.
GeometricBoundaryField< vector, fvPatchField, volMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Mesh data needed to do the Finite Volume discretisation.
virtual bool fixesValue() const
True if the patch field fixes a value.
virtual bool coupled() const
True if the patch field is coupled.
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
static void limitFace(typename outerProduct< vector, Type >::type &g, const Type &maxDelta, const Type &minDelta, const vector &dcf)
virtual tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > calcGrad(const GeometricField< Type, fvPatchField, volMesh > &vsf, const word &name) const
Return the gradient of the given field to the gradScheme::grad for optional caching.
static void correctBoundaryConditions(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > &)
Correct the boundary values of the gradient using the patchField snGrad functions.
const fvMesh & mesh() const
Return const reference to mesh.
A class for managing temporary objects.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
#define makeFvGradScheme(SS)
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< tensor, fvPatchField, volMesh > volTensorField
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Field< vector > vectorField
Specialisation of Field<T> for vector.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
UList< label > labelUList
A UList of labels.
fvPatchField< vector > fvPatchVectorField
#define forAll(list, i)
Loop across all elements in list.