All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
diameterModel.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) 2011-2020 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 "diameterModel.H"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32  defineTypeNameAndDebug(diameterModel, 0);
33  defineRunTimeSelectionTable(diameterModel, dictionary);
34 }
35 
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
40 {
41  if (!dPtr_.valid())
42  {
43  dPtr_.reset
44  (
45  new volScalarField
46  (
47  IOobject
48  (
49  IOobject::groupName("d", phase_.name()),
50  phase_.time().timeName(),
51  phase_.mesh()
52  ),
53  phase_.mesh(),
55  )
56  );
57  }
58 
59  return dPtr_();
60 }
61 
63 {
64  if (!aPtr_.valid())
65  {
66  aPtr_.reset
67  (
68  new volScalarField
69  (
70  IOobject
71  (
72  IOobject::groupName("a", phase_.name()),
73  phase_.time().timeName(),
74  phase_.mesh()
75  ),
76  phase_.mesh(),
78  )
79  );
80  }
81 
82  return aPtr_();
83 }
84 
85 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
86 
88 (
89  const dictionary& diameterProperties,
90  const phaseModel& phase
91 )
92 :
93  diameterProperties_(diameterProperties),
94  phase_(phase),
95  dPtr_(nullptr),
96  aPtr_(nullptr)
97 {
98  if (diameterProperties.lookupOrDefault("storeD", false))
99  {
100  dRef();
101  }
102  if (diameterProperties.lookupOrDefault("storeA", false))
103  {
104  aRef();
105  }
106 }
107 
108 
109 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
110 
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
118 {
119  if (dPtr_.valid())
120  {
121  return dPtr_();
122  }
123  else
124  {
125  return calcD();
126  }
127 }
128 
129 
131 {
132  if (aPtr_.valid())
133  {
134  return aPtr_();
135  }
136  else
137  {
138  return calcA();
139  }
140 }
141 
142 
144 {
145  if (dPtr_.valid())
146  {
147  tmp<volScalarField> td = calcD();
148  if (td.isTmp())
149  {
150  dPtr_() = td;
151  }
152  }
153  if (aPtr_.valid())
154  {
155  tmp<volScalarField> tA = calcA();
156  if (tA.isTmp())
157  {
158  aPtr_() = tA;
159  }
160  }
161 }
162 
163 
164 bool Foam::diameterModel::read(const dictionary& phaseProperties)
165 {
166  diameterProperties_ = phaseProperties.optionalSubDict(type() + "Coeffs");
167 
168  return true;
169 }
170 
171 
172 // ************************************************************************* //
diameterModel(const dictionary &diameterProperties, const phaseModel &phase)
virtual ~diameterModel()
Destructor.
tmp< volScalarField > a() const
Return the surface area per unit volume.
volScalarField & aRef()
Get a reference to the stored surface area per unit volume field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
static word groupName(Name name, const word &group)
virtual tmp< volScalarField > calcA() const
Get the diameter field.
volScalarField & dRef()
Get a reference to the stored diameter field.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
tmp< volScalarField > d() const
Return the diameter.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual tmp< volScalarField > calcD() const
Get the diameter field.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void correct()
Correct the diameter field.
Namespace for OpenFOAM.