Loading...
Searching...
No Matches
SLADelta.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) 2022 Upstream CFD GmbH
9 Copyright (C) 2022 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 "SLADelta.H"
30#include "wallDist.H"
31#include "maxDeltaxyz.H"
32#include "fvcGrad.H"
33#include "fvcCurl.H"
35
36#include "oneField.H"
38// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39
40namespace Foam
41{
42
43template<class Type, class WeightType = oneField>
45(
48 const WeightType& weights = oneField()
49)
50{
51 const fvMesh& mesh = field.mesh();
52
53 result == Zero;
54
55 auto tscaling = tmp<scalarField>::New(mesh.nCells(), Zero);
56 auto& scaling = tscaling.ref();
57
58 const auto& faceNeighbour = mesh.faceNeighbour();
59 const auto& faceOwner = mesh.faceOwner();
60
61 forAll(faceNeighbour, i)
62 {
63 const label own = faceOwner[i];
64 const label nbr = faceNeighbour[i];
65
66 result[own] += field[nbr];
67 result[nbr] += field[own];
68
69 scaling[own] = scaling[own] + weights[nbr];
70 scaling[nbr] = scaling[nbr] + weights[own];
71 }
72
73 const auto& pbm = mesh.boundaryMesh();
74
75 for (const polyPatch& pp : pbm)
76 {
77 const auto& pf = field.boundaryField()[pp.index()];
78 const labelUList& faceCells = pp.faceCells();
79
80 if (pf.coupled())
81 {
82 const scalarField nbrField(pf.patchNeighbourField());
83
84 forAll(faceCells, facei)
85 {
86 const label celli = faceCells[facei];
87 result[celli] += nbrField[facei];
88 scaling[celli] = scaling[celli] + weights[celli];
89 }
90 }
91 }
92
93 return tscaling;
94}
95
96} // End namespace Foam
97
99// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
100
101namespace Foam
102{
103namespace LESModels
104{
107}
108}
109
110
111// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
112
113void Foam::LESModels::SLADelta::calcDelta()
114{
115 const fvMesh& mesh = turbulenceModel_.mesh();
116
118 const volScalarField& nut = tnut.cref();
119
121 const volScalarField& nu = tnu.cref();
122
123 // Calculate vortex tilting measure, VTM
124 const volVectorField& U0 = turbulenceModel_.U();
125
126 tmp<volVectorField> tvorticity = fvc::curl(U0);
127 const volVectorField& vorticity = tvorticity.cref();
128
130 const volSymmTensorField& S = tS.cref();
131
132 const dimensionedScalar nuMin(nu.dimensions(), SMALL);
133 const dimensionedScalar eps(dimless/pow3(dimTime), SMALL);
134
136 (
137 max(scalar(1), 0.2*nu/max(nut, nuMin))
138 *sqrt(6.0)*mag((S & vorticity)^vorticity)
139 /max(magSqr(vorticity)*sqrt(3*magSqr(S) - sqr(tr(S))), eps)
140 );
141 vtm.correctBoundaryConditions();
142 tS.clear();
143
144 const dimensionedScalar vortMin(dimless/dimTime, SMALL);
145 const volVectorField nVecVort(vorticity/(max(mag(vorticity), vortMin)));
146 tvorticity.clear();
147
148
149 // Calculate averaged VTM
150 volScalarField vtmAve("vtmAve", vtm);
151 tmp<scalarField> tweights = sumNeighbours(vtm, vtmAve);
152
153 // Add cell centre values
154 vtmAve += vtm;
155
156 // Weights normalisation (add 1 for cell centres)
157 vtmAve.primitiveFieldRef() /= tweights + 1;
158
159
160 // Calculate DDES shielding function, fd
161 const volScalarField& y = wallDist::New(mesh).y();
162
163 const dimensionedScalar magGradUeps(dimless/dimTime, SMALL);
164 const dimensionedScalar nuEps(nu.dimensions(), ROOTSMALL);
165
166 tmp<volScalarField> tmagGradU = mag(fvc::grad(U0));
167
169 min
170 (
171 (nut + nu)/(max(tmagGradU, magGradUeps)*sqr(kappa_*y) + nuEps),
172 scalar(10)
173 );
174 tnut.clear();
175 tnu.clear();
176
177 const volScalarField fd(1.0 - tanh(pow(Cd1_*trd, Cd2_)));
178
179
180 // Assemble delta
181 const cellList& cells = mesh.cells();
182 const vectorField& cellCentres = mesh.cellCentres();
183 const vectorField& faceCentres = mesh.faceCentres();
184
185 forAll(cells, celli)
186 {
187 const labelList& cFaces = cells[celli];
188 const point& cc = cellCentres[celli];
189 const vector& nv = nVecVort[celli];
190
191 scalar deltaMaxTmp = 0;
192
193 for (const label facei : cFaces)
194 {
195 const point& fc = faceCentres[facei];
196 const scalar tmp = 2.0*mag(nv ^ (fc - cc));
197
198 if (tmp > deltaMaxTmp)
199 {
200 deltaMaxTmp = tmp;
201 }
202 }
203
204 scalar FKH =
205 max
206 (
207 FKHmin_,
208 min
209 (
210 FKHmax_,
211 FKHmin_
212 + (FKHmax_ - FKHmin_)*(vtmAve[celli] - a1_)/(a2_ - a1_)
213 )
214 );
215
216 if ((shielding_ == true) && (fd[celli] < (1.0 - epsilon_)))
217 {
218 FKH = 1;
219 }
220
221 delta_[celli] = deltaCoeff_*deltaMaxTmp*FKH;
222 }
223
224 const label nD = mesh.nGeometricD();
225
226 if (nD == 2)
227 {
229 << "Case is 2D, LES is not strictly applicable" << nl
230 << endl;
231 }
232 else if (nD != 3)
233 {
235 << "Case must be either 2D or 3D" << exit(FatalError);
236 }
237
238 // Handle coupled boundaries
239 delta_.correctBoundaryConditions();
240}
241
242
243// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
244
245Foam::LESModels::SLADelta::SLADelta
246(
247 const word& name,
249 const dictionary& dict
250)
251:
253 hmaxPtr_(nullptr),
254 deltaCoeff_
255 (
256 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
257 (
258 "deltaCoeff",
259 1.035
260 )
261 ),
262 requireUpdate_
263 (
264 dict.optionalSubDict(type() + "Coeffs").getOrDefault<bool>
265 (
266 "requireUpdate",
267 true
268 )
269 ),
270 shielding_
271 (
272 dict.optionalSubDict(type() + "Coeffs").getOrDefault<bool>
273 (
274 "shielding",
275 true
276 )
277 ),
278 FKHmin_
279 (
280 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
281 (
282 "FKHmin",
283 0.1
284 )
285 ),
286 FKHmax_
287 (
288 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
289 (
290 "FKHmax",
291 1.0
292 )
293 ),
294 a1_
295 (
296 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
297 (
298 "a1",
299 0.15
300 )
301 ),
302 a2_
303 (
304 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
305 (
306 "a2",
307 0.30
308 )
309 ),
310 epsilon_
311 (
312 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
313 (
314 "epsilon",
315 0.01
316 )
317 ),
318 kappa_
319 (
320 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
321 (
322 "kappa",
323 0.41
324 )
325 ),
326 Cd1_
327 (
328 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
329 (
330 "Cd1",
331 20
332 )
333 ),
334 Cd2_
335 (
336 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
337 (
338 "Cd2",
339 3
340 )
341 )
342{
343 if (dict.optionalSubDict(type() + "Coeffs").found("hmax"))
344 {
345 // User-defined hmax
346 hmaxPtr_ =
348 (
349 IOobject::groupName("hmax", turbulence.U().group()),
351 dict.optionalSubDict("hmaxCoeffs"),
352 "hmax"
353 );
354 }
355 else
356 {
357 Info<< "Employing " << maxDeltaxyz::typeName << " for hmax" << endl;
358
359 hmaxPtr_.reset
360 (
361 new maxDeltaxyz
362 (
363 IOobject::groupName("hmax", turbulence.U().group()),
365 dict.optionalSubDict("hmaxCoeffs")
366 )
367 );
368 }
369
370 if (mag(a2_ - a1_) < SMALL)
371 {
373 << "Model coefficients a1 = " << a1_
374 << ", and a2 = " << a2_ << " cannot be equal."
376 }
377}
378
379
380// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
381
383{
384 const dictionary& coeffsDict(dict.optionalSubDict(type() + "Coeffs"));
385
386 coeffsDict.readIfPresent<scalar>("deltaCoeff", deltaCoeff_);
387 coeffsDict.readIfPresent<bool>("requireUpdate", requireUpdate_);
388 coeffsDict.readIfPresent<scalar>("FKHmin", FKHmin_);
389 coeffsDict.readIfPresent<scalar>("FKHmax", FKHmax_);
390 coeffsDict.readIfPresent<scalar>("a1", a1_);
391 coeffsDict.readIfPresent<scalar>("a2", a2_);
392 coeffsDict.readIfPresent<scalar>("epsilon", epsilon_);
393 coeffsDict.readIfPresent<scalar>("kappa", kappa_);
394 coeffsDict.readIfPresent<scalar>("Cd1", Cd1_);
395 coeffsDict.readIfPresent<scalar>("Cd2", Cd2_);
396
397 calcDelta();
398}
399
400
402{
403 if (turbulenceModel_.mesh().changing() && requireUpdate_)
404 {
405 hmaxPtr_->correct();
406 }
407
408 calcDelta();
409}
410
411
412// ************************************************************************* //
scalar y
bool found
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
Generic GeometricField class.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Delta formulation that accounts for the orientation of the vorticity vector and a flow-sensitised fun...
Definition SLADelta.H:64
void read(const dictionary &)
Read the LESdelta dictionary.
Definition SLADelta.C:375
Abstract base class for LES deltas.
Definition LESdelta.H:50
LESdelta(const LESdelta &)=delete
No copy construct.
static autoPtr< LESdelta > New(const word &name, const turbulenceModel &turbulence, const dictionary &dict, const word &lookupName="delta")
Return a reference to the selected LES delta.
Definition LESdelta.C:61
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
Definition LESdelta.H:146
volScalarField delta_
Definition LESdelta.H:57
const turbulenceModel & turbulenceModel_
Definition LESdelta.H:55
static FOAM_NO_DANGLING_REFERENCE const wallDist & New(const fvMesh &mesh, Args &&... args)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Delta calculated by taking the maximum distance between the cell centre and any face centre....
A class representing the concept of a field of 1 used to avoid unnecessary manipulations for objects ...
Definition oneField.H:50
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
A class for managing temporary objects.
Definition tmp.H:75
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition tmpI.H:221
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
const fvMesh & mesh() const
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const volScalarField & y() const noexcept
Return reference to cached distance-to-wall field.
Definition wallDist.H:220
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
rDeltaTY field()
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
scalar nut
Calculate the curl of the given volField by constructing the Hodge-dual of the symmetric part of the ...
Calculate the gradient of the given field.
const cellShapeList & cells
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for LES SGS models.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > curl(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvcCurl.C:40
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
tmp< scalarField > sumNeighbours(const GeometricField< Type, fvPatchField, volMesh > &field, GeometricField< Type, fvPatchField, volMesh > &result, const WeightType &weights=oneField())
Definition SLADelta.C:38
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & nu
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299