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 | 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 "DispersionRASModel.H"
30 #include "turbulenceModel.H"
31 
32 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
33 
34 template<class CloudType>
37 {
38  const objectRegistry& obr = this->owner().mesh();
39  const word turbName =
40  IOobject::groupName
41  (
42  turbulenceModel::propertiesName,
43  this->owner().U().group()
44  );
45 
46  const turbulenceModel* turb = obr.findObject<turbulenceModel>(turbName);
47 
48  if (turb)
49  {
50  return turb->k();
51  }
52 
54  << "Turbulence model not found in mesh database" << nl
55  << "Database objects include: " << obr.sortedToc()
56  << abort(FatalError);
57 
58  return nullptr;
59 }
60 
61 
62 template<class CloudType>
65 {
66  const objectRegistry& obr = this->owner().mesh();
67  const word turbName =
68  IOobject::groupName
69  (
70  turbulenceModel::propertiesName,
71  this->owner().U().group()
72  );
73 
74  const turbulenceModel* turb = obr.findObject<turbulenceModel>(turbName);
75 
76  if (turb)
77  {
78  return turb->epsilon();
79  }
80 
82  << "Turbulence model not found in mesh database" << nl
83  << "Database objects include: " << obr.sortedToc()
84  << abort(FatalError);
85 
86  return nullptr;
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
91 
92 template<class CloudType>
94 (
95  const dictionary&,
96  CloudType& owner
97 )
98 :
99  DispersionModel<CloudType>(owner),
100  kPtr_(nullptr),
101  epsilonPtr_(nullptr),
102  ownK_(false),
103  ownEpsilon_(false)
104 {}
105 
106 
107 template<class CloudType>
109 (
111 )
112 :
114  kPtr_(dm.kPtr_),
115  epsilonPtr_(dm.epsilonPtr_),
116  ownK_(dm.ownK_),
117  ownEpsilon_(dm.ownEpsilon_)
118 {
119  dm.ownK_ = false;
120  dm.ownEpsilon_ = false;
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
125 
126 template<class CloudType>
128 {
129  cacheFields(false);
130 }
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134 
135 template<class CloudType>
137 {
138  if (store)
139  {
140  tmp<volScalarField> tk = this->kModel();
141  if (tk.movable())
142  {
143  // Take ownership
144  kPtr_ = tk.ptr();
145  ownK_ = true;
146  }
147  else
148  {
149  kPtr_ = &tk.cref();
150  ownK_ = false;
151  }
152 
153  tmp<volScalarField> tepsilon = this->epsilonModel();
154  if (tepsilon.movable())
155  {
156  // Take ownership
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_)
169  {
170  deleteDemandDrivenData(kPtr_);
171  ownK_ = false;
172  }
173  if (ownEpsilon_)
174  {
175  deleteDemandDrivenData(epsilonPtr_);
176  ownEpsilon_ = false;
177  }
178  }
179 }
180 
181 
182 template<class CloudType>
184 {
186 
187  os.writeEntry("ownK", ownK_);
188  os.writeEntry("ownEpsilon", ownEpsilon_);
189 }
190 
191 
192 // ************************************************************************* //
tmp< volScalarField > kModel() const
Return the k field from the turbulence model.
bool ownEpsilon_
Local ownership of the epsilon field.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
bool movable() const noexcept
True if this is a non-null managed pointer with a unique ref-count.
Definition: tmpI.H:214
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
compressible::turbulenceModel & turb
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition: tmpI.H:221
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
virtual void write(Ostream &os) const
Write to os.
virtual ~DispersionRASModel()
Destructor.
bool ownK_
Local ownership of the k field.
Abstract base class for turbulence models (RAS, LES and laminar).
constexpr const char *const group
Group name for atomic constants.
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 Foam::string.
Definition: word.H:63
virtual void write(Ostream &os) const
Write.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
OBJstream os(runTime.globalPath()/outputName)
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:37
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
The table of contents (the keys) in sorted order.
Definition: HashTable.C:139
Base class for dispersion modelling.
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
Definition: tmpI.H:256
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Registry of regIOobjects.
Base class for particle dispersion models based on RAS turbulence.
void deleteDemandDrivenData(DataPtr &dataPtr)
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67