Loading...
Searching...
No Matches
faceLimitedGrads.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "faceLimitedGrad.H"
30#include "gaussGrad.H"
31#include "fvMesh.H"
32#include "volMesh.H"
33#include "surfaceMesh.H"
34#include "volFields.H"
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39makeFvGradScheme(faceLimitedGrad)
40
41// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42
43template<>
46(
47 const volScalarField& vsf,
48 const word& name
49) const
50{
51 const fvMesh& mesh = vsf.mesh();
52
53 tmp<volVectorField> tGrad = basicGradScheme_().calcGrad(vsf, name);
54
55 if (k_ < SMALL)
56 {
57 return tGrad;
58 }
59
60 volVectorField& g = tGrad.ref();
61
62 const labelUList& owner = mesh.owner();
63 const labelUList& neighbour = mesh.neighbour();
64
65 const volVectorField& C = mesh.C();
66 const surfaceVectorField& Cf = mesh.Cf();
67
68 // create limiter
69 scalarField limiter(vsf.primitiveField().size(), 1.0);
70
71 const scalar rk = (1.0/k_ - 1.0);
72
73 forAll(owner, facei)
74 {
75 const label own = owner[facei];
76 const label nei = neighbour[facei];
77
78 const scalar vsfOwn = vsf[own];
79 const scalar vsfNei = vsf[nei];
80
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;
86
87 // owner side
88 limitFace
89 (
90 limiter[own],
91 maxFace - vsfOwn, minFace - vsfOwn,
92 (Cf[facei] - C[own]) & g[own]
93 );
94
95 // neighbour side
96 limitFace
97 (
98 limiter[nei],
99 maxFace - vsfNei, minFace - vsfNei,
100 (Cf[facei] - C[nei]) & g[nei]
101 );
102 }
103
104 const volScalarField::Boundary& bsf = vsf.boundaryField();
105
106 forAll(bsf, patchi)
107 {
108 const fvPatchScalarField& psf = bsf[patchi];
109
110 const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
111 const vectorField& pCf = Cf.boundaryField()[patchi];
112
113 if (psf.coupled())
114 {
115 const scalarField psfNei(psf.patchNeighbourField());
116
117 forAll(pOwner, pFacei)
118 {
119 const label own = pOwner[pFacei];
120
121 const scalar vsfOwn = vsf[own];
122 const scalar vsfNei = psfNei[pFacei];
123
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;
129
130 limitFace
131 (
132 limiter[own],
133 maxFace - vsfOwn, minFace - vsfOwn,
134 (pCf[pFacei] - C[own]) & g[own]
135 );
136 }
137 }
138 else if (psf.fixesValue())
139 {
140 forAll(pOwner, pFacei)
141 {
142 const label own = pOwner[pFacei];
143
144 const scalar vsfOwn = vsf[own];
145 const scalar vsfNei = psf[pFacei];
146
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;
152
153 limitFace
154 (
155 limiter[own],
156 maxFace - vsfOwn, minFace - vsfOwn,
157 (pCf[pFacei] - C[own]) & g[own]
158 );
159 }
160 }
161 }
162
163 if (fv::debug)
164 {
165 auto limits = gMinMax(limiter);
166 auto avg = gAverage(limiter);
167
168 Info<< "gradient limiter for: " << vsf.name()
169 << " min = " << limits.min()
170 << " max = " << limits.max()
171 << " average: " << avg << endl;
172 }
173
174 g.primitiveFieldRef() *= limiter;
175 g.correctBoundaryConditions();
176 gaussGrad<scalar>::correctBoundaryConditions(vsf, g);
178 return tGrad;
179}
180
181
182template<>
185(
186 const volVectorField& vvf,
187 const word& name
188) const
189{
190 const fvMesh& mesh = vvf.mesh();
191
192 tmp<volTensorField> tGrad = basicGradScheme_().calcGrad(vvf, name);
193
194 if (k_ < SMALL)
195 {
196 return tGrad;
197 }
198
199 volTensorField& g = tGrad.ref();
200
201 const labelUList& owner = mesh.owner();
202 const labelUList& neighbour = mesh.neighbour();
203
204 const volVectorField& C = mesh.C();
205 const surfaceVectorField& Cf = mesh.Cf();
206
207 // create limiter
209
210 const scalar rk = (1.0/k_ - 1.0);
211
212 forAll(owner, facei)
213 {
214 const label own = owner[facei];
215 const label nei = neighbour[facei];
216
217 const vector& vvfOwn = vvf[own];
218 const vector& vvfNei = vvf[nei];
219
220 // owner side
221 vector gradf((Cf[facei] - C[own]) & g[own]);
222
223 scalar vsfOwn = gradf & vvfOwn;
224 scalar vsfNei = gradf & vvfNei;
225
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;
231
232 limitFace
233 (
234 limiter[own],
235 maxFace - vsfOwn, minFace - vsfOwn,
236 magSqr(gradf)
237 );
238
239
240 // neighbour side
241 gradf = (Cf[facei] - C[nei]) & g[nei];
242
243 vsfOwn = gradf & vvfOwn;
244 vsfNei = gradf & vvfNei;
245
246 maxFace = max(vsfOwn, vsfNei);
247 minFace = min(vsfOwn, vsfNei);
248
249 limitFace
250 (
251 limiter[nei],
252 maxFace - vsfNei, minFace - vsfNei,
253 magSqr(gradf)
254 );
255 }
256
257
258 const volVectorField::Boundary& bvf = vvf.boundaryField();
259
260 forAll(bvf, patchi)
261 {
262 const fvPatchVectorField& psf = bvf[patchi];
263
264 const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
265 const vectorField& pCf = Cf.boundaryField()[patchi];
266
267 if (psf.coupled())
268 {
269 const vectorField psfNei(psf.patchNeighbourField());
270
271 forAll(pOwner, pFacei)
272 {
273 const label own = pOwner[pFacei];
274
275 const vector& vvfOwn = vvf[own];
276 const vector& vvfNei = psfNei[pFacei];
277
278 const vector gradf((pCf[pFacei] - C[own]) & g[own]);
279
280 const scalar vsfOwn = gradf & vvfOwn;
281 const scalar vsfNei = gradf & vvfNei;
282
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;
288
289 limitFace
290 (
291 limiter[own],
292 maxFace - vsfOwn, minFace - vsfOwn,
293 magSqr(gradf)
294 );
295 }
296 }
297 else if (psf.fixesValue())
298 {
299 forAll(pOwner, pFacei)
300 {
301 const label own = pOwner[pFacei];
302
303 const vector& vvfOwn = vvf[own];
304 const vector& vvfNei = psf[pFacei];
305
306 const vector gradf((pCf[pFacei] - C[own]) & g[own]);
307
308 const scalar vsfOwn = gradf & vvfOwn;
309 const scalar vsfNei = gradf & vvfNei;
310
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;
316
317 limitFace
318 (
319 limiter[own],
320 maxFace - vsfOwn, minFace - vsfOwn,
321 magSqr(gradf)
322 );
323 }
324 }
325 }
326
327 if (fv::debug)
328 {
329 auto limits = gMinMax(limiter);
330 auto avg = gAverage(limiter);
331
332 Info<< "gradient limiter for: " << vvf.name()
333 << " min = " << limits.min()
334 << " max = " << limits.max()
335 << " average: " << avg << endl;
336 }
337
338 g.primitiveFieldRef() *= limiter;
339 g.correctBoundaryConditions();
340 gaussGrad<vector>::correctBoundaryConditions(vvf, g);
341
342 return tGrad;
343}
344
345
346// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
const uniformDimensionedVectorField & g
Graphite solid properties.
Definition C.H:49
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.
Definition IOobjectI.H:205
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
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.
Definition gaussGrad.C:199
const fvMesh & mesh() const
Return const reference to mesh.
Definition gradScheme.H:138
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
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.
Definition word.H:66
auto limits
Definition setRDeltaT.H:186
dynamicFvMesh & mesh
auto & name
#define makeFvGradScheme(SS)
Definition gradScheme.H:238
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.
Definition hashSets.C:40
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.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
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.
Definition exprTraits.C:127
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
fvPatchField< vector > fvPatchVectorField
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition faNVDscheme.C:31
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299