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 volScalarField::New
91  (
92  scopedName("a"),
94  mesh_,
95  aRef_
96  );
97 }
98 
99 
101 Foam::functionObjects::proudmanAcousticPower::k() const
102 {
103  if (kName_ != "none")
104  {
105  return lookupObject<volScalarField>(kName_);
106  }
107 
108  const auto& turb =
109  lookupObject<turbulenceModel>(turbulenceModel::propertiesName);
110 
111  return turb.k();
112 }
113 
114 
116 Foam::functionObjects::proudmanAcousticPower::epsilon() const
117 {
118  if (epsilonName_ != "none")
119  {
120  return lookupObject<volScalarField>(epsilonName_);
121  }
122 
123  if (omegaName_ != "none")
124  {
125  // Construct epsilon on-the-fly
126  const auto& omega = lookupObject<volScalarField>(omegaName_);
127  const scalar betaStar = 0.09;
128  return betaStar*k()*omega;
129  }
130 
131  const auto& turb =
132  lookupObject<turbulenceModel>(turbulenceModel::propertiesName);
133 
134  return turb.epsilon();
135 }
136 
137 
138 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
139 
141 (
142  const word& name,
143  const Time& runTime,
144  const dictionary& dict
145 )
146 :
148  alphaEps_(0.1),
149  rhoInf_("0", dimDensity, -1),
150  aRef_(dimVelocity, Zero),
151  kName_("none"),
152  epsilonName_("none"),
153  omegaName_("none")
154 {
155  read(dict);
156 
157  volScalarField* PAPtr
158  (
159  new volScalarField
160  (
161  IOobject
162  (
163  scopedName("P_A"),
164  mesh_.time().timeName(),
165  mesh_,
169  ),
170  mesh_,
172  )
173  );
174 
175  PAPtr->store();
176 
177  volScalarField* LPPtr
178  (
179  new volScalarField
180  (
181  IOobject
182  (
183  scopedName("L_P"),
184  mesh_.time().timeName(),
185  mesh_,
189  ),
190  mesh_,
192  )
193  );
195  LPPtr->store();
196 }
197 
198 
199 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
200 
202 {
204  {
205  Info<< type() << " " << name() << nl;
206 
207  dict.readIfPresent("alphaEps", alphaEps_);
208  rhoInf_.readIfPresent("rhoInf", dict);
209  aRef_.readIfPresent("aRef", dict);
210 
211  if (dict.readIfPresent("k", kName_))
212  {
213  Info<< " k field: " << kName_ << endl;
214  }
215  else
216  {
217  Info<< " k field from turbulence model" << endl;
218  }
219 
220  if (dict.readIfPresent("epsilon", epsilonName_))
221  {
222  Info<< " epsilon field: " << epsilonName_ << endl;
223  }
224  else
225  {
226  Info<< " epsilon field from turbulence model (if needed)"
227  << endl;
228  }
229 
230  if (dict.readIfPresent("omega", omegaName_))
231  {
232  Info<< " omega field: " << omegaName_ << endl;
233  }
234  else
235  {
236  Info<< " omega field from turbulence model (if needed)" << endl;
237  }
238 
239  if (epsilonName_ != "none" && omegaName_ != "none")
240  {
242  << "either epsilon or omega field names can be set but not both"
243  << exit(FatalIOError);
244  }
245 
246  Info<< endl;
247 
248  return true;
249  }
250 
251  return false;
252 }
253 
254 
256 {
257  const volScalarField Mt(sqrt(2*k())/a());
258 
259  auto& P_A = mesh_.lookupObjectRef<volScalarField>(scopedName("P_A"));
260 
261  P_A = rhoScale(alphaEps_*epsilon()*pow5(Mt));
262 
263  auto& L_P = mesh_.lookupObjectRef<volScalarField>(scopedName("L_P"));
265  L_P = 10*log10(P_A/dimensionedScalar("PRef", dimPower/dimVolume, 1e-12));
266 
267  return true;
268 }
269 
270 
272 {
273  Log << type() << " " << name() << " write:" << nl;
274 
275  const auto& P_A = mesh_.lookupObject<volScalarField>(scopedName("P_A"));
276 
277  Log << " writing field " << P_A.name() << nl;
278 
279  P_A.write();
280 
281  const auto& L_P = mesh_.lookupObject<volScalarField>(scopedName("L_P"));
282 
283  Log << " writing field " << L_P.name() << nl;
284 
285  L_P.write();
286 
287  Log << endl;
288 
289  return true;
290 }
291 
292 
293 // ************************************************************************* //
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:608
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.
proudmanAcousticPower(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
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...
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:637
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:180
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)
Do not request registration (bool: false)
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