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 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 makeFaGradScheme(edgeLimitedGrad)
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace fa
44 {
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 template<class Type>
49 inline void edgeLimitedGrad<Type>::limitEdge
50 (
51  scalar& limiter,
52  const scalar maxDelta,
53  const scalar minDelta,
54  const scalar extrapolate
55 ) const
56 {
57  if (extrapolate > maxDelta + VSMALL)
58  {
59  limiter = min(limiter, maxDelta/extrapolate);
60  }
61  else if (extrapolate < minDelta - VSMALL)
62  {
63  limiter = min(limiter, minDelta/extrapolate);
64  }
65 }
66 
67 
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70 template<>
71 tmp<areaVectorField> edgeLimitedGrad<scalar>::calcGrad
72 (
73  const areaScalarField& vsf,
74  const word& name
75 ) const
76 {
77  const faMesh& mesh = vsf.mesh();
78 
79  tmp<areaVectorField> tGrad = basicGradScheme_().calcGrad(vsf, name);
80 
81  if (k_ < SMALL)
82  {
83  return tGrad;
84  }
85 
86  areaVectorField& g = tGrad.ref();
87 
88  const labelUList& owner = mesh.owner();
89  const labelUList& neighbour = mesh.neighbour();
90 
91  const areaVectorField& C = mesh.areaCentres();
92  const edgeVectorField& Cf = mesh.edgeCentres();
93 
94  // Create limiter field
95  scalarField limiter(vsf.internalField().size(), 1.0);
96 
97  const scalar rk = (1.0/k_ - 1.0);
98 
99  forAll(owner, edgei)
100  {
101  const label own = owner[edgei];
102  const label nei = neighbour[edgei];
103 
104  const scalar vsfOwn = vsf[own];
105  const scalar vsfNei = vsf[nei];
106 
107  scalar maxEdge = max(vsfOwn, vsfNei);
108  scalar minEdge = min(vsfOwn, vsfNei);
109  const scalar maxMinEdge = rk*(maxEdge - minEdge);
110  maxEdge += maxMinEdge;
111  minEdge -= maxMinEdge;
112 
113  // owner side
114  limitEdge
115  (
116  limiter[own],
117  maxEdge - vsfOwn,
118  minEdge - vsfOwn,
119  (Cf[edgei] - C[own]) & g[own]
120  );
121 
122  // neighbour side
123  limitEdge
124  (
125  limiter[nei],
126  maxEdge - vsfNei,
127  minEdge - vsfNei,
128  (Cf[edgei] - C[nei]) & g[nei]
129  );
130  }
131 
132  // Lambda expression to update limiter for boundary edges
133  auto updateLimiter = [&](const label patchi, const scalarField& fld) -> void
134  {
135  const labelUList& pOwner = mesh.boundary()[patchi].edgeFaces();
136  const vectorField& pCf = Cf.boundaryField()[patchi];
137 
138  forAll(pOwner, edgei)
139  {
140  const label own = pOwner[edgei];
141 
142  const scalar vsfOwn = vsf[own];
143  const scalar vsfNei = fld[edgei];
144 
145  scalar maxEdge = max(vsfOwn, vsfNei);
146  scalar minEdge = min(vsfOwn, vsfNei);
147  const scalar maxMinEdge = rk*(maxEdge - minEdge);
148  maxEdge += maxMinEdge;
149  minEdge -= maxMinEdge;
150 
151  limitEdge
152  (
153  limiter[own],
154  maxEdge - vsfOwn,
155  minEdge - vsfOwn,
156  (pCf[edgei] - C[own]) & g[own]
157  );
158  }
159  };
160 
161  const areaScalarField::Boundary& bsf = vsf.boundaryField();
162  forAll(bsf, patchi)
163  {
164  const faPatchScalarField& psf = bsf[patchi];
165 
166  if (psf.coupled())
167  {
168  updateLimiter(patchi, psf.patchNeighbourField());
169  }
170  else if (psf.fixesValue())
171  {
172  updateLimiter(patchi, psf);
173  }
174  }
175 
176  if (fa::debug)
177  {
178  Info<< "gradient limiter for: " << vsf.name()
179  << " max = " << gMax(limiter)
180  << " min = " << gMin(limiter)
181  << " average: " << gAverage(limiter) << endl;
182  }
183 
184  g.primitiveFieldRef() *= limiter;
185  g.correctBoundaryConditions();
187 
188  return tGrad;
189 }
190 
191 
192 template<>
193 tmp<areaTensorField> edgeLimitedGrad<vector>::calcGrad
194 (
195  const areaVectorField& vvf,
196  const word& name
197 ) const
198 {
199  const faMesh& mesh = vvf.mesh();
200 
201  tmp<areaTensorField> tGrad = basicGradScheme_().grad(vvf, name);
202 
203  if (k_ < SMALL)
204  {
205  return tGrad;
206  }
207 
208  areaTensorField& g = tGrad.ref();
209 
210  const labelUList& owner = mesh.owner();
211  const labelUList& neighbour = mesh.neighbour();
212 
213  const areaVectorField& C = mesh.areaCentres();
214  const edgeVectorField& Cf = mesh.edgeCentres();
215 
216  // Create limiter
217  scalarField limiter(vvf.internalField().size(), 1.0);
218 
219  const scalar rk = (1.0/k_ - 1.0);
220 
221  forAll(owner, edgei)
222  {
223  const label own = owner[edgei];
224  const label nei = neighbour[edgei];
225 
226  // owner side
227  vector gradf((Cf[edgei] - C[own]) & g[own]);
228 
229  scalar vsfOwn = gradf & vvf[own];
230  scalar vsfNei = gradf & vvf[nei];
231 
232  scalar maxEdge = max(vsfOwn, vsfNei);
233  scalar minEdge = min(vsfOwn, vsfNei);
234  const scalar maxMinEdge = rk*(maxEdge - minEdge);
235  maxEdge += maxMinEdge;
236  minEdge -= maxMinEdge;
237 
238  limitEdge
239  (
240  limiter[own],
241  maxEdge - vsfOwn,
242  minEdge - vsfOwn,
243  magSqr(gradf)
244  );
245 
246 
247  // neighbour side
248  gradf = (Cf[edgei] - C[nei]) & g[nei];
249 
250  vsfOwn = gradf & vvf[own];
251  vsfNei = gradf & vvf[nei];
252 
253  maxEdge = max(vsfOwn, vsfNei);
254  minEdge = min(vsfOwn, vsfNei);
255 
256  limitEdge
257  (
258  limiter[nei],
259  maxEdge - vsfNei,
260  minEdge - vsfNei,
261  magSqr(gradf)
262  );
263  }
264 
265 
266  // Lambda expression to update limiter for boundary edges
267  auto updateLimiter = [&](const label patchi, const vectorField& fld) -> void
268  {
269  const labelUList& pOwner = mesh.boundary()[patchi].edgeFaces();
270  const vectorField& pCf = Cf.boundaryField()[patchi];
271 
272  forAll(pOwner, edgei)
273  {
274  const label own = pOwner[edgei];
275 
276  const vector gradf((pCf[edgei] - C[own]) & g[own]);
277 
278  const scalar vsfOwn = gradf & vvf[own];
279  const scalar vsfNei = gradf & fld[edgei];
280 
281  scalar maxEdge = max(vsfOwn, vsfNei);
282  scalar minEdge = min(vsfOwn, vsfNei);
283  const scalar maxMinEdge = rk*(maxEdge - minEdge);
284  maxEdge += maxMinEdge;
285  minEdge -= maxMinEdge;
286 
287  limitEdge
288  (
289  limiter[own],
290  maxEdge - vsfOwn,
291  minEdge - vsfOwn,
292  magSqr(gradf)
293  );
294  }
295  };
296 
297 
298  const areaVectorField::Boundary& bvf = vvf.boundaryField();
299  forAll(bvf, patchi)
300  {
301  const faPatchVectorField& psf = bvf[patchi];
302 
303  if (psf.coupled())
304  {
305  updateLimiter(patchi, psf.patchNeighbourField());
306  }
307  else if (psf.fixesValue())
308  {
309  updateLimiter(patchi, psf);
310  }
311  }
312 
313  if (fa::debug)
314  {
315  Info<< "gradient limiter for: " << vvf.name()
316  << " max = " << gMax(limiter)
317  << " min = " << gMin(limiter)
318  << " average: " << gAverage(limiter) << endl;
319  }
320 
321  g.primitiveFieldRef() *= limiter;
322  g.correctBoundaryConditions();
324 
325  return tGrad;
326 }
327 
328 
329 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
330 
331 } // End namespace fa
332 
333 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
334 
335 } // End namespace Foam
336 
337 // ************************************************************************* //
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:60
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:88
const uniformDimensionedVectorField & g
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
volScalarField & C
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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:108
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:76
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Definition: areaFieldsFwd.H:72