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 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 basicThermo* thermoPtr =
59  getObjectPtr<basicThermo>(basicThermo::dictName);
60 
61  if (thermoPtr)
62  {
63  return fld*thermoPtr->rho();
64  }
65 
66  if (rhoInf_.value() < 0)
67  {
69  << type() << " " << name() << ": "
70  << "Incompressible calculation assumed, but no reference density "
71  << "set. Please set the entry 'rhoInf' to an appropriate value"
72  << nl
73  << exit(FatalError);
74  }
75 
76  return rhoInf_*fld;
77 }
78 
79 
81 Foam::functionObjects::proudmanAcousticPower::a() const
82 {
83  const basicThermo* thermoPtr =
84  getObjectPtr<basicThermo>(basicThermo::dictName);
85 
86  if (thermoPtr)
87  {
88  const basicThermo& thermo = *thermoPtr;
89  return sqrt(thermo.gamma()*thermo.p()/thermo.rho());
90  }
91 
92  return
94  (
95  IOobject
96  (
97  scopedName("a"),
98  mesh_.time().timeName(),
99  mesh_,
101  ),
102  mesh_,
103  aRef_
104  );
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
109 
111 (
112  const word& name,
113  const Time& runTime,
114  const dictionary& dict
115 )
116 :
118  alphaEps_(0.1),
119  rhoInf_("0", dimDensity, -1),
120  aRef_(dimVelocity, Zero)
121 {
122  read(dict);
123 
124  volScalarField* PAPtr
125  (
126  new volScalarField
127  (
128  IOobject
129  (
130  scopedName("P_A"),
131  mesh_.time().timeName(),
132  mesh_,
135  ),
136  mesh_,
138  )
139  );
140 
141  PAPtr->store();
142 
143  volScalarField* LPPtr
144  (
145  new volScalarField
146  (
147  IOobject
148  (
149  scopedName("L_P"),
150  mesh_.time().timeName(),
151  mesh_,
154  ),
155  mesh_,
157  )
158  );
159 
160  LPPtr->store();
161 }
162 
163 
164 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
165 
167 {
169  {
170  dict.readIfPresent("alphaEps", alphaEps_);
171  rhoInf_.readIfPresent("rhoInf", dict);
172  aRef_.readIfPresent("aRef", dict);
173 
174  return true;
175  }
176 
177  return false;
178 }
179 
180 
182 {
183  const turbulenceModel& turb =
184  lookupObject<turbulenceModel>(turbulenceModel::propertiesName);
185 
186  const volScalarField Mt(sqrt(2*turb.k())/a());
187 
188  volScalarField& P_A =
189  mesh_.lookupObjectRef<volScalarField>(scopedName("P_A"));
190 
191  P_A = rhoScale(alphaEps_*turb.epsilon()*pow5(Mt));
192 
193  volScalarField& L_P =
194  mesh_.lookupObjectRef<volScalarField>(scopedName("L_P"));
195 
196  L_P = 10*log10(P_A/dimensionedScalar("PRef", dimPower/dimVolume, 1e-12));
197 
198  return true;
199 }
200 
201 
203 {
204  Log << type() << " " << name() << " write:" << nl;
205 
206  const volScalarField& P_A =
207  mesh_.lookupObject<volScalarField>(scopedName("P_A"));
208 
209  Log << " writing field " << P_A.name() << nl;
210 
211  P_A.write();
212 
213  const volScalarField& L_P =
214  mesh_.lookupObject<volScalarField>(scopedName("L_P"));
215 
216  Log << " writing field " << L_P.name() << nl;
217 
218  L_P.write();
219 
220  Log << endl;
221 
222  return true;
223 }
224 
225 
226 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
volFields.H
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::functionObjects::proudmanAcousticPower::execute
virtual bool execute()
Calculate the Proudman acoustic power.
Definition: proudmanAcousticPower.C:181
basicThermo.H
Log
#define Log
Definition: PDRblock.C:35
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
Foam::functionObjects::proudmanAcousticPower::write
virtual bool write()
Write the Proudman acoustic power.
Definition: proudmanAcousticPower.C:202
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::log10
dimensionedScalar log10(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:263
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
fld
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))
Definition: gmvOutputLagrangian.H:23
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionary.H:458
Foam::dimPower
const dimensionSet dimPower
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::pow5
dimensionedScalar pow5(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:111
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::functionObjects::proudmanAcousticPower::proudmanAcousticPower
proudmanAcousticPower(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
Definition: proudmanAcousticPower.C:111
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::functionObject::name
const word & name() const
Return the name of this functionObject.
Definition: functionObject.C:131
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::functionObject::type
virtual const word & type() const =0
Runtime type information.
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
proudmanAcousticPower.H
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::GeometricField< scalar, fvPatchField, volMesh >
turb
compressible::turbulenceModel & turb
Definition: setRegionFluidFields.H:10
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::functionObjects::proudmanAcousticPower::read
virtual bool read(const dictionary &)
Read the Proudman acoustic power data.
Definition: proudmanAcousticPower.C:166
turbulenceModel.H