edgeLimitedFaGrads.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 "edgeLimitedFaGrad.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(edgeLimitedGrad)
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace fa
49 {
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 template<class Type>
54 inline void edgeLimitedGrad<Type>::limitEdge
55 (
56  scalar& limiter,
57  const scalar maxDelta,
58  const scalar minDelta,
59  const scalar extrapolate
60 ) const
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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74 
75 template<>
76 tmp<areaVectorField> edgeLimitedGrad<scalar>::calcGrad
77 (
78  const areaScalarField& vsf,
79  const word& name
80 ) const
81 {
82  const faMesh& mesh = vsf.mesh();
83 
84  tmp<areaVectorField> tGrad = basicGradScheme_().calcGrad(vsf, name);
85 
86  if (k_ < SMALL)
87  {
88  return tGrad;
89  }
90 
91  areaVectorField& g = tGrad.ref();
92 
93  const labelUList& owner = mesh.owner();
94  const labelUList& neighbour = mesh.neighbour();
95 
96  const areaVectorField& C = mesh.areaCentres();
97  const edgeVectorField& Cf = mesh.edgeCentres();
98 
99  // create limiter
100  scalarField limiter(vsf.internalField().size(), 1.0);
101 
102  scalar rk = (1.0/k_ - 1.0);
103 
104  forAll(owner, edgei)
105  {
106  label own = owner[edgei];
107  label nei = neighbour[edgei];
108 
109  scalar vsfOwn = vsf[own];
110  scalar vsfNei = vsf[nei];
111 
112  scalar maxEdge = max(vsfOwn, vsfNei);
113  scalar minEdge = min(vsfOwn, vsfNei);
114  scalar maxMinEdge = rk*(maxEdge - minEdge);
115  maxEdge += maxMinEdge;
116  minEdge -= maxMinEdge;
117 
118  // owner side
119  limitEdge
120  (
121  limiter[own],
122  maxEdge - vsfOwn, minEdge - vsfOwn,
123  (Cf[edgei] - C[own]) & g[own]
124  );
125 
126  // neighbour side
127  limitEdge
128  (
129  limiter[nei],
130  maxEdge - vsfNei, minEdge - vsfNei,
131  (Cf[edgei] - C[nei]) & g[nei]
132  );
133  }
134 
135  const areaScalarField::Boundary& bsf = vsf.boundaryField();
136 
137  forAll(bsf, patchi)
138  {
139  const faPatchScalarField& psf = bsf[patchi];
140 
141  const labelUList& pOwner = mesh.boundary()[patchi].edgeFaces();
142  const vectorField& pCf = Cf.boundaryField()[patchi];
143 
144  if (psf.coupled())
145  {
146  const scalarField psfNei(psf.patchNeighbourField());
147 
148  forAll(pOwner, pEdgei)
149  {
150  label own = pOwner[pEdgei];
151 
152  scalar vsfOwn = vsf[own];
153  scalar vsfNei = psfNei[pEdgei];
154 
155  scalar maxEdge = max(vsfOwn, vsfNei);
156  scalar minEdge = min(vsfOwn, vsfNei);
157  scalar maxMinEdge = rk*(maxEdge - minEdge);
158  maxEdge += maxMinEdge;
159  minEdge -= maxMinEdge;
160 
161  limitEdge
162  (
163  limiter[own],
164  maxEdge - vsfOwn, minEdge - vsfOwn,
165  (pCf[pEdgei] - C[own]) & g[own]
166  );
167  }
168  }
169  else if (psf.fixesValue())
170  {
171  forAll(pOwner, pEdgei)
172  {
173  label own = pOwner[pEdgei];
174 
175  scalar vsfOwn = vsf[own];
176  scalar vsfNei = psf[pEdgei];
177 
178  scalar maxEdge = max(vsfOwn, vsfNei);
179  scalar minEdge = min(vsfOwn, vsfNei);
180  scalar maxMinEdge = rk*(maxEdge - minEdge);
181  maxEdge += maxMinEdge;
182  minEdge -= maxMinEdge;
183 
184  limitEdge
185  (
186  limiter[own],
187  maxEdge - vsfOwn, minEdge - vsfOwn,
188  (pCf[pEdgei] - C[own]) & g[own]
189  );
190  }
191  }
192  }
193 
194  if (fa::debug)
195  {
196  Info<< "gradient limiter for: " << vsf.name()
197  << " max = " << gMax(limiter)
198  << " min = " << gMin(limiter)
199  << " average: " << gAverage(limiter) << endl;
200  }
201 
202  g.primitiveFieldRef() *= limiter;
203  g.correctBoundaryConditions();
205 
206  return tGrad;
207 }
208 
209 
210 template<>
211 tmp<areaTensorField> edgeLimitedGrad<vector>::calcGrad
212 (
213  const areaVectorField& vvf,
214  const word& name
215 ) const
216 {
217  const faMesh& mesh = vvf.mesh();
218 
219  tmp<areaTensorField> tGrad = basicGradScheme_().grad(vvf, name);
220 
221  if (k_ < SMALL)
222  {
223  return tGrad;
224  }
225 
226  areaTensorField& g = tGrad.ref();
227 
228  const labelUList& owner = mesh.owner();
229  const labelUList& neighbour = mesh.neighbour();
230 
231  const areaVectorField& C = mesh.areaCentres();
232  const edgeVectorField& Cf = mesh.edgeCentres();
233 
234  // create limiter
235  scalarField limiter(vvf.internalField().size(), 1.0);
236 
237  scalar rk = (1.0/k_ - 1.0);
238 
239  forAll(owner, edgei)
240  {
241  label own = owner[edgei];
242  label nei = neighbour[edgei];
243 
244  vector vvfOwn = vvf[own];
245  vector vvfNei = vvf[nei];
246 
247  // owner side
248  vector gradf = (Cf[edgei] - C[own]) & g[own];
249 
250  scalar vsfOwn = gradf & vvfOwn;
251  scalar vsfNei = gradf & vvfNei;
252 
253  scalar maxEdge = max(vsfOwn, vsfNei);
254  scalar minEdge = min(vsfOwn, vsfNei);
255  scalar maxMinEdge = rk*(maxEdge - minEdge);
256  maxEdge += maxMinEdge;
257  minEdge -= maxMinEdge;
258 
259  limitEdge
260  (
261  limiter[own],
262  maxEdge - vsfOwn, minEdge - vsfOwn,
263  magSqr(gradf)
264  );
265 
266 
267  // neighbour side
268  gradf = (Cf[edgei] - C[nei]) & g[nei];
269 
270  vsfOwn = gradf & vvfOwn;
271  vsfNei = gradf & vvfNei;
272 
273  maxEdge = max(vsfOwn, vsfNei);
274  minEdge = min(vsfOwn, vsfNei);
275 
276  limitEdge
277  (
278  limiter[nei],
279  maxEdge - vsfNei, minEdge - vsfNei,
280  magSqr(gradf)
281  );
282  }
283 
284 
285  const areaVectorField::Boundary& bvf = vvf.boundaryField();
286 
287  forAll(bvf, patchi)
288  {
289  const faPatchVectorField& psf = bvf[patchi];
290 
291  const labelUList& pOwner = mesh.boundary()[patchi].edgeFaces();
292  const vectorField& pCf = Cf.boundaryField()[patchi];
293 
294  if (psf.coupled())
295  {
296  const vectorField psfNei(psf.patchNeighbourField());
297 
298  forAll(pOwner, pEdgei)
299  {
300  label own = pOwner[pEdgei];
301 
302  vector vvfOwn = vvf[own];
303  vector vvfNei = psfNei[pEdgei];
304 
305  vector gradf = (pCf[pEdgei] - C[own]) & g[own];
306 
307  scalar vsfOwn = gradf & vvfOwn;
308  scalar vsfNei = gradf & vvfNei;
309 
310  scalar maxEdge = max(vsfOwn, vsfNei);
311  scalar minEdge = min(vsfOwn, vsfNei);
312  scalar maxMinEdge = rk*(maxEdge - minEdge);
313  maxEdge += maxMinEdge;
314  minEdge -= maxMinEdge;
315 
316  limitEdge
317  (
318  limiter[own],
319  maxEdge - vsfOwn, minEdge - vsfOwn,
320  magSqr(gradf)
321  );
322  }
323  }
324  else if (psf.fixesValue())
325  {
326  forAll(pOwner, pEdgei)
327  {
328  label own = pOwner[pEdgei];
329 
330  vector vvfOwn = vvf[own];
331  vector vvfNei = psf[pEdgei];
332 
333  vector gradf = (pCf[pEdgei] - C[own]) & g[own];
334 
335  scalar vsfOwn = gradf & vvfOwn;
336  scalar vsfNei = gradf & vvfNei;
337 
338  scalar maxEdge = max(vsfOwn, vsfNei);
339  scalar minEdge = min(vsfOwn, vsfNei);
340  scalar maxMinEdge = rk*(maxEdge - minEdge);
341  maxEdge += maxMinEdge;
342  minEdge -= maxMinEdge;
343 
344  limitEdge
345  (
346  limiter[own],
347  maxEdge - vsfOwn, minEdge - vsfOwn,
348  magSqr(gradf)
349  );
350  }
351  }
352  }
353 
354  if (fa::debug)
355  {
356  Info<< "gradient limiter for: " << vvf.name()
357  << " max = " << gMax(limiter)
358  << " min = " << gMin(limiter)
359  << " average: " << gAverage(limiter) << endl;
360  }
361 
362  g.primitiveFieldRef() *= limiter;
363  g.correctBoundaryConditions();
365 
366  return tGrad;
367 }
368 
369 
370 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
371 
372 } // End namespace fa
373 
374 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
375 
376 } // End namespace Foam
377 
378 // ************************************************************************* //
faPatchField< scalar > faPatchScalarField
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
GeometricBoundaryField< scalar, faPatchField, areaMesh > Boundary
Type of boundary fields.
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:31
makeFaGradScheme(edgeLimitedGrad)
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
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
const uniformDimensionedVectorField & g
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
volScalarField & C
Type gAverage(const FieldField< Field, Type > &f)
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.
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.
GeometricField< vector, faPatchField, areaMesh > areaVectorField
Definition: areaFieldsFwd.H:81
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Definition: areaFieldsFwd.H:80