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-------------------------------------------------------------------------------
10License
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
34template<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
62template<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
92template<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
107template<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
126template<class CloudType>
128{
129 cacheFields(false);
130}
131
132
133// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134
135template<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 {
169 ownK_ = false;
170 }
171 if (ownEpsilon_ && epsilonPtr_)
172 {
173 deleteDemandDrivenData(epsilonPtr_);
174 ownEpsilon_ = false;
175 }
176 }
177}
178
179
180template<class CloudType>
182{
184
185 os.writeEntry("ownK", ownK_);
186 os.writeEntry("ownEpsilon", ownEpsilon_);
187}
188
189
190// ************************************************************************* //
compressible::turbulenceModel & turb
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
Base class for dispersion modelling.
Base class for particle dispersion models based on RAS turbulence.
tmp< volScalarField > kModel() const
Return the k field from the turbulence model.
virtual void cacheFields(const bool store)
Cache carrier fields.
tmp< volScalarField > epsilonModel() const
Return the epsilon field from the turbulence model.
bool ownEpsilon_
Take ownership of the epsilon field.
virtual ~DispersionRASModel()
Destructor.
bool ownK_
Take ownership of the k field.
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:137
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool write()
Write the output fields.
Registry of regIOobjects.
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
A class for managing temporary objects.
Definition: tmp.H:65
bool isTmp() const noexcept
Identical to is_pointer()
Definition: tmp.H:295
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
Definition: tmpI.H:255
Abstract base class for turbulence models (RAS, LES and laminar).
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
void deleteDemandDrivenData(DataPtr &dataPtr)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53