47 const volScalarField& vsf,
51 const fvMesh&
mesh = vsf.mesh();
53 tmp<volVectorField> tGrad = basicGradScheme_().calcGrad(vsf,
name);
60 volVectorField&
g = tGrad.ref();
62 const labelUList& owner =
mesh.owner();
63 const labelUList& neighbour =
mesh.neighbour();
65 const volVectorField&
C =
mesh.C();
66 const surfaceVectorField& Cf =
mesh.Cf();
69 scalarField limiter(vsf.primitiveField().size(), 1.0);
71 const scalar rk = (1.0/k_ - 1.0);
75 const label own = owner[facei];
76 const label nei = neighbour[facei];
78 const scalar vsfOwn = vsf[own];
79 const scalar vsfNei = vsf[nei];
81 scalar maxFace = max(vsfOwn, vsfNei);
82 scalar minFace = min(vsfOwn, vsfNei);
83 const scalar maxMinFace = rk*(maxFace - minFace);
84 maxFace += maxMinFace;
85 minFace -= maxMinFace;
91 maxFace - vsfOwn, minFace - vsfOwn,
92 (Cf[facei] -
C[own]) &
g[own]
99 maxFace - vsfNei, minFace - vsfNei,
100 (Cf[facei] -
C[nei]) &
g[nei]
104 const volScalarField::Boundary& bsf = vsf.boundaryField();
108 const fvPatchScalarField& psf = bsf[patchi];
110 const labelUList& pOwner =
mesh.boundary()[patchi].faceCells();
111 const vectorField& pCf = Cf.boundaryField()[patchi];
115 const scalarField psfNei(psf.patchNeighbourField());
119 const label own = pOwner[pFacei];
121 const scalar vsfOwn = vsf[own];
122 const scalar vsfNei = psfNei[pFacei];
124 scalar maxFace = max(vsfOwn, vsfNei);
125 scalar minFace = min(vsfOwn, vsfNei);
126 const scalar maxMinFace = rk*(maxFace - minFace);
127 maxFace += maxMinFace;
128 minFace -= maxMinFace;
133 maxFace - vsfOwn, minFace - vsfOwn,
134 (pCf[pFacei] -
C[own]) &
g[own]
138 else if (psf.fixesValue())
142 const label own = pOwner[pFacei];
144 const scalar vsfOwn = vsf[own];
145 const scalar vsfNei = psf[pFacei];
147 scalar maxFace = max(vsfOwn, vsfNei);
148 scalar minFace = min(vsfOwn, vsfNei);
149 const scalar maxMinFace = rk*(maxFace - minFace);
150 maxFace += maxMinFace;
151 minFace -= maxMinFace;
156 maxFace - vsfOwn, minFace - vsfOwn,
157 (pCf[pFacei] -
C[own]) &
g[own]
165 auto limits = gMinMax(limiter);
166 auto avg = gAverage(limiter);
168 Info<<
"gradient limiter for: " << vsf.name()
169 <<
" min = " <<
limits.min()
170 <<
" max = " <<
limits.max()
171 <<
" average: " << avg << endl;
174 g.primitiveFieldRef() *= limiter;
175 g.correctBoundaryConditions();
176 gaussGrad<scalar>::correctBoundaryConditions(vsf,
g);
210 const scalar rk = (1.0/k_ - 1.0);
214 const label own = owner[facei];
215 const label nei = neighbour[facei];
217 const vector& vvfOwn = vvf[own];
218 const vector& vvfNei = vvf[nei];
221 vector gradf((Cf[facei] -
C[own]) &
g[own]);
223 scalar vsfOwn = gradf & vvfOwn;
224 scalar vsfNei = gradf & vvfNei;
226 scalar maxFace =
max(vsfOwn, vsfNei);
227 scalar minFace =
min(vsfOwn, vsfNei);
228 const scalar maxMinFace = rk*(maxFace - minFace);
229 maxFace += maxMinFace;
230 minFace -= maxMinFace;
235 maxFace - vsfOwn, minFace - vsfOwn,
241 gradf = (Cf[facei] -
C[nei]) &
g[nei];
243 vsfOwn = gradf & vvfOwn;
244 vsfNei = gradf & vvfNei;
246 maxFace =
max(vsfOwn, vsfNei);
247 minFace =
min(vsfOwn, vsfNei);
252 maxFace - vsfNei, minFace - vsfNei,
262 const fvPatchVectorField& psf = bvf[patchi];
264 const labelUList& pOwner =
mesh.boundary()[patchi].faceCells();
269 const vectorField psfNei(psf.patchNeighbourField());
273 const label own = pOwner[pFacei];
275 const vector& vvfOwn = vvf[own];
276 const vector& vvfNei = psfNei[pFacei];
278 const vector gradf((pCf[pFacei] -
C[own]) &
g[own]);
280 const scalar vsfOwn = gradf & vvfOwn;
281 const scalar vsfNei = gradf & vvfNei;
283 scalar maxFace = max(vsfOwn, vsfNei);
284 scalar minFace = min(vsfOwn, vsfNei);
285 const scalar maxMinFace = rk*(maxFace - minFace);
286 maxFace += maxMinFace;
287 minFace -= maxMinFace;
292 maxFace - vsfOwn, minFace - vsfOwn,
297 else if (psf.fixesValue())
301 const label own = pOwner[pFacei];
303 const vector& vvfOwn = vvf[own];
304 const vector& vvfNei = psf[pFacei];
306 const vector gradf((pCf[pFacei] -
C[own]) &
g[own]);
308 const scalar vsfOwn = gradf & vvfOwn;
309 const scalar vsfNei = gradf & vvfNei;
311 scalar maxFace =
max(vsfOwn, vsfNei);
312 scalar minFace =
min(vsfOwn, vsfNei);
313 const scalar maxMinFace = rk*(maxFace - minFace);
314 maxFace += maxMinFace;
315 minFace -= maxMinFace;
320 maxFace - vsfOwn, minFace - vsfOwn,
332 Info<<
"gradient limiter for: " << vvf.
name()
333 <<
" min = " <<
limits.min()
334 <<
" max = " <<
limits.max()
335 <<
" average: " << avg <<
endl;
339 g.correctBoundaryConditions();
340 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 Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
const word & name() const noexcept
Return the object name.
void size(const label n)
Older name for setAddressableSize.
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.
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
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
GeometricField< vector, fvPatchField, volMesh > volVectorField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
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.
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
UList< label > labelUList
A UList of labels.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
fvPatchField< vector > fvPatchVectorField
tmp< areaScalarField > limiter(const areaScalarField &phi)
#define forAll(list, i)
Loop across all elements in list.