Loading...
Searching...
No Matches
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"
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(
51) const
52{
53 if (phi.dimensions() == dimMass/dimTime)
54 {
56 }
57
58 return phi;
59}
60
61
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 "nuEff",
86 mesh_,
88 );
89 }
90 else
91 {
93 << "Unable to determine the viscosity"
94 << exit(FatalError);
95 }
96
97
98 const surfaceScalarField& phi =
99 mesh_.lookupObject<surfaceScalarField>(fieldName_);
100
101 return store
102 (
103 resultName_,
104 mag(rhoScale(phi))
105 /(
106 mesh_.magSf()
107 *mesh_.surfaceInterpolation::deltaCoeffs()
108 *fvc::interpolate(nuEff)
109 )
110 );
111 }
113 return false;
114}
115
116
117// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
118
120(
121 const word& name,
122 const Time& runTime,
123 const dictionary& dict
124)
125:
127 rhoName_("rho")
128{
129 setResultName("Pe", "phi");
130 read(dict);
131}
132
133
134// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
135
137{
138 rhoName_ = dict.getOrDefault<word>("rho", "rho");
139
140 return true;
141}
142
143
144// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base-class for Time/database function objects.
Computes the Peclet number as a surfaceScalarField.
Definition PecletNo.H:132
PecletNo(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
Definition PecletNo.C:113
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
Definition PecletNo.C:129
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
fieldExpression(const word &name, const Time &runTime, const dictionary &dict, const word &fieldName=word::null, const word &resultName=word::null)
Construct from name, Time and dictionary.
virtual bool calc()=0
Calculate the components of the field and return true if successful.
void setResultName(const word &typeName, const word &defaultArg)
Set the name of result field.
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
A class for managing temporary objects.
Definition tmp.H:75
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
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
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
dictionary dict