Loading...
Searching...
No Matches
DeltaOmegaTildeDelta.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
30#include "fvcCurl.H"
31#include "maxDeltaxyz.H"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace LESModels
39{
42}
43}
44
45
46// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47
48void Foam::LESModels::DeltaOmegaTildeDelta::calcDelta()
49{
51
52 const volVectorField& U0 = turbulenceModel_.U();
53 tmp<volVectorField> tvorticity = fvc::curl(U0);
54 const volVectorField& vorticity = tvorticity.cref();
55 const volVectorField nvecvort
56 (
57 vorticity
58 / (
59 max
60 (
61 mag(vorticity),
63 )
64 )
65 );
66 tvorticity.clear();
67
68 const cellList& cells = mesh.cells();
69 const vectorField& cellCentres = mesh.cellCentres();
70 const vectorField& faceCentres = mesh.faceCentres();
71
72 forAll(cells, celli)
73 {
74 const labelList& cFaces = cells[celli];
75 const point& cc = cellCentres[celli];
76 const vector& nv = nvecvort[celli];
77
78 scalar deltaMaxTmp = 0;
79
80 for (const label facei : cFaces)
81 {
82 const point& fc = faceCentres[facei];
83 const scalar tmp = 2.0*mag(nv ^ (fc - cc));
84
85 if (tmp > deltaMaxTmp)
86 {
87 deltaMaxTmp = tmp;
88 }
89 }
90
91 delta_[celli] = deltaCoeff_*deltaMaxTmp;
92 }
93
94 const label nD = mesh.nGeometricD();
95
96 if (nD == 2)
97 {
99 << "Case is 2D, LES is not strictly applicable" << nl
100 << endl;
101 }
102 else if (nD != 3)
103 {
105 << "Case must be either 2D or 3D" << exit(FatalError);
106 }
107
108 // Handle coupled boundaries
109 delta_.correctBoundaryConditions();
110}
111
112
113// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
114
115Foam::LESModels::DeltaOmegaTildeDelta::DeltaOmegaTildeDelta
116(
117 const word& name,
119 const dictionary& dict
120)
121:
123 hmaxPtr_(nullptr),
124 deltaCoeff_
125 (
126 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
127 (
128 "deltaCoeff",
129 1.035
130 )
131 ),
132 requireUpdate_
133 (
134 dict.optionalSubDict(type() + "Coeffs").getOrDefault<bool>
135 (
136 "requireUpdate",
137 true
138 )
139 )
140{
141 if (dict.optionalSubDict(type() + "Coeffs").found("hmax"))
142 {
143 // User-defined hmax
144 hmaxPtr_ =
146 (
147 IOobject::groupName("hmax", turbulence.U().group()),
149 dict.optionalSubDict("hmaxCoeffs"),
150 "hmax"
151 );
152 }
153 else
154 {
155 Info<< "Employing " << maxDeltaxyz::typeName << " for hmax" << endl;
156
157 hmaxPtr_.reset
158 (
159 new maxDeltaxyz
160 (
161 IOobject::groupName("hmax", turbulence.U().group()),
163 dict.optionalSubDict("hmaxCoeffs")
164 )
165 );
166 }
168 calcDelta();
169}
170
171
172// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
173
175{
176 const dictionary& coeffsDict = dict.optionalSubDict(type() + "Coeffs");
177
178 coeffsDict.readIfPresent<scalar>("deltaCoeff", deltaCoeff_);
179 coeffsDict.readIfPresent<bool>("requireUpdate", requireUpdate_);
180
181 calcDelta();
182}
183
184
186{
187 if (turbulenceModel_.mesh().changing() && requireUpdate_)
188 {
189 hmaxPtr_->correct();
190 }
191
192 calcDelta();
193}
194
195
196// ************************************************************************* //
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.
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. In "2D-regions" (i....
virtual void read(const dictionary &)
Read the LESdelta dictionary.
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
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...
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 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
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
const fvMesh & mesh() const
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Calculate the curl of the given volField by constructing the Hodge-dual of the symmetric part of the ...
const cellShapeList & cells
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for LES SGS models.
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
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
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
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299