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 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace functionObjects
40 {
41  defineTypeNameAndDebug(CourantNo, 0);
42  addToRunTimeSelectionTable(functionObject, CourantNo, dictionary);
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
50 Foam::functionObjects::CourantNo::byRho
51 (
52  const tmp<volScalarField::Internal>& Co
53 ) const
54 {
55  if (Co().dimensions() == dimDensity)
56  {
57  return Co/obr_.lookupObject<volScalarField>(rhoName_);
58  }
59 
60  return Co;
61 }
62 
63 
64 bool Foam::functionObjects::CourantNo::calc()
65 {
66  if (foundObject<surfaceScalarField>(fieldName_))
67  {
68  const surfaceScalarField& phi =
69  lookupObject<surfaceScalarField>(fieldName_);
70 
71  tmp<volScalarField::Internal> Coi
72  (
73  byRho
74  (
75  (0.5*mesh_.time().deltaT())
76  *fvc::surfaceSum(mag(phi))()()
77  /mesh_.V()
78  )
79  );
80 
81  auto* resultPtr = getObjectPtr<volScalarField>(resultName_);
82 
83  if (!resultPtr)
84  {
85  resultPtr = new volScalarField
86  (
87  IOobject
88  (
89  resultName_,
90  mesh_.time().timeName(),
91  mesh_,
95  ),
96  mesh_,
99  );
100  mesh_.objectRegistry::store(resultPtr);
101  }
102  auto& Co = *resultPtr;
103 
104  Co.internalFieldRef() = Coi;
105  Co.correctBoundaryConditions();
106 
107  return true;
108  }
109 
110  return false;
111 }
112 
113 
114 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
115 
117 (
118  const word& name,
119  const Time& runTime,
120  const dictionary& dict
121 )
122 :
123  fieldExpression(name, runTime, dict, "phi"),
124  rhoName_("rho")
125 {
126  setResultName("Co", "phi");
127  read(dict);
128 }
129 
130 
131 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
132 
134 {
136 
137  rhoName_ = dict.getOrDefault<word>("rho", "rho");
138 
139  return true;
140 }
141 
142 
143 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
Foam::surfaceFields.
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
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...
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
engineTime & runTime
virtual bool read(const dictionary &dict)
Read the fieldExpression data.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
CourantNo(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: CourantNo.C:110
Macros for easy insertion into run-time selection tables.
virtual bool read(const dictionary &)
Read the CourantNo data.
Definition: CourantNo.C:126
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
void setResultName(const word &typeName, const word &defaultArg)
Set the name of result field.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
Internal & internalFieldRef(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
const dimensionSet dimDensity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Intermediate class for handling field expression function objects (e.g. blendingFactor etc...
const objectRegistry & obr_
Reference to the region objectRegistry.
Nothing to be read.
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Request registration (bool: true)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127