Loading...
Searching...
No Matches
DSMCParcel.H
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-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::DSMCParcel
29
30Description
31 DSMC parcel class
32
33SourceFiles
34 DSMCParcelI.H
35 DSMCParcel.C
36 DSMCParcelIO.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef DSMCParcel_H
41#define DSMCParcel_H
42
43#include "particle.H"
44#include "IOstream.H"
45#include "contiguous.H"
46#include "DSMCCloud.H"
47
48// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49
50namespace Foam
51{
52
53template<class ParcelType>
55
56// Forward declaration of friend functions
57
58template<class ParcelType>
59Ostream& operator<<
60(
61 Ostream&,
63);
64
65/*---------------------------------------------------------------------------*\
66 Class DSMCParcel Declaration
67\*---------------------------------------------------------------------------*/
68
69template<class ParcelType>
70class DSMCParcel
71:
72 public ParcelType
73{
74public:
75
76 //- Size in bytes of the fields
77 static const std::size_t sizeofFields;
78
79
80 //- Class to hold DSMC particle constant properties
82 {
83 // Private data
84
85 //- Particle mass [kg] (constant)
86 scalar mass_;
87
88 //- Particle hard sphere diameter [m] (constant)
89 scalar d_;
90
91 //- Internal degrees of freedom
92 direction internalDegreesOfFreedom_;
93
94 //- Viscosity index
95 scalar omega_;
96
97
98 public:
99
100 // Constructors
101
102 //- Null constructor, allows List of constantProperties to be
103 // created before the contents is initialised
104 inline constantProperties();
105
106 //- Constructor from dictionary
107 inline constantProperties(const dictionary& dict);
108
109
110 // Member functions
111
112 //- Return const access to the particle mass [kg]
113 inline scalar mass() const;
114
115 //- Return const access to the hard sphere diameter [m]
116 inline scalar d() const;
117
118 //- Return the reference total collision cross section
119 inline scalar sigmaT() const;
120
121 //- Return the internalDegreesOfFreedom
122 inline direction internalDegreesOfFreedom() const;
123
124 //- Return the viscosity index
125 inline scalar omega() const;
126
127 };
128
129
130 //- Use base tracking data
131 typedef typename ParcelType::trackingData trackingData;
132
133
134protected:
135
136 // Protected member data
137
138 // Parcel properties
139
140 //- Velocity of Parcel [m/s]
141 vector U_;
142
143 //- Internal energy of the Parcel, covering all non-translational
144 // degrees of freedom [J]
145 scalar Ei_;
146
147 //- Parcel type id
148 label typeId_;
149
150
151public:
152
153 //- Runtime type information
154 TypeName("DSMCParcel");
156 friend class Cloud<ParcelType>;
157
158
159 // Constructors
160
161 //- Construct from components
162 inline DSMCParcel
163 (
164 const polyMesh& mesh,
166 const label celli,
167 const label tetFacei,
168 const label tetPti,
169 const vector& U,
170 const scalar Ei,
171 const label typeId
172 );
173
174 //- Construct from a position and a cell, searching for the rest of the
175 // required topology
176 inline DSMCParcel
177 (
178 const polyMesh& mesh,
180 const label celli,
181 const vector& U,
182 const scalar Ei,
183 const label typeId
184 );
185
186 //- Construct from Istream
188 (
189 const polyMesh& mesh,
190 Istream& is,
191 bool readFields = true,
192 bool newFormat = true
193 );
194
195 //- Return a clone
196 virtual autoPtr<particle> clone() const
198 return particle::Clone(*this);
199 }
200
201 //- Factory class to read-construct particles (for parallel transfer)
202 class iNew
203 {
204 const polyMesh& mesh_;
205
206 public:
207
208 iNew(const polyMesh& mesh)
209 :
210 mesh_(mesh)
211 {}
212
213 autoPtr<DSMCParcel<ParcelType>> operator()(Istream& is) const
216 (
217 new DSMCParcel<ParcelType>(mesh_, is, true)
218 );
219 }
220 };
221
222
223 // Member Functions
224
225 // Access
226
227 //- Return type id
228 inline label typeId() const;
229
230 //- Return const access to velocity
231 inline const vector& U() const;
232
233 //- Return const access to internal energy
234 inline scalar Ei() const;
235
236
237 // Edit
239 //- Return access to velocity
240 inline vector& U();
241
242 //- Return access to internal energy
243 inline scalar& Ei();
244
245
246 // Main calculation loop
247
248 // Tracking
249
250 //- Move the parcel
251 template<class TrackCloudType>
252 bool move
253 (
254 TrackCloudType& cloud,
256 const scalar trackTime
257 );
258
259
260 // Patch interactions
261
262 //- Overridable function to handle the particle hitting a patch
263 // Executed before other patch-hitting functions
264 template<class TrackCloudType>
265 bool hitPatch(TrackCloudType&, trackingData&);
266
267 //- Overridable function to handle the particle hitting a
268 // processorPatch
269 template<class TrackCloudType>
270 void hitProcessorPatch(TrackCloudType&, trackingData&);
271
272 //- Overridable function to handle the particle hitting a wallPatch
273 template<class TrackCloudType>
274 void hitWallPatch(TrackCloudType&, trackingData&);
275
276 //- Transform the physical properties of the particle
277 // according to the given transformation tensor
278 virtual void transformProperties(const tensor& T);
280 //- Transform the physical properties of the particle
281 // according to the given separation vector
282 virtual void transformProperties(const vector& separation);
283
285 // I-O
286
287 static void readFields(Cloud<DSMCParcel<ParcelType>>& c);
288
289 static void writeFields(const Cloud<DSMCParcel<ParcelType>>& c);
290
291
292 // Ostream Operator
293
294 friend Ostream& operator<< <ParcelType>
295 (
296 Ostream&,
298 );
299};
300
301
302// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303
304} // End namespace Foam
305
306// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
307
309
310// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
311
312#ifdef NoRepository
313 #include "DSMCParcel.C"
314#endif
315
316// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
317
318#endif
319
320// ************************************************************************* //
Base cloud calls templated on particle type.
Definition Cloud.H:64
scalar d() const
Return const access to the hard sphere diameter [m].
Definition DSMCParcelI.H:95
constantProperties()
Null constructor, allows List of constantProperties to be.
Definition DSMCParcelI.H:26
scalar omega() const
Return the viscosity index.
scalar mass() const
Return const access to the particle mass [kg].
Definition DSMCParcelI.H:88
scalar sigmaT() const
Return the reference total collision cross section.
direction internalDegreesOfFreedom() const
Return the internalDegreesOfFreedom.
autoPtr< DSMCParcel< ParcelType > > operator()(Istream &is) const
Definition DSMCParcel.H:257
iNew(const polyMesh &mesh)
Definition DSMCParcel.H:252
DSMC parcel class.
Definition DSMCParcel.H:68
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition DSMCParcel.C:181
const vector & U() const
virtual void transformProperties(const vector &separation)
Transform the physical properties of the particle.
Definition DSMCParcel.C:190
void hitProcessorPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
Definition DSMCParcel.C:79
virtual autoPtr< particle > clone() const
Return a clone.
Definition DSMCParcel.H:238
static const std::size_t sizeofFields
Definition DSMCParcel.H:74
bool hitPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a patch.
Definition DSMCParcel.C:70
particle::trackingData trackingData
Definition DSMCParcel.H:155
static void readFields(Cloud< DSMCParcel< particle > > &c)
static void writeFields(const Cloud< DSMCParcel< ParcelType > > &c)
DSMCParcel(const polyMesh &mesh, const vector &position, const label celli, const vector &U, const scalar Ei, const label typeId)
Construct from a position and a cell, searching for the rest of the.
Definition DSMCParcelI.H:68
DSMCParcel(const polyMesh &mesh, Istream &is, bool readFields=true, bool newFormat=true)
Construct from Istream.
bool move(TrackCloudType &cloud, trackingData &td, const scalar trackTime)
Move the parcel.
Definition DSMCParcel.C:29
vector & U()
Return access to velocity.
DSMCParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const vector &U, const scalar Ei, const label typeId)
Construct from components.
Definition DSMCParcelI.H:48
void hitWallPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
Definition DSMCParcel.C:91
scalar & Ei()
Return access to internal energy.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
vector position() const
Return current particle position.
Definition particleI.H:283
static autoPtr< particle > Clone(const Derived &p)
Clone a particle.
Definition particle.H:552
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition particleI.H:110
const barycentric & coordinates() const noexcept
Return current particle coordinates.
Definition particleI.H:116
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
dynamicFvMesh & mesh
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace for OpenFOAM.
Tensor< scalar > tensor
Definition symmTensor.H:57
uint8_t direction
Definition direction.H:49
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
Vector< scalar > vector
Definition vector.H:57
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68