Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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 "faceLimitedFaGrad.H"
30#include "gaussFaGrad.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34makeFaGradScheme(faceLimitedGrad)
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38namespace Foam
39{
40
41// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42
43namespace fa
44{
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48template<>
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
68template<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
92template<>
93tmp<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 auto limits = gMinMax(limiter);
225 auto avg = gAverage(limiter);
226
227 Info<< "gradient limiter for: " << vsf.name()
228 << " min = " << limits.min()
229 << " max = " << limits.max()
230 << " average: " << avg << endl;
231 }
232
233 g.primitiveFieldRef() *= limiter;
234 g.correctBoundaryConditions();
236
237 return tGrad;
238}
239
240
241template<>
242tmp<areaTensorField> faceLimitedGrad<vector>::calcGrad
243(
244 const areaVectorField& vsf,
245 const word& name
246) const
247{
248 const faMesh& mesh = vsf.mesh();
249
250 tmp<areaTensorField> tGrad = basicGradScheme_().grad(vsf, name);
251
252 if (k_ < SMALL)
253 {
254 return tGrad;
255 }
256
257 areaTensorField& g = tGrad.ref();
258
259 const labelUList& owner = mesh.owner();
260 const labelUList& neighbour = mesh.neighbour();
261
262 const areaVectorField& C = mesh.areaCentres();
263 const edgeVectorField& Cf = mesh.edgeCentres();
264
265 vectorField maxVsf(vsf.internalField());
266 vectorField minVsf(vsf.internalField());
267
268 forAll(owner, facei)
269 {
270 const label own = owner[facei];
271 const label nei = neighbour[facei];
272
273 const vector& vsfOwn = vsf[own];
274 const vector& vsfNei = vsf[nei];
275
276 maxVsf[own] = max(maxVsf[own], vsfNei);
277 minVsf[own] = min(minVsf[own], vsfNei);
278
279 maxVsf[nei] = max(maxVsf[nei], vsfOwn);
280 minVsf[nei] = min(minVsf[nei], vsfOwn);
281 }
282
283
284 // Lambda expression to update limiter for boundary edges
285 auto updateLimiter = [&](const label patchi, const vectorField& fld) -> void
286 {
287 const labelUList& pOwner = mesh.boundary()[patchi].edgeFaces();
288
289 forAll(pOwner, facei)
290 {
291 const label own = pOwner[facei];
292 const vector& vsf = fld[facei];
293
294 maxVsf[own] = max(maxVsf[own], vsf);
295 minVsf[own] = min(minVsf[own], vsf);
296 }
297 };
298
299 const areaVectorField::Boundary& bsf = vsf.boundaryField();
300 forAll(bsf, patchi)
301 {
302 const faPatchVectorField& psf = bsf[patchi];
303
304 if (psf.coupled())
305 {
306 updateLimiter(patchi, psf.patchNeighbourField());
307 }
308 else
309 {
310 updateLimiter(patchi, psf);
311 }
312 }
313
314 maxVsf -= vsf;
315 minVsf -= vsf;
316
317 if (k_ < 1.0)
318 {
319 const vectorField maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
320 maxVsf += maxMinVsf;
321 minVsf -= maxMinVsf;
322
323 //maxVsf *= 1.0/k_;
324 //minVsf *= 1.0/k_;
325 }
326
327
328 // Create limiter
329 vectorField limiter(vsf.internalField().size(), vector::one);
330
331 forAll(owner, facei)
332 {
333 const label own = owner[facei];
334 const label nei = neighbour[facei];
335
336 // owner side
337 limitEdge
338 (
339 limiter[own],
340 maxVsf[own],
341 minVsf[own],
342 (Cf[facei] - C[own]) & g[own]
343 );
344
345 // neighbour side
346 limitEdge
347 (
348 limiter[nei],
349 maxVsf[nei],
350 minVsf[nei],
351 (Cf[facei] - C[nei]) & g[nei]
352 );
353 }
354
355 forAll(bsf, patchi)
356 {
357 const labelUList& pOwner = mesh.boundary()[patchi].edgeFaces();
358 const vectorField& pCf = Cf.boundaryField()[patchi];
359
360 forAll(pOwner, pFacei)
361 {
362 const label own = pOwner[pFacei];
363
364 limitEdge
365 (
366 limiter[own],
367 maxVsf[own],
368 minVsf[own],
369 ((pCf[pFacei] - C[own]) & g[own])
370 );
371 }
372 }
373
374 if (fa::debug)
375 {
376 auto limits = gMinMax(limiter);
377 auto avg = gAverage(limiter);
378
379 Info<< "gradient limiter for: " << vsf.name()
380 << " min = " << limits.min()
381 << " max = " << limits.max()
382 << " average: " << avg << endl;
383 }
384
385 tensorField& gIf = g.primitiveFieldRef();
386
387 forAll(gIf, celli)
388 {
389 gIf[celli] = tensor
390 (
391 cmptMultiply(limiter[celli], gIf[celli].x()),
392 cmptMultiply(limiter[celli], gIf[celli].y()),
393 cmptMultiply(limiter[celli], gIf[celli].z())
394 );
395 }
396
397 g.correctBoundaryConditions();
399
400 return tGrad;
401}
402
403
404// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
405
406} // End namespace fa
407
408// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
409
410} // End namespace Foam
411
412// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
scalar y
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
const uniformDimensionedVectorField & g
GeometricBoundaryField< scalar, faPatchField, areaMesh > Boundary
static const Form one
static void limitEdge(Type &limiter, const Type &maxDelta, const Type &minDelta, const Type &extrapolate)
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.
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.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
auto limits
Definition setRDeltaT.H:186
dynamicFvMesh & mesh
#define makeFaGradScheme(SS)
auto & name
Namespace for OpenFOAM.
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, faePatchField, edgeMesh > edgeVectorField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< tensor, faPatchField, areaMesh > areaTensorField
Tensor< scalar > tensor
Definition symmTensor.H:57
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
faPatchField< vector > faPatchVectorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< vector, faPatchField, areaMesh > areaVectorField
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
faPatchField< scalar > faPatchScalarField
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition faNVDscheme.C:31
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299