consumptionSpeed.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) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "consumptionSpeed.H"
29 
30 /* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */
31 
32 namespace Foam
33 {
35 }
36 
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
40 Foam::consumptionSpeed::consumptionSpeed
41 (
42  const dictionary& dict
43 )
44 : omega0_(dict.get<scalar>("omega0")),
45  eta_(dict.get<scalar>("eta")),
46  sigmaExt_(dict.get<scalar>("sigmaExt")),
47  omegaMin_(dict.get<scalar>("omegaMin"))
48 {}
49 
50 
51 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
52 
54 {}
55 
56 
57 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
58 
59 Foam::scalar Foam::consumptionSpeed::omega0Sigma
60 (
61  scalar sigma,
62  scalar a
63 ) const
64 {
65  scalar omega0 = 0.0;
66 
67  if (sigma < sigmaExt_)
68  {
69  omega0 = max
70  (
71  a*omega0_*(1.0 - exp(eta_*max(sigma, 0.0))),
72  omegaMin_
73  ) ;
74  }
75 
76  return omega0;
77 }
78 
79 
80 Foam::tmp<Foam::volScalarField> Foam::consumptionSpeed::omega0Sigma
81 (
82  const volScalarField& sigma
83 )
84 {
85  auto tomega0 = volScalarField::New
86  (
87  "omega0",
89  sigma.mesh(),
90  dimensionedScalar(dimensionSet(1, -2, -1, 0, 0, 0, 0), Zero)
91  );
92  auto& omega0 = tomega0.ref();
93 
94  volScalarField::Internal& iomega0 = omega0;
95 
96  forAll(iomega0, celli)
97  {
98  iomega0[celli] = omega0Sigma(sigma[celli], 1.0);
99  }
100 
101  volScalarField::Boundary& bomega0 = omega0.boundaryFieldRef();
102 
103  forAll(bomega0, patchi)
104  {
105  forAll(bomega0[patchi], facei)
106  {
107  bomega0[patchi][facei] =
108  omega0Sigma
109  (
110  sigma.boundaryField()[patchi][facei],
111  1.0
112  );
113  }
114  }
115 
116  return tomega0;
117 }
118 
119 
120 void Foam::consumptionSpeed::read(const dictionary& dict)
121 {
122  dict.readEntry("omega0", omega0_);
123  dict.readEntry("eta", eta_);
124  dict.readEntry("sigmaExt", sigmaExt_);
125  dict.readEntry("omegaMin", omegaMin_);
126 }
127 
128 
129 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
DimensionedField< scalar, volMesh > Internal
The internal field type from which this GeometricField is derived.
void read(const dictionary &dict)
Update properties.
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
dimensionedScalar exp(const dimensionedScalar &ds)
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
defineTypeNameAndDebug(combustionModel, 0)
virtual ~consumptionSpeed()
Destructor.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Do not request registration (bool: false)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127