faceLimitedFaGrads.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) 2016-2017 Wikki Ltd
9  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "faceLimitedFaGrad.H"
30 #include "gaussFaGrad.H"
31 #include "faMesh.H"
32 #include "areaFaMesh.H"
33 #include "edgeFaMesh.H"
34 #include "areaFields.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 makeFaGradScheme(faceLimitedGrad)
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace fa
49 {
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 template<>
55 (
56  scalar& limiter,
57  const scalar& maxDelta,
58  const scalar& minDelta,
59  const scalar& extrapolate
60 )
61 {
62  if (extrapolate > maxDelta + VSMALL)
63  {
64  limiter = min(limiter, maxDelta/extrapolate);
65  }
66  else if (extrapolate < minDelta - VSMALL)
67  {
68  limiter = min(limiter, minDelta/extrapolate);
69  }
70 }
71 
72 
73 template<class Type>
75 (
76  Type& limiter,
77  const Type& maxDelta,
78  const Type& minDelta,
79  const Type& extrapolate
80 )
81 {
82  for (direction cmpt = 0; cmpt < Type::nComponents; ++cmpt)
83  {
85  (
86  limiter.component(cmpt),
87  maxDelta.component(cmpt),
88  minDelta.component(cmpt),
89  extrapolate.component(cmpt)
90  );
91  }
92 }
93 
94 
95 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
96 
97 template<>
98 tmp<areaVectorField> faceLimitedGrad<scalar>::calcGrad
99 (
100  const areaScalarField& vsf,
101  const word& name
102 ) const
103 {
104  const faMesh& mesh = vsf.mesh();
105 
106  tmp<areaVectorField> tGrad = basicGradScheme_().calcGrad(vsf, name);
107 
108  if (k_ < SMALL)
109  {
110  return tGrad;
111  }
112 
113  areaVectorField& g = tGrad.ref();
114 
115  const labelUList& owner = mesh.owner();
116  const labelUList& neighbour = mesh.neighbour();
117 
118  const areaVectorField& C = mesh.areaCentres();
119  const edgeVectorField& Cf = mesh.edgeCentres();
120 
121  scalarField maxVsf(vsf.internalField());
122  scalarField minVsf(vsf.internalField());
123 
124  forAll(owner, facei)
125  {
126  label own = owner[facei];
127  label nei = neighbour[facei];
128 
129  scalar vsfOwn = vsf[own];
130  scalar vsfNei = vsf[nei];
131 
132  maxVsf[own] = max(maxVsf[own], vsfNei);
133  minVsf[own] = min(minVsf[own], vsfNei);
134 
135  maxVsf[nei] = max(maxVsf[nei], vsfOwn);
136  minVsf[nei] = min(minVsf[nei], vsfOwn);
137  }
138 
139 
140  const areaScalarField::Boundary& bsf = vsf.boundaryField();
141 
142  forAll(bsf, patchi)
143  {
144  const faPatchScalarField& psf = bsf[patchi];
145 
146  const labelUList& pOwner = mesh.boundary()[patchi].edgeFaces();
147 
148  if (psf.coupled())
149  {
150  scalarField psfNei(psf.patchNeighbourField());
151 
152  forAll(pOwner, pFacei)
153  {
154  label own = pOwner[pFacei];
155  scalar vsfNei = psfNei[pFacei];
156 
157  maxVsf[own] = max(maxVsf[own], vsfNei);
158  minVsf[own] = min(minVsf[own], vsfNei);
159  }
160  }
161  else
162  {
163  forAll(pOwner, pFacei)
164  {
165  label own = pOwner[pFacei];
166  scalar vsfNei = psf[pFacei];
167 
168  maxVsf[own] = max(maxVsf[own], vsfNei);
169  minVsf[own] = min(minVsf[own], vsfNei);
170  }
171  }
172  }
173 
174  maxVsf -= vsf;
175  minVsf -= vsf;
176 
177  if (k_ < 1.0)
178  {
179  scalarField maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
180  maxVsf += maxMinVsf;
181  minVsf -= maxMinVsf;
182  }
183 
184 
185  // create limiter
186  scalarField limiter(vsf.internalField().size(), 1.0);
187 
188  forAll(owner, facei)
189  {
190  label own = owner[facei];
191  label nei = neighbour[facei];
192 
193  // owner side
194  limitEdge
195  (
196  limiter[own],
197  maxVsf[own],
198  minVsf[own],
199  (Cf[facei] - C[own]) & g[own]
200  );
201 
202  // neighbour side
203  limitEdge
204  (
205  limiter[nei],
206  maxVsf[nei],
207  minVsf[nei],
208  (Cf[facei] - C[nei]) & g[nei]
209  );
210  }
211 
212  forAll(bsf, patchi)
213  {
214  const labelUList& pOwner = mesh.boundary()[patchi].edgeFaces();
215  const vectorField& pCf = Cf.boundaryField()[patchi];
216 
217  forAll(pOwner, pFacei)
218  {
219  label own = pOwner[pFacei];
220 
221  limitEdge
222  (
223  limiter[own],
224  maxVsf[own],
225  minVsf[own],
226  (pCf[pFacei] - C[own]) & g[own]
227  );
228  }
229  }
230 
231  if (fa::debug)
232  {
233  Info<< "gradient limiter for: " << vsf.name()
234  << " max = " << gMax(limiter)
235  << " min = " << gMin(limiter)
236  << " average: " << gAverage(limiter) << endl;
237  }
238 
239  g.primitiveFieldRef() *= limiter;
240  g.correctBoundaryConditions();
242 
243  return tGrad;
244 }
245 
246 
247 template<>
248 tmp<areaTensorField> faceLimitedGrad<vector>::calcGrad
249 (
250  const areaVectorField& vsf,
251  const word& name
252 ) const
253 {
254  const faMesh& mesh = vsf.mesh();
255 
256  tmp<areaTensorField> tGrad = basicGradScheme_().grad(vsf, name);
257 
258  if (k_ < SMALL)
259  {
260  return tGrad;
261  }
262 
263  areaTensorField& g = tGrad.ref();
264 
265  const labelUList& owner = mesh.owner();
266  const labelUList& neighbour = mesh.neighbour();
267 
268  const areaVectorField& C = mesh.areaCentres();
269  const edgeVectorField& Cf = mesh.edgeCentres();
270 
271  vectorField maxVsf(vsf.internalField());
272  vectorField minVsf(vsf.internalField());
273 
274  forAll(owner, facei)
275  {
276  label own = owner[facei];
277  label nei = neighbour[facei];
278 
279  const vector& vsfOwn = vsf[own];
280  const vector& vsfNei = vsf[nei];
281 
282  maxVsf[own] = max(maxVsf[own], vsfNei);
283  minVsf[own] = min(minVsf[own], vsfNei);
284 
285  maxVsf[nei] = max(maxVsf[nei], vsfOwn);
286  minVsf[nei] = min(minVsf[nei], vsfOwn);
287  }
288 
289 
290  const areaVectorField::Boundary& bsf = vsf.boundaryField();
291 
292  forAll(bsf, patchi)
293  {
294  const faPatchVectorField& psf = bsf[patchi];
295  const labelUList& pOwner = mesh.boundary()[patchi].edgeFaces();
296 
297  if (psf.coupled())
298  {
299  vectorField psfNei(psf.patchNeighbourField());
300 
301  forAll(pOwner, pFacei)
302  {
303  label own = pOwner[pFacei];
304  const vector& vsfNei = psfNei[pFacei];
305 
306  maxVsf[own] = max(maxVsf[own], vsfNei);
307  minVsf[own] = min(minVsf[own], vsfNei);
308  }
309  }
310  else
311  {
312  forAll(pOwner, pFacei)
313  {
314  label own = pOwner[pFacei];
315  const vector& vsfNei = psf[pFacei];
316 
317  maxVsf[own] = max(maxVsf[own], vsfNei);
318  minVsf[own] = min(minVsf[own], vsfNei);
319  }
320  }
321  }
322 
323  maxVsf -= vsf;
324  minVsf -= vsf;
325 
326  if (k_ < 1.0)
327  {
328  vectorField maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
329  maxVsf += maxMinVsf;
330  minVsf -= maxMinVsf;
331 
332  //maxVsf *= 1.0/k_;
333  //minVsf *= 1.0/k_;
334  }
335 
336 
337  // create limiter
338  vectorField limiter(vsf.internalField().size(), vector::one);
339 
340  forAll(owner, facei)
341  {
342  label own = owner[facei];
343  label nei = neighbour[facei];
344 
345  // owner side
346  limitEdge
347  (
348  limiter[own],
349  maxVsf[own],
350  minVsf[own],
351  (Cf[facei] - C[own]) & g[own]
352  );
353 
354  // neighbour side
355  limitEdge
356  (
357  limiter[nei],
358  maxVsf[nei],
359  minVsf[nei],
360  (Cf[facei] - C[nei]) & g[nei]
361  );
362  }
363 
364  forAll(bsf, patchi)
365  {
366  const labelUList& pOwner = mesh.boundary()[patchi].edgeFaces();
367  const vectorField& pCf = Cf.boundaryField()[patchi];
368 
369  forAll(pOwner, pFacei)
370  {
371  label own = pOwner[pFacei];
372 
373  limitEdge
374  (
375  limiter[own],
376  maxVsf[own],
377  minVsf[own],
378  ((pCf[pFacei] - C[own]) & g[own])
379  );
380  }
381  }
382 
383  if (fa::debug)
384  {
385  Info<< "gradient limiter for: " << vsf.name()
386  << " max = " << gMax(limiter)
387  << " min = " << gMin(limiter)
388  << " average: " << gAverage(limiter) << endl;
389  }
390 
391  tensorField& gIf = g.primitiveFieldRef();
392 
393  forAll(gIf, celli)
394  {
395  gIf[celli] = tensor
396  (
397  cmptMultiply(limiter[celli], gIf[celli].x()),
398  cmptMultiply(limiter[celli], gIf[celli].y()),
399  cmptMultiply(limiter[celli], gIf[celli].z())
400  );
401  }
402 
403  g.correctBoundaryConditions();
405 
406  return tGrad;
407 }
408 
409 
410 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
411 
412 } // End namespace fa
413 
414 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
415 
416 } // End namespace Foam
417 
418 // ************************************************************************* //
faPatchField< scalar > faPatchScalarField
uint8_t direction
Definition: direction.H:46
virtual tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh > > calcGrad(const GeometricField< Type, faPatchField, areaMesh > &vsf, const word &name) const
Return the gradient of the given field to the gradScheme::grad for optional caching.
Type gMin(const FieldField< Field, Type > &f)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Tensor< scalar > tensor
Definition: symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
GeometricBoundaryField< scalar, faPatchField, areaMesh > Boundary
Type of boundary fields.
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:31
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
scalar y
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
GeometricField< vector, faePatchField, edgeMesh > edgeVectorField
Definition: edgeFieldsFwd.H:58
Vector< scalar > vector
Definition: vector.H:57
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
GeometricField< tensor, faPatchField, areaMesh > areaTensorField
Definition: areaFieldsFwd.H:85
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
const uniformDimensionedVectorField & g
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
volScalarField & C
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
Type gAverage(const FieldField< Field, Type > &f)
makeFaGradScheme(faceLimitedGrad)
static void correctBoundaryConditions(const GeometricField< Type, faPatchField, areaMesh > &, GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh > &)
Correct the boundary values of the gradient using the patchField.
Definition: gaussFaGrad.C:110
faPatchField< vector > faPatchVectorField
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< vector > vectorField
Specialisation of Field<T> for vector.
GeometricField< vector, faPatchField, areaMesh > areaVectorField
Definition: areaFieldsFwd.H:81
static void limitEdge(Type &limiter, const Type &maxDelta, const Type &minDelta, const Type &extrapolate)
Namespace for OpenFOAM.
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Definition: areaFieldsFwd.H:80