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 -------------------------------------------------------------------------------
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 "DeltaOmegaTildeDelta.H"
30 #include "fvcCurl.H"
31 #include "maxDeltaxyz.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace LESModels
39 {
40  defineTypeNameAndDebug(DeltaOmegaTildeDelta, 0);
41  addToRunTimeSelectionTable(LESdelta, DeltaOmegaTildeDelta, dictionary);
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::LESModels::DeltaOmegaTildeDelta::calcDelta()
49 {
50  const fvMesh& mesh = turbulenceModel_.mesh();
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
110 }
111 
112 
113 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
114 
115 Foam::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()),
148  turbulence,
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()),
162  turbulence,
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 // ************************************************************************* //
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:300
const fvMesh & mesh() const
dictionary dict
List< cell > cellList
List of cell.
Definition: cellListFwd.H:39
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
virtual void read(const dictionary &)
Read the LESdelta dictionary.
const volVectorField & U() const
Access function to velocity field.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
Definition: LESdelta.H:146
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Abstract base class for LES deltas.
Definition: LESdelta.H:49
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const dimensionSet dimless
Dimensionless.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
tmp< GeometricField< Type, fvPatchField, volMesh > > curl(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcCurl.C:40
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
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
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
Calculate the curl of the given volField by constructing the Hodge-dual of the symmetric part of the ...
A class for handling words, derived from Foam::string.
Definition: word.H:63
Vector< scalar > vector
Definition: vector.H:57
volScalarField delta_
Definition: LESdelta.H:57
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...
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
defineTypeNameAndDebug(cubeRootVolDelta, 0)
vector point
Point is a vector.
Definition: point.H:37
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
Delta calculated by taking the maximum distance between the cell centre and any face centre...
const turbulenceModel & turbulenceModel_
Definition: LESdelta.H:55
Namespace for OpenFOAM.