PecletNo.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2022 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 "momentumTransportModel.H"
28 #include "surfaceInterpolate.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace functionObjects
36 {
38 
40  (
42  PecletNo,
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  (
58  );
59 
60  const surfaceScalarField& phi =
62 
63  return store
64  (
66  mag(phi)
67  /(
68  mesh_.magSf()
69  *mesh_.surfaceInterpolation::deltaCoeffs()
71  )
72  );
73  }
74  else
75  {
76  cannotFindObject<surfaceScalarField>(fieldName_);
77 
78  return false;
79  }
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84 
86 (
87  const word& name,
88  const Time& runTime,
89  const dictionary& dict
90 )
91 :
92  fieldExpression(name, runTime, dict, "Pe", "phi")
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
97 
99 {}
100 
101 
102 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Abstract base-class for Time/database functionObjects.
Calculates and outputs the Peclet number as a surfaceScalarField.
Definition: PecletNo.H:57
virtual ~PecletNo()
Destructor.
Definition: PecletNo.C:98
PecletNo(const word &name, const Time &runTime, const dictionary &dict)
Construct for given objectRegistry and dictionary.
Definition: PecletNo.C:86
word resultName_
Name of result field.
const word fieldName_
Name of field to process.
const fvMesh & mesh_
Reference to the fvMesh.
bool store(const tmp< ObjectType > &tfield)
Store the given field in the objectRegistry.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
const Type & lookupType(const word &group=word::null) const
Lookup and return the object of the given Type.
A class for handling words, derived from string.
Definition: word.H:62
const scalar nuEff
compressibleMomentumTransportModel momentumTransportModel
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Namespace for OpenFOAM.
SurfaceField< scalar > surfaceScalarField
dimensioned< scalar > mag(const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dictionary dict