73 const auto* thermoPtr =
77 return thermoPtr->nu();
80 const auto* laminarPtr =
81 mesh_.cfindObject<transportModel>(
"transportProperties");
84 return laminarPtr->nu();
87 const auto* dictPtr = mesh_.cfindObject<
dictionary>(
"transportProperties");
100 <<
"No valid model for laminar viscosity"
112 const word& modelType,
138 coeffs_.readIfPresent(
"nut", nutName_);
139 coeffs_.readIfPresent(
"c", c_);
141 fieldNames_.resize(1, nutName_);
150 if (canWriteHeader())
152 writeFileHeader(file());
162 const auto& tnu = this->
nu();
163 const auto&
nu = tnu();
168 for (
const label celli : cells_)
170 const scalar nutLim = c_*
nu[celli];
172 if (
nut[celli] > nutLim)
179 reduce(nAboveMax, sumOp<label>());
182 const auto percent = [](scalar num, label denom) -> scalar
184 return (denom ? 1
e-2*round(1e4*num/denom) : 0);
187 const scalar nAboveMaxPercent = percent(nAboveMax, nTotCells);
189 Info<<
type() <<
' ' << name_ <<
" limited " << nAboveMax <<
" ("
190 << nAboveMaxPercent <<
"%) of cells" <<
endl;
192 if (canWriteToFile())
195 << mesh_.time().timeOutputValue() <<
token::TAB
208 for (
auto& nutp : nutbf)
212 const scalarField& nup = nubf[nutp.patch().index()];
216 scalar nutLim = c_*nup[facei];
217 if (nutp[facei] > nutLim)
219 nutp[facei] = nutLim;
237 nut.boundaryFieldRef().evaluateCoupled<processorFvPatch>();
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,...
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
static bool & parRun() noexcept
Test if this a parallel run.
static const word dictName
The dictionary name ("thermophysicalProperties").
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Fundamental fluid thermodynamic properties.
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
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.
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
bool writtenHeader_
Flag to identify whether the header has been written.
virtual OFstream & file()
Return access to the file (if only 1).
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
virtual void resetFile(const word &name)
Reset internal file pointer to new file with new name.
virtual bool canWriteHeader() const
Flag to allow writing the header.
virtual bool canResetFile() const
Flag to allow resetting the file.
Mesh data needed to do the Finite Volume discretisation.
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.
scalar c_
Limiting coefficient [].
Base abstract class for handling finite volume options (i.e. fvOption).
const word & name() const noexcept
Return const access to the source name.
const fvMesh & mesh_
Reference to the mesh database.
wordList fieldNames_
Field names to apply source to - populated by derived models.
dictionary coeffs_
Dictionary containing source coefficients.
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
const word name_
Source name.
A class for managing temporary objects.
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.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
Namespace for finite-volume.
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.
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.
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.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define forAll(list, i)
Loop across all elements in list.