turbulenceIntensity.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2018-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "turbulenceIntensity.H"
27 #include "momentumTransportModel.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace functionObjects
35 {
36  defineTypeNameAndDebug(turbulenceIntensity, 0);
37 
39  (
40  functionObject,
41  turbulenceIntensity,
42  dictionary
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void Foam::functionObjects::turbulenceIntensity::writeFileHeader(const label i)
51 {
52  writeHeader(file(), "I ()");
53  writeCommented(file(), "Time");
54  writeTabbed(file(), "min");
55  writeTabbed(file(), "max");
56  writeTabbed(file(), "average");
57  file() << endl;
58 }
59 
60 
61 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
62 
64 (
65  const word& name,
66  const Time& runTime,
67  const dictionary& dict
68 )
69 :
70  fvMeshFunctionObject(name, runTime, dict),
71  logFiles(obr_, name),
72  writeLocalObjects(obr_, log)
73 {
74  read(dict);
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
79 
81 {}
82 
83 
84 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
85 
87 {
90 
91  resetName(typeName);
92  resetLocalObjectName("I");
93 
94  return true;
95 }
96 
97 
99 {
100  if
101  (
102  mesh_.foundObject<momentumTransportModel>
103  (
104  momentumTransportModel::typeName
105  )
106  )
107  {
108  const momentumTransportModel& turbModel =
109  mesh_.lookupObject<momentumTransportModel>
110  (
111  momentumTransportModel::typeName
112  );
113 
114  volScalarField uPrime(sqrt((2.0/3.0)*turbModel.k()));
115 
116  word name("I");
117 
118  return
119  store
120  (
121  name,
122  uPrime
123  /max
124  (
125  max(uPrime, mag(turbModel.U())),
127  )
128  );
129  }
130  else
131  {
133  << "Unable to find turbulence model in the "
134  << "database" << exit(FatalError);
135 
136  return false;
137  }
138 }
139 
140 
142 {
143  Log << type() << " " << name() << " write:" << nl;
144 
146 
147  logFiles::write();
148 
150  mesh_.lookupObject<volScalarField>("I");
151 
152  const scalar minTurbulenceIntensity = min(turbulenceIntensity).value();
153  const scalar maxTurbulenceIntensity = max(turbulenceIntensity).value();
154  const scalar avgTurbulenceIntensity = turbulenceIntensity.average().value();
155 
156  if (Pstream::master())
157  {
158  Log << " I : min = " << minTurbulenceIntensity
159  << ", max = " << maxTurbulenceIntensity
160  << ", average = " << avgTurbulenceIntensity << nl;
161 
162  writeTime(file());
163  file()
164  << tab << minTurbulenceIntensity
165  << tab << maxTurbulenceIntensity
166  << tab << avgTurbulenceIntensity
167  << endl;
168  }
169 
170  Log << endl;
171 
172  return true;
173 }
174 
175 
176 // ************************************************************************* //
virtual bool write()
Write function.
Calculates the natural logarithm of the specified scalar field.
Definition: log.H:103
virtual bool write()
Write function.
Definition: logFiles.C:167
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
virtual bool write()
Write the turbulenceIntensity field.
static const char tab
Definition: Ostream.H:259
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
const ObjectType & lookupObject(const word &fieldName) const
Lookup object from the objectRegistry.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
void writeHeader(std::ostream &, const bool isBinary, const std::string &title)
Write header.
virtual bool execute()
Calculate the turbulenceIntensity field.
A class for handling words, derived from string.
Definition: word.H:59
turbulenceIntensity(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Evaluates and writes the turbulence intensity field &#39;I&#39;.
static const char nl
Definition: Ostream.H:260
const dimensionSet dimVelocity
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
defineTypeNameAndDebug(Qdot, 0)
dimensioned< Type > average() const
Calculate and return arithmetic average.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
virtual bool read(const dictionary &)
Read the turbulenceIntensity data.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
FunctionObject base class for managing a list of objects on behalf of the inheriting function object...
virtual bool read(const dictionary &)
Read the list of objects to be written.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
#define Log
Report write to Foam::Info if the local log switch is true.
dimensioned< scalar > mag(const dimensioned< Type > &)
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Namespace for OpenFOAM.
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:57