DispersionRASModel.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) 2011-2015 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 "DispersionRASModel.H"
27 #include "demandDrivenData.H"
28 #include "turbulenceModel.H"
29 
30 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
31 
32 template<class CloudType>
35 {
36  const objectRegistry& obr = this->owner().mesh();
37  const word turbName =
38  IOobject::groupName
39  (
40  turbulenceModel::propertiesName,
41  this->owner().U().group()
42  );
43 
44  if (obr.foundObject<turbulenceModel>(turbName))
45  {
46  const turbulenceModel& model =
47  obr.lookupObject<turbulenceModel>(turbName);
48  return model.k();
49  }
50  else
51  {
53  (
54  "Foam::tmp<Foam::volScalarField>"
55  "Foam::DispersionRASModel<CloudType>::kModel() const"
56  )
57  << "Turbulence model not found in mesh database" << nl
58  << "Database objects include: " << obr.sortedToc()
59  << abort(FatalError);
60 
61  return tmp<volScalarField>(NULL);
62  }
63 }
64 
65 
66 template<class CloudType>
69 {
70  const objectRegistry& obr = this->owner().mesh();
71  const word turbName =
72  IOobject::groupName
73  (
74  turbulenceModel::propertiesName,
75  this->owner().U().group()
76  );
77 
78  if (obr.foundObject<turbulenceModel>(turbName))
79  {
80  const turbulenceModel& model =
81  obr.lookupObject<turbulenceModel>(turbName);
82  return model.epsilon();
83  }
84  else
85  {
87  (
88  "Foam::tmp<Foam::volScalarField>"
89  "Foam::DispersionRASModel<CloudType>::epsilonModel() const"
90  )
91  << "Turbulence model not found in mesh database" << nl
92  << "Database objects include: " << obr.sortedToc()
93  << abort(FatalError);
94 
95  return tmp<volScalarField>(NULL);
96  }
97 }
98 
99 
100 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
101 
102 template<class CloudType>
104 (
105  const dictionary&,
106  CloudType& owner
107 )
108 :
110  kPtr_(NULL),
111  ownK_(false),
112  epsilonPtr_(NULL),
113  ownEpsilon_(false)
114 {}
115 
116 
117 template<class CloudType>
119 (
121 )
122 :
124  kPtr_(dm.kPtr_),
125  ownK_(dm.ownK_),
126  epsilonPtr_(dm.epsilonPtr_),
127  ownEpsilon_(dm.ownEpsilon_)
128 {
129  dm.ownK_ = false;
130  dm.ownEpsilon_ = false;
131 }
132 
133 
134 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
135 
136 template<class CloudType>
138 {
139  cacheFields(false);
140 }
141 
142 
143 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
144 
145 template<class CloudType>
147 {
148  if (store)
149  {
150  tmp<volScalarField> tk = this->kModel();
151  if (tk.isTmp())
152  {
153  kPtr_ = tk.ptr();
154  ownK_ = true;
155  }
156  else
157  {
158  kPtr_ = tk.operator->();
159  ownK_ = false;
160  }
161 
162  tmp<volScalarField> tepsilon = this->epsilonModel();
163  if (tepsilon.isTmp())
164  {
165  epsilonPtr_ = tepsilon.ptr();
166  ownEpsilon_ = true;
167  }
168  else
169  {
170  epsilonPtr_ = tepsilon.operator->();
171  ownEpsilon_ = false;
172  }
173  }
174  else
175  {
176  if (ownK_ && kPtr_)
177  {
178  deleteDemandDrivenData(kPtr_);
179  ownK_ = false;
180  }
181  if (ownEpsilon_ && epsilonPtr_)
182  {
183  deleteDemandDrivenData(epsilonPtr_);
184  ownEpsilon_ = false;
185  }
186  }
187 }
188 
189 
190 template<class CloudType>
192 {
194 
195  os.writeKeyword("ownK") << ownK_ << token::END_STATEMENT << endl;
196  os.writeKeyword("ownEpsilon") << ownEpsilon_ << token::END_STATEMENT
197  << endl;
198 }
199 
200 
201 // ************************************************************************* //
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:148
virtual void cacheFields(const bool store)
Cache carrier fields.
const volScalarField * kPtr_
Turbulence k.
virtual ~DispersionRASModel()
Destructor.
virtual void write(Ostream &os) const
Write.
tmp< volScalarField > epsilonModel() const
Return the epsilon field from the turbulence model.
Base class for particle dispersion models based on RAS turbulence.
bool ownK_
Take ownership of the k field.
bool foundObject(const word &name) const
Is the named Type found?
const char *const group
Group name for atomic constants.
bool isTmp() const
Return true if this is really a temporary object.
Definition: tmpI.H:127
List< Key > sortedToc() const
Return the table of contents as a sorted list.
Definition: HashTable.C:216
void deleteDemandDrivenData(DataPtr &dataPtr)
A class for handling words, derived from string.
Definition: word.H:59
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
runTime write()
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< volScalarField > kModel() const
Return the k field from the turbulence model.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Template functions to aid in the implementation of demand driven data.
const volScalarField * epsilonPtr_
Turbulence epsilon.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
Abstract base class for turbulence models (RAS, LES and laminar).
error FatalError
Registry of regIOobjects.
bool ownEpsilon_
Take ownership of the epsilon field.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
const fvMesh & mesh() const
Return refernce to the mesh.
Definition: DSMCCloudI.H:41
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
DispersionRASModel(const dictionary &dict, CloudType &owner)
Construct from components.
A class for managing temporary objects.
Definition: PtrList.H:118
U
Definition: pEqn.H:82