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 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 makeFaGradScheme(faceLimitedGrad)
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace fa
44 {
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 template<>
50 (
51  scalar& limiter,
52  const scalar& maxDelta,
53  const scalar& minDelta,
54  const scalar& extrapolate
55 )
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 template<class Type>
70 (
71  Type& limiter,
72  const Type& maxDelta,
73  const Type& minDelta,
74  const Type& extrapolate
75 )
76 {
77  for (direction cmpt = 0; cmpt < Type::nComponents; ++cmpt)
78  {
80  (
81  limiter.component(cmpt),
82  maxDelta.component(cmpt),
83  minDelta.component(cmpt),
84  extrapolate.component(cmpt)
85  );
86  }
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
91 
92 template<>
93 tmp<areaVectorField> faceLimitedGrad<scalar>::calcGrad
94 (
95  const areaScalarField& vsf,
96  const word& name
97 ) const
98 {
99  const faMesh& mesh = vsf.mesh();
100 
101  tmp<areaVectorField> tGrad = basicGradScheme_().calcGrad(vsf, name);
102 
103  if (k_ < SMALL)
104  {
105  return tGrad;
106  }
107 
108  areaVectorField& g = tGrad.ref();
109 
110  const labelUList& owner = mesh.owner();
111  const labelUList& neighbour = mesh.neighbour();
112 
113  const areaVectorField& C = mesh.areaCentres();
114  const edgeVectorField& Cf = mesh.edgeCentres();
115 
116  scalarField maxVsf(vsf.internalField());
117  scalarField minVsf(vsf.internalField());
118 
119  forAll(owner, facei)
120  {
121  const label own = owner[facei];
122  const label nei = neighbour[facei];
123 
124  const scalar vsfOwn = vsf[own];
125  const scalar vsfNei = vsf[nei];
126 
127  maxVsf[own] = max(maxVsf[own], vsfNei);
128  minVsf[own] = min(minVsf[own], vsfNei);
129 
130  maxVsf[nei] = max(maxVsf[nei], vsfOwn);
131  minVsf[nei] = min(minVsf[nei], vsfOwn);
132  }
133 
134 
135  // Lambda expression to update limiter for boundary edges
136  auto updateLimiter = [&](const label patchi, const scalarField& fld) -> void
137  {
138  const labelUList& pOwner = mesh.boundary()[patchi].edgeFaces();
139 
140  forAll(pOwner, facei)
141  {
142  const label own = pOwner[facei];
143  const scalar vsf = fld[facei];
144 
145  maxVsf[own] = max(maxVsf[own], vsf);
146  minVsf[own] = min(minVsf[own], vsf);
147  }
148  };
149 
150  const areaScalarField::Boundary& bsf = vsf.boundaryField();
151  forAll(bsf, patchi)
152  {
153  const faPatchScalarField& psf = bsf[patchi];
154 
155  if (psf.coupled())
156  {
157  updateLimiter(patchi, psf.patchNeighbourField());
158  }
159  else
160  {
161  updateLimiter(patchi, psf);
162  }
163  }
164 
165  maxVsf -= vsf;
166  minVsf -= vsf;
167 
168  if (k_ < 1.0)
169  {
170  const scalarField maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
171  maxVsf += maxMinVsf;
172  minVsf -= maxMinVsf;
173  }
174 
175 
176  // Create limiter
177  scalarField limiter(vsf.internalField().size(), 1.0);
178 
179  forAll(owner, facei)
180  {
181  const label own = owner[facei];
182  const label nei = neighbour[facei];
183 
184  // owner side
185  limitEdge
186  (
187  limiter[own],
188  maxVsf[own],
189  minVsf[own],
190  (Cf[facei] - C[own]) & g[own]
191  );
192 
193  // neighbour side
194  limitEdge
195  (
196  limiter[nei],
197  maxVsf[nei],
198  minVsf[nei],
199  (Cf[facei] - C[nei]) & g[nei]
200  );
201  }
202 
203  forAll(bsf, patchi)
204  {
205  const labelUList& pOwner = mesh.boundary()[patchi].edgeFaces();
206  const vectorField& pCf = Cf.boundaryField()[patchi];
207 
208  forAll(pOwner, pFacei)
209  {
210  const label own = pOwner[pFacei];
211 
212  limitEdge
213  (
214  limiter[own],
215  maxVsf[own],
216  minVsf[own],
217  (pCf[pFacei] - C[own]) & g[own]
218  );
219  }
220  }
221 
222  if (fa::debug)
223  {
224  Info<< "gradient limiter for: " << vsf.name()
225  << " max = " << gMax(limiter)
226  << " min = " << gMin(limiter)
227  << " average: " << gAverage(limiter) << endl;
228  }
229 
230  g.primitiveFieldRef() *= limiter;
231  g.correctBoundaryConditions();
233 
234  return tGrad;
235 }
236 
237 
238 template<>
239 tmp<areaTensorField> faceLimitedGrad<vector>::calcGrad
240 (
241  const areaVectorField& vsf,
242  const word& name
243 ) const
244 {
245  const faMesh& mesh = vsf.mesh();
246 
247  tmp<areaTensorField> tGrad = basicGradScheme_().grad(vsf, name);
248 
249  if (k_ < SMALL)
250  {
251  return tGrad;
252  }
253 
254  areaTensorField& g = tGrad.ref();
255 
256  const labelUList& owner = mesh.owner();
257  const labelUList& neighbour = mesh.neighbour();
258 
259  const areaVectorField& C = mesh.areaCentres();
260  const edgeVectorField& Cf = mesh.edgeCentres();
261 
262  vectorField maxVsf(vsf.internalField());
263  vectorField minVsf(vsf.internalField());
264 
265  forAll(owner, facei)
266  {
267  const label own = owner[facei];
268  const label nei = neighbour[facei];
269 
270  const vector& vsfOwn = vsf[own];
271  const vector& vsfNei = vsf[nei];
272 
273  maxVsf[own] = max(maxVsf[own], vsfNei);
274  minVsf[own] = min(minVsf[own], vsfNei);
275 
276  maxVsf[nei] = max(maxVsf[nei], vsfOwn);
277  minVsf[nei] = min(minVsf[nei], vsfOwn);
278  }
279 
280 
281  // Lambda expression to update limiter for boundary edges
282  auto updateLimiter = [&](const label patchi, const vectorField& fld) -> void
283  {
284  const labelUList& pOwner = mesh.boundary()[patchi].edgeFaces();
285 
286  forAll(pOwner, facei)
287  {
288  const label own = pOwner[facei];
289  const vector& vsf = fld[facei];
290 
291  maxVsf[own] = max(maxVsf[own], vsf);
292  minVsf[own] = min(minVsf[own], vsf);
293  }
294  };
295 
296  const areaVectorField::Boundary& bsf = vsf.boundaryField();
297  forAll(bsf, patchi)
298  {
299  const faPatchVectorField& psf = bsf[patchi];
300 
301  if (psf.coupled())
302  {
303  updateLimiter(patchi, psf.patchNeighbourField());
304  }
305  else
306  {
307  updateLimiter(patchi, psf);
308  }
309  }
310 
311  maxVsf -= vsf;
312  minVsf -= vsf;
313 
314  if (k_ < 1.0)
315  {
316  const vectorField maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
317  maxVsf += maxMinVsf;
318  minVsf -= maxMinVsf;
319 
320  //maxVsf *= 1.0/k_;
321  //minVsf *= 1.0/k_;
322  }
323 
324 
325  // Create limiter
326  vectorField limiter(vsf.internalField().size(), vector::one);
327 
328  forAll(owner, facei)
329  {
330  const label own = owner[facei];
331  const label nei = neighbour[facei];
332 
333  // owner side
334  limitEdge
335  (
336  limiter[own],
337  maxVsf[own],
338  minVsf[own],
339  (Cf[facei] - C[own]) & g[own]
340  );
341 
342  // neighbour side
343  limitEdge
344  (
345  limiter[nei],
346  maxVsf[nei],
347  minVsf[nei],
348  (Cf[facei] - C[nei]) & g[nei]
349  );
350  }
351 
352  forAll(bsf, patchi)
353  {
354  const labelUList& pOwner = mesh.boundary()[patchi].edgeFaces();
355  const vectorField& pCf = Cf.boundaryField()[patchi];
356 
357  forAll(pOwner, pFacei)
358  {
359  const label own = pOwner[pFacei];
360 
361  limitEdge
362  (
363  limiter[own],
364  maxVsf[own],
365  minVsf[own],
366  ((pCf[pFacei] - C[own]) & g[own])
367  );
368  }
369  }
370 
371  if (fa::debug)
372  {
373  Info<< "gradient limiter for: " << vsf.name()
374  << " max = " << gMax(limiter)
375  << " min = " << gMin(limiter)
376  << " average: " << gAverage(limiter) << endl;
377  }
378 
379  tensorField& gIf = g.primitiveFieldRef();
380 
381  forAll(gIf, celli)
382  {
383  gIf[celli] = tensor
384  (
385  cmptMultiply(limiter[celli], gIf[celli].x()),
386  cmptMultiply(limiter[celli], gIf[celli].y()),
387  cmptMultiply(limiter[celli], gIf[celli].z())
388  );
389  }
390 
391  g.correctBoundaryConditions();
393 
394  return tGrad;
395 }
396 
397 
398 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
399 
400 } // End namespace fa
401 
402 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
403 
404 } // End namespace Foam
405 
406 // ************************************************************************* //
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: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
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
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))
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:108
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:76
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:72