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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "PecletNo.H"
30 #include "turbulenceModel.H"
31 #include "surfaceInterpolate.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace functionObjects
39 {
40  defineTypeNameAndDebug(PecletNo, 0);
41  addToRunTimeSelectionTable(functionObject, PecletNo, dictionary);
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 Foam::tmp<Foam::surfaceScalarField> Foam::functionObjects::PecletNo::rhoScale
49 (
50  const surfaceScalarField& phi
51 ) const
52 {
53  if (phi.dimensions() == dimMass/dimTime)
54  {
55  return phi/fvc::interpolate(lookupObject<volScalarField>(rhoName_));
56  }
57 
58  return phi;
59 }
60 
61 
62 bool Foam::functionObjects::PecletNo::calc()
63 {
64  if (foundObject<surfaceScalarField>(fieldName_))
65  {
66  tmp<volScalarField> nuEff;
67  if (mesh_.foundObject<turbulenceModel>(turbulenceModel::propertiesName))
68  {
69  const turbulenceModel& model =
70  lookupObject<turbulenceModel>
71  (
73  );
74 
75  nuEff = model.nuEff();
76  }
77  else if (mesh_.foundObject<dictionary>("transportProperties"))
78  {
79  const dictionary& model =
80  mesh_.lookupObject<dictionary>("transportProperties");
81 
83  (
84  IOobject
85  (
86  "nuEff",
87  mesh_.time().timeName(),
88  mesh_,
91  ),
92  mesh_,
93  dimensionedScalar("nu", dimViscosity, model)
94  );
95  }
96  else
97  {
99  << "Unable to determine the viscosity"
100  << exit(FatalError);
101  }
102 
103 
104  const surfaceScalarField& phi =
105  mesh_.lookupObject<surfaceScalarField>(fieldName_);
106 
107  return store
108  (
109  resultName_,
110  mag(rhoScale(phi))
111  /(
112  mesh_.magSf()
113  *mesh_.surfaceInterpolation::deltaCoeffs()
114  *fvc::interpolate(nuEff)
115  )
116  );
117  }
118 
119  return false;
120 }
121 
122 
123 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
124 
126 (
127  const word& name,
128  const Time& runTime,
129  const dictionary& dict
130 )
131 :
132  fieldExpression(name, runTime, dict, "phi"),
133  rhoName_("rho")
134 {
135  setResultName("Pe", "phi");
136  read(dict);
137 }
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
143 {
144  rhoName_ = dict.getOrDefault<word>("rho", "rho");
145 
146  return true;
147 }
148 
149 
150 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
runTime
engineTime & runTime
Definition: createEngineTime.H:13
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:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
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
PecletNo.H
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::functionObjects::PecletNo::PecletNo
PecletNo(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: PecletNo.C:126
Foam::compressible::turbulenceModel
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Definition: turbulentFluidThermoModel.H:63
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dimViscosity
const dimensionSet dimViscosity
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:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::fvc::interpolate
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.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::functionObjects::PecletNo::read
virtual bool read(const dictionary &)
Read the PecletNo data.
Definition: PecletNo.C:142
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::functionObjects::fieldExpression
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
Definition: fieldExpression.H:120
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::IOobject::NO_READ
Definition: IOobject.H:188
turbulenceModel.H