Loading...
Searching...
No Matches
CourantNo.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) 2020-2025 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 "CourantNo.H"
30#include "surfaceFields.H"
31#include "fvcSurfaceIntegrate.H"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace functionObjects
39{
42}
43}
44
45
46// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47
49{
51 {
53
55 (
56 (0.5*mesh_.time().deltaT())
57 * fvc::surfaceSum(Foam::mag(phi))().internalField()
58 / mesh_.V()
59 );
60
61 if (tCo().dimensions() == dimDensity)
62 {
63 tCo.ref() /= obr_.lookupObject<volScalarField>(rhoName_);
64 }
65
66
68
69 if (!resultPtr)
70 {
71 resultPtr = new volScalarField
72 (
73 IOobject
74 (
77 mesh_.thisDb(),
81 ),
82 mesh_,
85 );
86 regIOobject::store(resultPtr);
87 }
88 auto& result = *resultPtr;
89
90 result.internalFieldRef() = tCo;
91 result.correctBoundaryConditions();
92
93 return true;
94 }
96 return false;
97}
98
99
100// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
101
103(
104 const word& name,
105 const Time& runTime,
106 const dictionary& dict
107)
108:
110 rhoName_("rho")
111{
112 setResultName("Co", "phi");
113 read(dict);
114}
115
116
117// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118
120{
122
123 rhoName_ = dict.getOrDefault<word>("rho", "rho");
124
125 return true;
126}
127
128
129// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
dimensionedScalar deltaT() const
Return time step.
Definition TimeStateI.H:61
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
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 Courant number field for time-variant simulations.
Definition CourantNo.H:132
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
Definition CourantNo.C:112
CourantNo(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
Definition CourantNo.C:96
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
bool foundObject(const word &name, const bool verbose=true) const
Return true if required objects are found.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
word resultName_
Name of result field.
word fieldName_
Name of field to process.
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 fvMesh & mesh_
Reference to the fvMesh.
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
ObjectType * getObjectPtr(const word &fieldName) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
const objectRegistry & obr_
Reference to the region objectRegistry.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition fvMesh.H:376
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
bool store()
Register object with its registry and transfer ownership to the registry.
A class for managing temporary objects.
Definition tmp.H:75
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
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
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.
const dimensionSet dimDensity
dictionary dict
Foam::surfaceFields.