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 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 Class
28  Foam::DSMCParcel
29 
30 Description
31  DSMC parcel class
32 
33 SourceFiles
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 "autoPtr.H"
46 #include "contiguous.H"
47 #include "DSMCCloud.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 
54 template<class ParcelType>
55 class DSMCParcel;
56 
57 // Forward declaration of friend functions
58 
59 template<class ParcelType>
60 Ostream& operator<<
61 (
62  Ostream&,
64 );
65 
66 /*---------------------------------------------------------------------------*\
67  Class DSMCParcel Declaration
68 \*---------------------------------------------------------------------------*/
69 
70 template<class ParcelType>
71 class DSMCParcel
72 :
73  public ParcelType
74 {
75 public:
76 
77  //- Size in bytes of the fields
78  static const std::size_t sizeofFields;
79 
80 
81  //- Class to hold DSMC particle constant properties
82  class constantProperties
83  {
84  // Private data
85 
86  //- Particle mass [kg] (constant)
87  scalar mass_;
88 
89  //- Particle hard sphere diameter [m] (constant)
90  scalar d_;
91 
92  //- Internal degrees of freedom
93  direction internalDegreesOfFreedom_;
94 
95  //- Viscosity index
96  scalar omega_;
97 
98 
99  public:
100 
101  // Constructors
102 
103  //- Null constructor, allows List of constantProperties to be
104  // created before the contents is initialised
105  inline constantProperties();
106 
107  //- Constructor from dictionary
108  inline constantProperties(const dictionary& dict);
109 
110 
111  // Member functions
112 
113  //- Return const access to the particle mass [kg]
114  inline scalar mass() const;
115 
116  //- Return const access to the hard sphere diameter [m]
117  inline scalar d() const;
118 
119  //- Return the reference total collision cross section
120  inline scalar sigmaT() const;
121 
122  //- Return the internalDegreesOfFreedom
123  inline direction internalDegreesOfFreedom() const;
124 
125  //- Return the viscosity index
126  inline scalar omega() const;
127 
128  };
129 
130 
131  //- Use base tracking data
132  typedef typename ParcelType::trackingData trackingData;
133 
134 
135 protected:
136 
137  // Protected member data
138 
139  // Parcel properties
140 
141  //- Velocity of Parcel [m/s]
142  vector U_;
143 
144  //- Internal energy of the Parcel, covering all non-translational
145  // degrees of freedom [J]
146  scalar Ei_;
147 
148  //- Parcel type id
149  label typeId_;
150 
151 
152 public:
153 
154  //- Runtime type information
155  TypeName("DSMCParcel");
156 
157  friend class Cloud<ParcelType>;
158 
159 
160  // Constructors
161 
162  //- Construct from components
163  inline DSMCParcel
164  (
165  const polyMesh& mesh,
166  const barycentric& coordinates,
167  const label celli,
168  const label tetFacei,
169  const label tetPti,
170  const vector& U,
171  const scalar Ei,
172  const label typeId
173  );
174 
175  //- Construct from a position and a cell, searching for the rest of the
176  // required topology
177  inline DSMCParcel
178  (
179  const polyMesh& mesh,
180  const vector& position,
181  const label celli,
182  const vector& U,
183  const scalar Ei,
184  const label typeId
185  );
186 
187  //- Construct from Istream
188  DSMCParcel
189  (
190  const polyMesh& mesh,
191  Istream& is,
192  bool readFields = true,
193  bool newFormat = true
194  );
195 
196  //- Construct and return a clone
197  virtual autoPtr<particle> clone() const
198  {
199  return autoPtr<particle>(new DSMCParcel<ParcelType>(*this));
200  }
201 
202  //- Factory class to read-construct particles used for
203  // parallel transfer
204  class iNew
205  {
206  const polyMesh& mesh_;
207 
208  public:
209 
210  iNew(const polyMesh& mesh)
211  :
212  mesh_(mesh)
213  {}
214 
216  {
218  (
219  new DSMCParcel<ParcelType>(mesh_, is, true)
220  );
221  }
222  };
223 
224 
225  // Member Functions
226 
227  // Access
228 
229  //- Return type id
230  inline label typeId() const;
231 
232  //- Return const access to velocity
233  inline const vector& U() const;
234 
235  //- Return const access to internal energy
236  inline scalar Ei() const;
237 
238 
239  // Edit
240 
241  //- Return access to velocity
242  inline vector& U();
243 
244  //- Return access to internal energy
245  inline scalar& Ei();
246 
247 
248  // Main calculation loop
249 
250  // Tracking
251 
252  //- Move the parcel
253  template<class TrackCloudType>
254  bool move
255  (
256  TrackCloudType& cloud,
257  trackingData& td,
258  const scalar trackTime
259  );
260 
261 
262  // Patch interactions
263 
264  //- Overridable function to handle the particle hitting a patch
265  // Executed before other patch-hitting functions
266  template<class TrackCloudType>
267  bool hitPatch(TrackCloudType&, trackingData&);
268 
269  //- Overridable function to handle the particle hitting a
270  // processorPatch
271  template<class TrackCloudType>
272  void hitProcessorPatch(TrackCloudType&, trackingData&);
273 
274  //- Overridable function to handle the particle hitting a wallPatch
275  template<class TrackCloudType>
276  void hitWallPatch(TrackCloudType&, trackingData&);
277 
278  //- Transform the physical properties of the particle
279  // according to the given transformation tensor
280  virtual void transformProperties(const tensor& T);
281 
282  //- Transform the physical properties of the particle
283  // according to the given separation vector
284  virtual void transformProperties(const vector& separation);
285 
286 
287  // I-O
288 
289  static void readFields(Cloud<DSMCParcel<ParcelType>>& c);
290 
291  static void writeFields(const Cloud<DSMCParcel<ParcelType>>& c);
292 
293 
294  // Ostream Operator
295 
296  friend Ostream& operator<< <ParcelType>
297  (
298  Ostream&,
300  );
301 };
302 
303 
304 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
305 
306 } // End namespace Foam
307 
308 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
309 
310 #include "DSMCParcelI.H"
311 
312 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
313 
314 #ifdef NoRepository
315  #include "DSMCParcel.C"
316 #endif
317 
318 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319 
320 #endif
321 
322 // ************************************************************************* //
Foam::DSMCParcel::U
const vector & U() const
Return const access to velocity.
Definition: DSMCParcelI.H:143
Foam::Tensor< scalar >
Foam::DSMCParcel
DSMC parcel class.
Definition: DSMCParcel.H:54
DSMCParcelI.H
Foam::DSMCParcel::DSMCParcel
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:55
Foam::DSMCParcel::hitPatch
bool hitPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a patch.
Definition: DSMCParcel.C:77
Foam::DSMCParcel::constantProperties
Class to hold DSMC particle constant properties.
Definition: DSMCParcel.H:81
Foam::DSMCParcel::hitProcessorPatch
void hitProcessorPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
Definition: DSMCParcel.C:86
Foam::DSMCParcel::TypeName
TypeName("DSMCParcel")
Runtime type information.
Foam::DSMCParcel::iNew
Factory class to read-construct particles used for.
Definition: DSMCParcel.H:203
Foam::DSMCParcel::typeId
label typeId() const
Return type id.
Definition: DSMCParcelI.H:136
DSMCParcel.C
Foam::DSMCParcel::constantProperties::constantProperties
constantProperties()
Null constructor, allows List of constantProperties to be.
Definition: DSMCParcelI.H:33
Foam::DSMCParcel::iNew::operator()
autoPtr< DSMCParcel< ParcelType > > operator()(Istream &is) const
Definition: DSMCParcel.H:214
Foam::DSMCParcel::typeId_
label typeId_
Parcel type id.
Definition: DSMCParcel.H:148
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::DSMCParcel::readFields
static void readFields(Cloud< DSMCParcel< ParcelType >> &c)
Definition: DSMCParcelIO.C:88
Foam::DSMCParcel::constantProperties::sigmaT
scalar sigmaT() const
Return the reference total collision cross section.
Definition: DSMCParcelI.H:110
Foam::DSMCParcel::constantProperties::internalDegreesOfFreedom
direction internalDegreesOfFreedom() const
Return the internalDegreesOfFreedom.
Definition: DSMCParcelI.H:118
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
coordinates
PtrList< coordinateSystem > coordinates(solidRegions.size())
IOstream.H
Foam::DSMCParcel::constantProperties::omega
scalar omega() const
Return the viscosity index.
Definition: DSMCParcelI.H:127
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::Barycentric< scalar >
Foam::DSMCParcel::U_
vector U_
Velocity of Parcel [m/s].
Definition: DSMCParcel.H:141
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::DSMCParcel::clone
virtual autoPtr< particle > clone() const
Construct and return a clone.
Definition: DSMCParcel.H:196
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::DSMCParcel::Ei_
scalar Ei_
Internal energy of the Parcel, covering all non-translational.
Definition: DSMCParcel.H:145
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::DSMCParcel::iNew::iNew
iNew(const polyMesh &mesh)
Definition: DSMCParcel.H:209
Foam::DSMCParcel::constantProperties::mass
scalar mass() const
Return const access to the particle mass [kg].
Definition: DSMCParcelI.H:95
Foam::DSMCParcel::move
bool move(TrackCloudType &cloud, trackingData &td, const scalar trackTime)
Move the parcel.
Definition: DSMCParcel.C:36
Foam::DSMCParcel::hitWallPatch
void hitWallPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
Definition: DSMCParcel.C:98
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:57
Foam::Vector< scalar >
contiguous.H
Foam::DSMCParcel::writeFields
static void writeFields(const Cloud< DSMCParcel< ParcelType >> &c)
Definition: DSMCParcelIO.C:120
Foam::Cloud
Base cloud calls templated on particle type.
Definition: Cloud.H:55
DSMCCloud.H
Foam::DSMCParcel::trackingData
ParcelType::trackingData trackingData
Use base tracking data.
Definition: DSMCParcel.H:131
Foam::direction
uint8_t direction
Definition: direction.H:52
particle.H
Foam::DSMCParcel::constantProperties::d
scalar d() const
Return const access to the hard sphere diameter [m].
Definition: DSMCParcelI.H:102
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::DSMCParcel::transformProperties
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition: DSMCParcel.C:188
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::DSMCParcel::Ei
scalar Ei() const
Return const access to internal energy.
Definition: DSMCParcelI.H:150
Foam::DSMCParcel::sizeofFields
static const std::size_t sizeofFields
Size in bytes of the fields.
Definition: DSMCParcel.H:77
autoPtr.H