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 |
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 
44 Foam::Peclet::Peclet
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  rhoName_("rho")
57 {
58  // Check if the available mesh is an fvMesh, otherwise deactivate
59  if (!isA<fvMesh>(obr_))
60  {
61  active_ = false;
62  WarningIn
63  (
64  "Peclet::Peclet"
65  "("
66  "const word&, "
67  "const objectRegistry&, "
68  "const dictionary&, "
69  "const bool"
70  ")"
71  ) << "No fvMesh available, deactivating " << name_ << nl
72  << endl;
73  }
74 
75  read(dict);
76 
77  if (active_)
78  {
79  const fvMesh& mesh = refCast<const fvMesh>(obr_);
80 
81  surfaceScalarField* PecletPtr
82  (
84  (
85  IOobject
86  (
87  type(),
88  mesh.time().timeName(),
89  mesh,
92  ),
93  mesh,
94  dimensionedScalar("0", dimless, 0.0)
95  )
96  );
97 
98  mesh.objectRegistry::store(PecletPtr);
99  }
100 }
101 
102 
103 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
104 
106 {}
107 
108 
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110 
112 {
113  if (active_)
114  {
115  phiName_ = dict.lookupOrDefault<word>("phiName", "phi");
116  rhoName_ = dict.lookupOrDefault<word>("rhoName", "rho");
117  }
118 }
119 
120 
122 {
123  typedef compressible::turbulenceModel cmpTurbModel;
124  typedef incompressible::turbulenceModel icoTurbModel;
125 
126  if (active_)
127  {
128  const fvMesh& mesh = refCast<const fvMesh>(obr_);
129 
130  tmp<volScalarField> nuEff;
131  if (mesh.foundObject<cmpTurbModel>(turbulenceModel::propertiesName))
132  {
133  const cmpTurbModel& model =
134  mesh.lookupObject<cmpTurbModel>
135  (
137  );
138 
139  const volScalarField& rho =
140  mesh.lookupObject<volScalarField>(rhoName_);
141 
142  nuEff = model.muEff()/rho;
143  }
144  else if
145  (
146  mesh.foundObject<icoTurbModel>(turbulenceModel::propertiesName)
147  )
148  {
149  const icoTurbModel& model =
150  mesh.lookupObject<icoTurbModel>
151  (
153  );
154 
155  nuEff = model.nuEff();
156  }
157  else if (mesh.foundObject<dictionary>("transportProperties"))
158  {
159  const dictionary& model =
160  mesh.lookupObject<dictionary>("transportProperties");
161 
162  nuEff =
164  (
165  new volScalarField
166  (
167  IOobject
168  (
169  "nuEff",
170  mesh.time().timeName(),
171  mesh,
174  ),
175  mesh,
176  dimensionedScalar(model.lookup("nu"))
177  )
178  );
179  }
180  else
181  {
182  FatalErrorIn("void Foam::Peclet::write()")
183  << "Unable to determine the viscosity"
184  << exit(FatalError);
185  }
186 
187  const surfaceScalarField& phi =
188  mesh.lookupObject<surfaceScalarField>(phiName_);
189 
191  const_cast<surfaceScalarField&>
192  (
194  );
195 
196  Peclet =
197  mag(phi)
198  /(
199  mesh.magSf()
200  *mesh.surfaceInterpolation::deltaCoeffs()
201  *fvc::interpolate(nuEff)
202  );
203  }
204 }
205 
206 
208 {
209  if (active_)
210  {
211  execute();
212  }
213 }
214 
216 {
217  // Do nothing
218 }
219 
220 
222 {
223  if (active_)
224  {
225  const surfaceScalarField& Peclet =
226  obr_.lookupObject<surfaceScalarField>(type());
227 
228  Info<< type() << " " << name_ << " output:" << nl
229  << " writing field " << Peclet.name() << nl
230  << endl;
231 
232  Peclet.write();
233  }
234 }
235 
236 
237 // ************************************************************************* //
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
This function object calculates and outputs the Peclet number as a surfaceScalarField.
Definition: Peclet.H:62
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::surfaceFields.
virtual void read(const dictionary &)
Read the Peclet data.
Definition: Peclet.C:111
dimensioned< scalar > mag(const dimensioned< Type > &)
bool foundObject(const word &name) const
Is the named Type found?
Templated abstract base class for single-phase incompressible turbulence models.
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
Definition: Peclet.C:215
A class for handling words, derived from string.
Definition: word.H:59
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
virtual void end()
Execute at the final time-loop, currently does nothing.
Definition: Peclet.C:207
virtual void write()
Calculate the Peclet and write.
Definition: Peclet.C:221
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Surface Interpolation.
virtual ~Peclet()
Destructor.
Definition: Peclet.C:105
messageStream Info
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
virtual bool write() const
Write using setting from DB.
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
#define WarningIn(functionName)
Report a warning using Foam::Warning.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
const word & name() const
Return name.
Definition: IOobject.H:260
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
virtual void execute()
Execute, currently does nothing.
Definition: Peclet.C:121
surfaceScalarField & phi
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
error FatalError
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Registry of regIOobjects.
static const word propertiesName
Default name of the turbulence properties dictionary.
bool read(const char *, int32_t &)
Definition: int32IO.C:87
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:118
defineTypeNameAndDebug(combustionModel, 0)