Peclet.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 | Copyright (C) 2013-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 "Peclet.H"
27 #include "volFields.H"
28 #include "dictionary.H"
29 #include "surfaceFields.H"
32 #include "surfaceInterpolate.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(Peclet, 0);
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const word& name,
47  const objectRegistry& obr,
48  const dictionary& dict,
49  const bool loadFromFiles
50 )
51 :
52  name_(name),
53  obr_(obr),
54  active_(true),
55  phiName_("phi"),
56  resultName_(name),
57  log_(true)
58 {
59  // Check if the available mesh is an fvMesh, otherwise deactivate
60  if (!isA<fvMesh>(obr_))
61  {
62  active_ = false;
64  << "No fvMesh available, deactivating " << name_ << nl
65  << endl;
66  }
67 
68  read(dict);
69 
70  if (active_)
71  {
72  const fvMesh& mesh = refCast<const fvMesh>(obr_);
73 
74  surfaceScalarField* PecletPtr
75  (
77  (
78  IOobject
79  (
80  resultName_,
81  mesh.time().timeName(),
82  mesh,
83  IOobject::NO_READ,
84  IOobject::NO_WRITE
85  ),
86  mesh,
87  dimensionedScalar("0", dimless, 0.0)
88  )
89  );
90 
91  mesh.objectRegistry::store(PecletPtr);
92  }
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
97 
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
105 {
106  if (active_)
107  {
108  log_.readIfPresent("log", dict);
109  dict.readIfPresent("phiName", phiName_);
110  dict.readIfPresent("resultName", resultName_);
111  }
112 }
113 
114 
116 {
117  typedef compressible::turbulenceModel cmpTurbModel;
118  typedef incompressible::turbulenceModel icoTurbModel;
119 
120  if (active_)
121  {
122  const fvMesh& mesh = refCast<const fvMesh>(obr_);
123 
124  // Obtain nuEff of muEff. Assumes that a compressible flux is present
125  // when using a compressible turbulence model, and an incompressible
126  // flux when using an incompressible turbulence model
127  tmp<volScalarField> nuOrMuEff;
129  {
130  const cmpTurbModel& model =
131  mesh.lookupObject<cmpTurbModel>
132  (
134  );
135 
136  nuOrMuEff = model.muEff();
137  }
138  else if
139  (
141  )
142  {
143  const icoTurbModel& model =
144  mesh.lookupObject<icoTurbModel>
145  (
147  );
148 
149  nuOrMuEff = model.nuEff();
150  }
151  else if (mesh.foundObject<dictionary>("transportProperties"))
152  {
153  const dictionary& model =
154  mesh.lookupObject<dictionary>("transportProperties");
155 
156  nuOrMuEff =
158  (
159  new volScalarField
160  (
161  IOobject
162  (
163  "nuEff",
164  mesh.time().timeName(),
165  mesh,
168  ),
169  mesh,
170  dimensionedScalar(model.lookup("nu"))
171  )
172  );
173  }
174  else
175  {
177  << "Unable to determine the viscosity"
178  << exit(FatalError);
179  }
180 
181  // Note: dimensions of phi will change depending on whether this is
182  // applied to an incompressible or compressible case
183  const surfaceScalarField& phi =
185 
187  const_cast<surfaceScalarField&>
188  (
189  mesh.lookupObject<surfaceScalarField>(resultName_)
190  );
191 
192  Peclet =
193  mag(phi)
194  /(
195  mesh.magSf()
196  *mesh.surfaceInterpolation::deltaCoeffs()
197  *fvc::interpolate(nuOrMuEff)
198  );
199  }
200 }
201 
202 
204 {
205  // Do nothing
206 }
207 
208 
210 {
211  // Do nothing
212 }
213 
214 
216 {
217  if (active_)
218  {
219  const surfaceScalarField& Peclet =
220  obr_.lookupObject<surfaceScalarField>(resultName_);
221 
222  if (log_) Info
223  << type() << " " << name_ << " output:" << nl
224  << " writing field " << Peclet.name() << nl
225  << endl;
226 
227  Peclet.write();
228  }
229 }
230 
231 
232 // ************************************************************************* //
Foam::Peclet
This function object calculates and outputs the Peclet number as a surfaceScalarField.
Definition: Peclet.H:106
volFields.H
Foam::Peclet::execute
virtual void execute()
Execute, currently does nothing.
Definition: Peclet.C:115
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::Peclet::read
virtual void read(const dictionary &)
Read the Peclet data.
Definition: Peclet.C:104
Foam::Peclet::timeSet
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
Definition: Peclet.C:209
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dictionary::readIfPresent
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Definition: dictionaryTemplates.C:94
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::fvc::interpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Definition: surfaceInterpolate.C:110
turbulentTransportModel.H
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::Peclet::~Peclet
virtual ~Peclet()
Destructor.
Definition: Peclet.C:98
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::fvMesh::magSf
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Definition: fvMeshGeometry.C:358
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::Peclet::Peclet
Peclet(const Peclet &)
Disallow default bitwise copy construct.
Foam::Peclet::end
virtual void end()
Execute at the final time-loop, currently does nothing.
Definition: Peclet.C:203
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::Peclet::name
virtual const word & name() const
Return name of the set of Peclet.
Definition: Peclet.H:164
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
surfaceInterpolate.H
Surface Interpolation.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::Peclet::write
virtual void write()
Calculate the Peclet and write.
Definition: Peclet.C:215
Foam::objectRegistry::foundObject
bool foundObject(const word &name) const
Is the named Type found?
Definition: objectRegistryTemplates.C:142
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: ThermalDiffusivity.H:48
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
dictionary.H
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::IncompressibleTurbulenceModel
Templated abstract base class for single-phase incompressible turbulence models.
Definition: IncompressibleTurbulenceModel.H:52
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Peclet.H
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
turbulentFluidThermoModel.H