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 | Copyright (C) 2013-2016 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 {
38  defineTypeNameAndDebug(CourantNo, 0);
39 
41  (
42  functionObject,
43  CourantNo,
44  dictionary
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::Internal> Coi
77  (
78  byRho
79  (
80  (0.5*mesh_.time().deltaT())
81  *fvc::surfaceSum(mag(phi))()()
82  /mesh_.V()
83  )
84  );
85 
86  if (foundObject<volScalarField>(resultName_))
87  {
88  volScalarField& Co
89  (
90  const_cast<volScalarField&>
91  (
92  lookupObject<volScalarField>(resultName_)
93  )
94  );
95 
96  Co.ref() = Coi();
97  Co.correctBoundaryConditions();
98  }
99  else
100  {
101  tmp<volScalarField> tCo
102  (
103  new volScalarField
104  (
105  IOobject
106  (
107  resultName_,
108  mesh_.time().timeName(),
109  mesh_,
112  ),
113  mesh_,
114  dimensionedScalar("0", dimless, 0.0),
115  zeroGradientFvPatchScalarField::typeName
116  )
117  );
118  tCo.ref().ref() = Coi();
119  tCo.ref().correctBoundaryConditions();
120  mesh_.objectRegistry::store(tCo.ptr());
121  }
122 
123  return true;
124  }
125  else
126  {
127  return false;
128  }
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
133 
135 (
136  const word& name,
137  const Time& runTime,
138  const dictionary& dict
139 )
140 :
141  fieldExpression(name, runTime, dict, "phi")
142 {
143  setResultName("Co", "phi");
144  read(dict);
145 }
146 
147 
148 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
149 
151 {}
152 
153 
154 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
155 
157 {
158  fieldExpression::read(dict);
159 
160  rhoName_ = dict.lookupOrDefault<word>("rho", "rho");
161 
162  return true;
163 }
164 
165 
166 // ************************************************************************* //
Foam::surfaceFields.
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
virtual bool read(const dictionary &)
Read the CourantNo data.
Definition: CourantNo.C:156
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
bool read(const char *, int32_t &)
Definition: int32IO.C:85
A class for handling words, derived from string.
Definition: word.H:59
virtual bool read(const dictionary &)
Read the fieldExpression data.
const dimensionSet dimDensity
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
addToRunTimeSelectionTable(functionObject, blendingFactor, dictionary)
virtual ~CourantNo()
Destructor.
Definition: CourantNo.C:150
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
A class for managing temporary objects.
Definition: PtrList.H:54
CourantNo(const word &name, const Time &, const dictionary &)
Construct from Time and dictionary.
Definition: CourantNo.C:135
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Namespace for OpenFOAM.