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"
29 #include "demandDrivenData.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 :
100  kPtr_(nullptr),
101  ownK_(false),
102  epsilonPtr_(nullptr),
103  ownEpsilon_(false)
104 {}
105 
106 
107 template<class CloudType>
109 (
111 )
112 :
114  kPtr_(dm.kPtr_),
115  ownK_(dm.ownK_),
116  epsilonPtr_(dm.epsilonPtr_),
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.isTmp())
142  {
143  kPtr_ = tk.ptr();
144  ownK_ = true;
145  }
146  else
147  {
148  kPtr_ = &tk();
149  ownK_ = false;
150  }
151 
152  tmp<volScalarField> tepsilon = this->epsilonModel();
153  if (tepsilon.isTmp())
154  {
155  epsilonPtr_ = tepsilon.ptr();
156  ownEpsilon_ = true;
157  }
158  else
159  {
160  epsilonPtr_ = &tepsilon();
161  ownEpsilon_ = false;
162  }
163  }
164  else
165  {
166  if (ownK_ && kPtr_)
167  {
168  deleteDemandDrivenData(kPtr_);
169  ownK_ = false;
170  }
171  if (ownEpsilon_ && epsilonPtr_)
172  {
173  deleteDemandDrivenData(epsilonPtr_);
174  ownEpsilon_ = false;
175  }
176  }
177 }
178 
179 
180 template<class CloudType>
182 {
184 
185  os.writeEntry("ownK", ownK_);
186  os.writeEntry("ownEpsilon", ownEpsilon_);
187 }
188 
189 
190 // ************************************************************************* //
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::constant::atomic::group
constexpr const char *const group
Group name for atomic constants.
Definition: atomicConstants.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::DispersionModel
Base class for dispersion modelling.
Definition: KinematicCloud.H:86
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::DispersionRASModel::epsilonPtr_
const volScalarField * epsilonPtr_
Turbulence epsilon.
Definition: DispersionRASModel.H:66
Foam::DispersionRASModel::write
virtual void write(Ostream &os) const
Write.
Definition: DispersionRASModel.C:181
Foam::DispersionRASModel::DispersionRASModel
DispersionRASModel(const dictionary &dict, CloudType &owner)
Construct from components.
Definition: DispersionRASModel.C:94
Foam::tmp::isTmp
bool isTmp() const noexcept
Identical to is_pointer()
Definition: tmp.H:181
Foam::DispersionRASModel::kPtr_
const volScalarField * kPtr_
Turbulence k.
Definition: DispersionRASModel.H:60
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::DispersionRASModel::cacheFields
virtual void cacheFields(const bool store)
Cache carrier fields.
Definition: DispersionRASModel.C:136
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::DispersionRASModel::kModel
tmp< volScalarField > kModel() const
Return the k field from the turbulence model.
Definition: DispersionRASModel.C:36
Foam::tmp::ptr
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
Definition: tmpI.H:259
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::HashTable::sortedToc
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:136
Foam::DispersionRASModel::ownEpsilon_
bool ownEpsilon_
Take ownership of the epsilon field.
Definition: DispersionRASModel.H:69
U
U
Definition: pEqn.H:72
Foam::DispersionRASModel::~DispersionRASModel
virtual ~DispersionRASModel()
Destructor.
Definition: DispersionRASModel.C:127
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::objectRegistry::findObject
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
Definition: objectRegistryTemplates.C:401
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::DispersionRASModel::ownK_
bool ownK_
Take ownership of the k field.
Definition: DispersionRASModel.H:63
DispersionRASModel.H
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:35
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
turb
compressible::turbulenceModel & turb
Definition: setRegionFluidFields.H:10
turbulenceModel.H
Foam::DispersionRASModel
Base class for particle dispersion models based on RAS turbulence.
Definition: DispersionRASModel.H:49
Foam::DispersionRASModel::epsilonModel
tmp< volScalarField > epsilonModel() const
Return the epsilon field from the turbulence model.
Definition: DispersionRASModel.C:64