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-------------------------------------------------------------------------------
11License
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
36namespace Foam
37{
38namespace functionObjects
39{
42}
43}
44
45
46// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47
48Foam::tmp<Foam::surfaceScalarField> Foam::functionObjects::PecletNo::rhoScale
49(
50 const surfaceScalarField& phi
51) const
52{
54 {
55 return phi/fvc::interpolate(lookupObject<volScalarField>(rhoName_));
56 }
57
58 return phi;
59}
60
61
62bool 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_,
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:
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// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
surfaceScalarField & phi
const dimensionSet & dimensions() const
Return dimensions.
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base-class for Time/database function objects.
Computes the Peclet number as a surfaceScalarField.
Definition: PecletNo.H:148
virtual bool read(const dictionary &)
Read the PecletNo data.
Definition: PecletNo.C:142
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
void setResultName(const word &typeName, const word &defaultArg)
Set the name of result field.
A class for managing temporary objects.
Definition: tmp.H:65
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
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.
Namespace for OpenFOAM.
const dimensionSet dimViscosity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
dictionary dict