Loading...
Searching...
No Matches
limitTurbulenceViscosity.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) 2023 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29#include "fvMesh.H"
30#include "transportModel.H"
31#include "fluidThermo.H"
32#include "turbulenceModel.H"
33#include "processorFvPatch.H"
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
40namespace fv
41{
44}
45}
46
47
48// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49
51{
53 writeHeaderValue(os, "c", c_);
54 writeCommented(os, "Time");
55 writeTabbed(os, "nLimitedCells_[count]");
56 writeTabbed(os, "nLimitedCells_[%]");
58 os << endl;
59
60 writtenHeader_ = true;
61}
62
63
65{
66 const auto* turbPtr =
68 if (turbPtr)
69 {
70 return turbPtr->nu();
71 }
72
73 const auto* thermoPtr =
74 mesh_.cfindObject<fluidThermo>(fluidThermo::dictName);
75 if (thermoPtr)
76 {
77 return thermoPtr->nu();
78 }
79
80 const auto* laminarPtr =
81 mesh_.cfindObject<transportModel>("transportProperties");
82 if (laminarPtr)
83 {
84 return laminarPtr->nu();
85 }
86
87 const auto* dictPtr = mesh_.cfindObject<dictionary>("transportProperties");
88 if (dictPtr)
89 {
91 (
92 "nu",
94 mesh_,
95 dimensionedScalar("nu", dimViscosity, *dictPtr)
96 );
97 }
98
100 << "No valid model for laminar viscosity"
101 << exit(FatalError);
103 return nullptr;
104}
105
106
107// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
108
110(
111 const word& name,
112 const word& modelType,
113 const dictionary& dict,
114 const fvMesh& mesh
115)
116:
117 fv::cellSetOption(name, modelType, dict, mesh),
118 writeFile(mesh, name, typeName, dict, false),
119 nutName_("nut"),
120 c_(1e5)
121{
122 if (isActive())
123 {
125 }
126}
127
128
129// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
130
132{
133 if (!(fv::cellSetOption::read(dict) && writeFile::read(dict)))
134 {
135 return false;
136 }
137
138 coeffs_.readIfPresent("nut", nutName_);
139 coeffs_.readIfPresent("c", c_);
140
141 fieldNames_.resize(1, nutName_);
142
144
145 if (canResetFile())
146 {
147 resetFile(typeName);
148 }
149
150 if (canWriteHeader())
151 {
152 writeFileHeader(file());
154
155
156 return true;
157}
158
159
161{
162 const auto& tnu = this->nu();
163 const auto& nu = tnu();
164
165 const label nTotCells(returnReduce(cells_.size(), sumOp<label>()));
166
167 label nAboveMax = 0;
168 for (const label celli : cells_)
169 {
170 const scalar nutLim = c_*nu[celli];
171
172 if (nut[celli] > nutLim)
173 {
174 nut[celli] = nutLim;
175 ++nAboveMax;
176 }
177 }
178
179 reduce(nAboveMax, sumOp<label>());
180
181 // Percent, max 2 decimal places
182 const auto percent = [](scalar num, label denom) -> scalar
183 {
184 return (denom ? 1e-2*round(1e4*num/denom) : 0);
185 };
186
187 const scalar nAboveMaxPercent = percent(nAboveMax, nTotCells);
188
189 Info<< type() << ' ' << name_ << " limited " << nAboveMax << " ("
190 << nAboveMaxPercent << "%) of cells" << endl;
191
192 if (canWriteToFile())
193 {
194 file()
195 << mesh_.time().timeOutputValue() << token::TAB
196 << nAboveMax << token::TAB
197 << nAboveMaxPercent
198 << endl;
199 }
200
201
202 // Handle boundaries in the case of 'all'
204 {
205 const volScalarField::Boundary& nubf = nu.boundaryField();
206 volScalarField::Boundary& nutbf = nut.boundaryFieldRef();
207
208 for (auto& nutp : nutbf)
209 {
210 // if (!nutp.fixesValue())
211 {
212 const scalarField& nup = nubf[nutp.patch().index()];
213
214 forAll(nutp, facei)
215 {
216 scalar nutLim = c_*nup[facei];
217 if (nutp[facei] > nutLim)
218 {
219 nutp[facei] = nutLim;
220 }
221 }
222 }
223 }
224 }
225
226
227 if (nAboveMax)
228 {
229 // We've changed internal values so give boundary conditions
230 // opportunity to correct
231
232 // Note: calling nut.correctBoundaryConditions() will re-evaluate wall
233 // functions and potentially (re-)introduce out-of-bounds values
234 //nut.correctBoundaryConditions();
235 if (UPstream::parRun())
236 {
237 nut.boundaryFieldRef().evaluateCoupled<processorFvPatch>();
238 }
239 }
240}
241
242
243// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
@ NO_REGISTER
Do not request registration (bool: false).
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
static const word dictName
The dictionary name ("thermophysicalProperties").
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Fundamental fluid thermodynamic properties.
Definition fluidThermo.H:52
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition writeFile.C:334
void writeHeaderValue(Ostream &os, const string &property, const Type &value) const
Write a (commented) header property and value pair.
virtual bool canWriteToFile() const
Flag to allow writing to the file.
Definition writeFile.C:292
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
Definition writeFile.C:200
bool writtenHeader_
Flag to identify whether the header has been written.
Definition writeFile.H:157
virtual OFstream & file()
Return access to the file (if only 1).
Definition writeFile.C:270
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition writeFile.C:318
virtual void resetFile(const word &name)
Reset internal file pointer to new file with new name.
Definition writeFile.C:165
virtual bool canWriteHeader() const
Flag to allow writing the header.
Definition writeFile.C:304
virtual bool canResetFile() const
Flag to allow resetting the file.
Definition writeFile.C:298
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Intermediate abstract class for handling cell-set options for the derived fvOptions.
labelList cells_
Set of cells to apply source to.
virtual bool read(const dictionary &dict)
Read source dictionary.
cellSetOption(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
bool useSubMesh() const noexcept
True if sub-selection should be used.
virtual bool isActive()
Is the source active?
Corrects turbulence viscosity field (e.g. nut) within a specified region by applying a maximum limit,...
limitTurbulenceViscosity(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
virtual void correct(volScalarField &nut)
Correct the energy field.
word nutName_
Name of turbulence viscosity field.
tmp< volScalarField > nu() const
Return laminar viscosity.
virtual bool read(const dictionary &dict)
Read dictionary.
void writeFileHeader(Ostream &os)
Write file header information.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition fvOption.H:124
const word & name() const noexcept
Return const access to the source name.
Definition fvOptionI.H:24
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition fvOption.H:157
dictionary coeffs_
Dictionary containing source coefficients.
Definition fvOption.H:152
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition fvOption.C:41
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition fvOptionI.H:30
const word name_
Source name.
Definition fvOption.H:132
A class for managing temporary objects.
Definition tmp.H:75
@ TAB
Tab [isspace].
Definition token.H:142
Base-class for all transport models used by the incompressible turbulence models.
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
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
scalar nut
OBJstream os(runTime.globalPath()/outputName)
Namespace for finite-volume.
Namespace for OpenFOAM.
const dimensionSet dimViscosity
GeometricField< scalar, fvPatchField, volMesh > volScalarField
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
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.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
volScalarField & nu
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299