PecletNo.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-2016 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 "PecletNo.H"
27 #include "turbulenceModel.H"
28 #include "surfaceInterpolate.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace functionObjects
36 {
37  defineTypeNameAndDebug(PecletNo, 0);
38 
40  (
41  functionObject,
42  PecletNo,
43  dictionary
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 bool Foam::functionObjects::PecletNo::calc()
52 {
53  if (foundObject<surfaceScalarField>(fieldName_))
54  {
55  tmp<volScalarField> nuEff
56  (
57  mesh_.lookupObject<turbulenceModel>
58  (
60  ).nuEff()
61  );
62 
63  const surfaceScalarField& phi =
64  mesh_.lookupObject<surfaceScalarField>(fieldName_);
65 
66  return store
67  (
68  resultName_,
69  mag(phi)
70  /(
71  mesh_.magSf()
72  *mesh_.surfaceInterpolation::deltaCoeffs()
73  *fvc::interpolate(nuEff)
74  )
75  );
76  }
77  else
78  {
79  return false;
80  }
81 }
82 
83 
84 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
85 
87 (
88  const word& name,
89  const Time& runTime,
90  const dictionary& dict
91 )
92 :
93  fieldExpression(name, runTime, dict, "phi")
94 {
95  setResultName("Pe", "phi");
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
100 
102 {}
103 
104 
105 // ************************************************************************* //
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
PecletNo(const word &name, const Time &runTime, const dictionary &dict)
Construct for given objectRegistry and dictionary.
Definition: PecletNo.C:87
virtual ~PecletNo()
Destructor.
Definition: PecletNo.C:101
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
addToRunTimeSelectionTable(functionObject, add, dictionary)
Namespace for OpenFOAM.