Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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 "edgeLimitedFaGrad.H"
30#include "gaussFaGrad.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34makeFaGradScheme(edgeLimitedGrad)
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38namespace Foam
39{
40
41// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42
43namespace fa
44{
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48template<class Type>
49inline 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
70template<>
71tmp<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 auto limits = gMinMax(limiter);
179 auto avg = gAverage(limiter);
180
181 Info<< "gradient limiter for: " << vsf.name()
182 << " min = " << limits.min()
183 << " max = " << limits.max()
184 << " average: " << avg << endl;
185 }
186
187 g.primitiveFieldRef() *= limiter;
188 g.correctBoundaryConditions();
190
191 return tGrad;
192}
193
194
195template<>
196tmp<areaTensorField> edgeLimitedGrad<vector>::calcGrad
197(
198 const areaVectorField& vvf,
199 const word& name
200) const
201{
202 const faMesh& mesh = vvf.mesh();
203
204 tmp<areaTensorField> tGrad = basicGradScheme_().grad(vvf, name);
205
206 if (k_ < SMALL)
207 {
208 return tGrad;
209 }
210
211 areaTensorField& g = tGrad.ref();
212
213 const labelUList& owner = mesh.owner();
214 const labelUList& neighbour = mesh.neighbour();
215
216 const areaVectorField& C = mesh.areaCentres();
217 const edgeVectorField& Cf = mesh.edgeCentres();
218
219 // Create limiter
220 scalarField limiter(vvf.internalField().size(), 1.0);
221
222 const scalar rk = (1.0/k_ - 1.0);
223
224 forAll(owner, edgei)
225 {
226 const label own = owner[edgei];
227 const label nei = neighbour[edgei];
228
229 // owner side
230 vector gradf((Cf[edgei] - C[own]) & g[own]);
231
232 scalar vsfOwn = gradf & vvf[own];
233 scalar vsfNei = gradf & vvf[nei];
234
235 scalar maxEdge = max(vsfOwn, vsfNei);
236 scalar minEdge = min(vsfOwn, vsfNei);
237 const scalar maxMinEdge = rk*(maxEdge - minEdge);
238 maxEdge += maxMinEdge;
239 minEdge -= maxMinEdge;
240
241 limitEdge
242 (
243 limiter[own],
244 maxEdge - vsfOwn,
245 minEdge - vsfOwn,
246 magSqr(gradf)
247 );
248
249
250 // neighbour side
251 gradf = (Cf[edgei] - C[nei]) & g[nei];
252
253 vsfOwn = gradf & vvf[own];
254 vsfNei = gradf & vvf[nei];
255
256 maxEdge = max(vsfOwn, vsfNei);
257 minEdge = min(vsfOwn, vsfNei);
258
259 limitEdge
260 (
261 limiter[nei],
262 maxEdge - vsfNei,
263 minEdge - vsfNei,
264 magSqr(gradf)
265 );
266 }
267
268
269 // Lambda expression to update limiter for boundary edges
270 auto updateLimiter = [&](const label patchi, const vectorField& fld) -> void
271 {
272 const labelUList& pOwner = mesh.boundary()[patchi].edgeFaces();
273 const vectorField& pCf = Cf.boundaryField()[patchi];
274
275 forAll(pOwner, edgei)
276 {
277 const label own = pOwner[edgei];
278
279 const vector gradf((pCf[edgei] - C[own]) & g[own]);
280
281 const scalar vsfOwn = gradf & vvf[own];
282 const scalar vsfNei = gradf & fld[edgei];
283
284 scalar maxEdge = max(vsfOwn, vsfNei);
285 scalar minEdge = min(vsfOwn, vsfNei);
286 const scalar maxMinEdge = rk*(maxEdge - minEdge);
287 maxEdge += maxMinEdge;
288 minEdge -= maxMinEdge;
289
290 limitEdge
291 (
292 limiter[own],
293 maxEdge - vsfOwn,
294 minEdge - vsfOwn,
295 magSqr(gradf)
296 );
297 }
298 };
299
300
301 const areaVectorField::Boundary& bvf = vvf.boundaryField();
302 forAll(bvf, patchi)
303 {
304 const faPatchVectorField& psf = bvf[patchi];
305
306 if (psf.coupled())
307 {
308 updateLimiter(patchi, psf.patchNeighbourField());
309 }
310 else if (psf.fixesValue())
311 {
312 updateLimiter(patchi, psf);
313 }
314 }
315
316 if (fa::debug)
317 {
318 auto limits = gMinMax(limiter);
319 auto avg = gAverage(limiter);
320
321 Info<< "gradient limiter for: " << vvf.name()
322 << " min = " << limits.min()
323 << " max = " << limits.max()
324 << " average: " << avg << endl;
325 }
326
327 g.primitiveFieldRef() *= limiter;
328 g.correctBoundaryConditions();
330
331 return tGrad;
332}
333
334
335// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
336
337} // End namespace fa
338
339// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
340
341} // End namespace Foam
342
343// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
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
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
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.
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
UList< label > labelUList
A UList of labels.
Definition UList.H:75
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
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