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