proudmanAcousticPower.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) 2019-2021 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 
28 #include "proudmanAcousticPower.H"
29 #include "volFields.H"
30 #include "basicThermo.H"
31 #include "turbulenceModel.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace functionObjects
39 {
40  defineTypeNameAndDebug(proudmanAcousticPower, 0);
42  (
43  functionObject,
44  proudmanAcousticPower,
45  dictionary
46  );
47 }
48 }
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
53 Foam::functionObjects::proudmanAcousticPower::rhoScale
54 (
55  const tmp<volScalarField>& fld
56 ) const
57 {
58  const auto* thermoPtr = getObjectPtr<basicThermo>(basicThermo::dictName);
59 
60  if (thermoPtr)
61  {
62  return fld*thermoPtr->rho();
63  }
64 
65  if (rhoInf_.value() < 0)
66  {
68  << type() << " " << name() << ": "
69  << "Incompressible calculation assumed, but no reference density "
70  << "set. Please set the entry 'rhoInf' to an appropriate value"
71  << nl
72  << exit(FatalError);
73  }
74 
75  return rhoInf_*fld;
76 }
77 
78 
80 Foam::functionObjects::proudmanAcousticPower::a() const
81 {
82  const auto* thermoPtr = getObjectPtr<basicThermo>(basicThermo::dictName);
83 
84  if (thermoPtr)
85  {
86  const basicThermo& thermo = *thermoPtr;
87  return sqrt(thermo.gamma()*thermo.p()/thermo.rho());
88  }
89 
90  return
92  (
93  IOobject
94  (
95  scopedName("a"),
96  mesh_.time().timeName(),
97  mesh_
98  ),
99  mesh_,
100  aRef_
101  );
102 }
103 
104 
106 Foam::functionObjects::proudmanAcousticPower::k() const
107 {
108  if (kName_ != "none")
109  {
110  return lookupObject<volScalarField>(kName_);
111  }
112 
113  const auto& turb =
114  lookupObject<turbulenceModel>(turbulenceModel::propertiesName);
115 
116  return turb.k();
117 }
118 
119 
121 Foam::functionObjects::proudmanAcousticPower::epsilon() const
122 {
123  if (epsilonName_ != "none")
124  {
125  return lookupObject<volScalarField>(epsilonName_);
126  }
127 
128  if (omegaName_ != "none")
129  {
130  // Construct epsilon on-the-fly
131  const auto& omega = lookupObject<volScalarField>(omegaName_);
132  const scalar betaStar = 0.09;
133  return betaStar*k()*omega;
134  }
135 
136  const auto& turb =
137  lookupObject<turbulenceModel>(turbulenceModel::propertiesName);
138 
139  return turb.epsilon();
140 }
141 
142 
143 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
144 
146 (
147  const word& name,
148  const Time& runTime,
149  const dictionary& dict
150 )
151 :
153  alphaEps_(0.1),
154  rhoInf_("0", dimDensity, -1),
155  aRef_(dimVelocity, Zero),
156  kName_("none"),
157  epsilonName_("none"),
158  omegaName_("none")
159 {
160  read(dict);
161 
162  volScalarField* PAPtr
163  (
164  new volScalarField
165  (
166  IOobject
167  (
168  scopedName("P_A"),
169  mesh_.time().timeName(),
170  mesh_,
174  ),
175  mesh_,
177  )
178  );
179 
180  PAPtr->store();
181 
182  volScalarField* LPPtr
183  (
184  new volScalarField
185  (
186  IOobject
187  (
188  scopedName("L_P"),
189  mesh_.time().timeName(),
190  mesh_,
194  ),
195  mesh_,
197  )
198  );
200  LPPtr->store();
201 }
202 
203 
204 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
205 
207 {
209  {
210  Info<< type() << " " << name() << nl;
211 
212  dict.readIfPresent("alphaEps", alphaEps_);
213  rhoInf_.readIfPresent("rhoInf", dict);
214  aRef_.readIfPresent("aRef", dict);
215 
216  if (dict.readIfPresent("k", kName_))
217  {
218  Info<< " k field: " << kName_ << endl;
219  }
220  else
221  {
222  Info<< " k field from turbulence model" << endl;
223  }
224 
225  if (dict.readIfPresent("epsilon", epsilonName_))
226  {
227  Info<< " epsilon field: " << epsilonName_ << endl;
228  }
229  else
230  {
231  Info<< " epsilon field from turbulence model (if needed)"
232  << endl;
233  }
234 
235  if (dict.readIfPresent("omega", omegaName_))
236  {
237  Info<< " omega field: " << omegaName_ << endl;
238  }
239  else
240  {
241  Info<< " omega field from turbulence model (if needed)" << endl;
242  }
243 
244  if (epsilonName_ != "none" && omegaName_ != "none")
245  {
247  << "either epsilon or omega field names can be set but not both"
248  << exit(FatalIOError);
249  }
250 
251  Info<< endl;
252 
253  return true;
254  }
255 
256  return false;
257 }
258 
259 
261 {
262  const volScalarField Mt(sqrt(2*k())/a());
263 
264  auto& P_A = mesh_.lookupObjectRef<volScalarField>(scopedName("P_A"));
265 
266  P_A = rhoScale(alphaEps_*epsilon()*pow5(Mt));
267 
268  auto& L_P = mesh_.lookupObjectRef<volScalarField>(scopedName("L_P"));
270  L_P = 10*log10(P_A/dimensionedScalar("PRef", dimPower/dimVolume, 1e-12));
271 
272  return true;
273 }
274 
275 
277 {
278  Log << type() << " " << name() << " write:" << nl;
279 
280  const auto& P_A = mesh_.lookupObject<volScalarField>(scopedName("P_A"));
281 
282  Log << " writing field " << P_A.name() << nl;
283 
284  P_A.write();
285 
286  const auto& L_P = mesh_.lookupObject<volScalarField>(scopedName("L_P"));
287 
288  Log << " writing field " << L_P.name() << nl;
289 
290  L_P.write();
291 
292  Log << endl;
293 
294  return true;
295 }
296 
297 
298 // ************************************************************************* //
virtual bool read(const dictionary &)
Read the Proudman acoustic power data.
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
const Type & value() const noexcept
Return const reference to value.
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
virtual bool execute()
Calculate the Proudman acoustic power.
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
compressible::turbulenceModel & turb
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual bool write()
Write the Proudman acoustic power.
bool store()
Register object with its registry and transfer ownership to the registry.
Definition: regIOobjectI.H:36
label k
Boltzmann constant.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
dimensionedScalar pow5(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
const word & name() const noexcept
Return the name of this functionObject.
psiReactionThermo & thermo
Definition: createFields.H:28
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
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual const word & type() const =0
Runtime type information.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
proudmanAcousticPower(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
const dimensionSet dimPower
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
const dimensionSet dimDensity
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
scalar epsilon
Nothing to be read.
#define Log
Definition: PDRblock.C:28
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual bool read(const dictionary &dict)
Read optional controls.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Request registration (bool: true)
word scopedName(const word &name) const
Return a scoped (prefixed) name.
const fvMesh & mesh_
Reference to the fvMesh.
dimensionedScalar log10(const dimensionedScalar &ds)
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity