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 -------------------------------------------------------------------------------
10 License
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"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace fv
41 {
42  defineTypeNameAndDebug(limitTurbulenceViscosity, 0);
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49 
51 {
52  writeHeaderValue(os, "Nut", nutName_);
53  writeHeaderValue(os, "c", c_);
54  writeCommented(os, "Time");
55  writeTabbed(os, "nLimitedCells_[count]");
56  writeTabbed(os, "nLimitedCells_[%]");
57 
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  {
90  return volScalarField::New
91  (
92  "nu",
94  mesh_,
95  dimensionedScalar("nu", dimViscosity, *dictPtr)
96  );
97  }
98 
100  << "No valid model for laminar viscosity"
101  << exit(FatalError);
102 
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 {
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());
153  }
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 // ************************************************************************* //
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
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
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
Tab [isspace].
Definition: token.H:129
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
const dimensionSet dimViscosity
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1061
virtual bool read(const dictionary &dict)
Read source dictionary.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
bool writtenHeader_
Flag to identify whether the header has been written.
Definition: writeFile.H:157
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
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
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:313
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
virtual bool isActive()
Is the source active?
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:48
void writeFileHeader(Ostream &os)
Write file header information.
limitTurbulenceViscosity(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
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())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
void writeHeaderValue(Ostream &os, const string &property, const Type &value) const
Write a (commented) header property and value pair.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
OBJstream os(runTime.globalPath()/outputName)
scalar c_
Limiting coefficient [].
tmp< volScalarField > nu() const
Return laminar viscosity.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
Definition: fluidThermo.C:96
messageStream Info
Information stream (stdout output on master, null elsewhere)
Corrects turbulence viscosity field (e.g. nut) within a specified region by applying a maximum limit...
Intermediate abstract class for handling cell-set options for the derived fvOptions.
word nutName_
Name of turbulence viscosity field.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual bool read(const dictionary &dict)
Read dictionary.
bool useSubMesh() const noexcept
True if sub-selection should be used.
volScalarField & nu
virtual void correct(volScalarField &nut)
Correct the energy field.
Do not request registration (bool: false)
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
scalar nut
Namespace for OpenFOAM.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:123
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:329