CourantNo.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-2025 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 "CourantNo.H"
27 #include "surfaceFields.H"
28 #include "fvcSurfaceIntegrate.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace functionObjects
37 {
39 
41  (
43  CourantNo,
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
53 Foam::functionObjects::CourantNo::byRho
54 (
55  const tmp<volScalarField::Internal>& Co
56 ) const
57 {
58  if (Co().dimensions() == dimDensity)
59  {
60  return Co/obr_.lookupObject<volScalarField>(rhoName_);
61  }
62  else
63  {
64  return Co;
65  }
66 }
67 
68 
69 bool Foam::functionObjects::CourantNo::calc()
70 {
71  if (foundObject<surfaceScalarField>(fieldName_))
72  {
73  const surfaceScalarField& phi =
74  lookupObject<surfaceScalarField>(fieldName_);
75 
76  tmp<volScalarField> tCo
77  (
79  (
80  resultName_,
81  mesh_,
83  zeroGradientFvPatchScalarField::typeName
84  )
85  );
86 
87  tCo->internalFieldRef() =
88  byRho
89  (
90  (0.5*time_.deltaT())
91  *fvc::surfaceSum(mag(phi))/mesh_.V()
92  );
93 
94  tCo->correctBoundaryConditions();
95 
96  store(resultName_, tCo);
97 
98  return true;
99  }
100  else
101  {
102  cannotFindObject<surfaceScalarField>(fieldName_);
103 
104  return false;
105  }
106 }
107 
108 
109 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
110 
112 (
113  const word& name,
114  const Time& runTime,
115  const dictionary& dict
116 )
117 :
118  fieldExpression(name, runTime, dict, "Co", "phi")
119 {
120  read(dict);
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
125 
127 {}
128 
129 
130 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131 
133 {
135 
136  rhoName_ = dict.lookupOrDefault<word>("rho", "rho");
137 
138  return true;
139 }
140 
141 
142 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Calculates and outputs the maximum Courant Number.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base-class for Time/database functionObjects.
Calculates and outputs the Courant number as a volScalarField. The field is stored on the mesh databa...
Definition: CourantNo.H:60
CourantNo(const word &name, const Time &, const dictionary &)
Construct from Time and dictionary.
Definition: CourantNo.C:112
virtual ~CourantNo()
Destructor.
Definition: CourantNo.C:126
virtual bool read(const dictionary &)
Read the CourantNo data.
Definition: CourantNo.C:132
const objectRegistry & obr_
Reference to the objectRegistry.
virtual bool read(const dictionary &)
Read optional controls.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
tmp< VolInternalField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
Namespace for OpenFOAM.
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimDensity
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
Foam::surfaceFields.