CollidingParcel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::CollidingParcel
26 
27 Description
28  Wrapper around parcel types to add collision modelling
29 
30 SourceFiles
31  CollidingParcelI.H
32  CollidingParcel.C
33  CollidingParcelIO.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef CollidingParcel_H
38 #define CollidingParcel_H
39 
40 #include "particle.H"
41 #include "demandDrivenEntry.H"
42 #include "CollisionRecordList.H"
43 #include "labelFieldIOField.H"
44 #include "vectorFieldIOField.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
54 
55 template<class ParcelType>
56 class CollidingParcel;
57 
58 // Forward declaration of friend functions
59 
60 template<class ParcelType>
61 Ostream& operator<<
62 (
63  Ostream&,
65 );
66 
67 /*---------------------------------------------------------------------------*\
68  Class CollidingParcelName Declaration
69 \*---------------------------------------------------------------------------*/
70 
72 
73 
74 /*---------------------------------------------------------------------------*\
75  Class CollidingParcel Declaration
76 \*---------------------------------------------------------------------------*/
77 
78 template<class ParcelType>
79 class CollidingParcel
80 :
81  public ParcelType,
82  public CollidingParcelName
83 {
84  // Private member data
85 
86  //- Size in bytes of the fields
87  static const std::size_t sizeofFields_;
88 
89 
90 public:
91 
92  //- Class to hold thermo particle constant properties
93  class constantProperties
94  :
95  public ParcelType::constantProperties
96  {
97 
98  // Private Data
99 
100  //- Young's modulus [N/m^2]
101  demandDrivenEntry<scalar> youngsModulus_;
102 
103  //- Poisson's ratio
104  demandDrivenEntry<scalar> poissonsRatio_;
105 
106 
107  public:
108 
109  // Constructors
110 
111  //- Null constructor
113 
114  //- Copy constructor
116 
117  //- Construct from dictionary
118  constantProperties(const dictionary& parentDict);
119 
120 
121  // Member Functions
122 
123  //- Return const access to Young's Modulus
124  inline scalar youngsModulus() const;
125 
126  //- Return const access to Poisson's ratio
127  inline scalar poissonsRatio() const;
128  };
129 
130 
131  //- Class to hold temporary data during tracking
132  class trackingData
133  :
134  public ParcelType::trackingData
135  {
136  public:
138  enum trackPart
139  {
140  tpVelocityHalfStep,
141  tpLinearTrack,
142  tpRotationalTrack
143  };
144 
145 
146  private:
147 
148  // Private Data
149 
150  //- Which part of the integration algorithm is taking place
151  trackPart part_;
152 
153 
154  public:
155 
156  // Constructors
157 
158  //- Construct from components
159  template <class TrackCloudType>
160  inline trackingData(const TrackCloudType& cloud);
161 
162 
163  // Member Functions
164 
165  //- Return the part of the tracking operation taking place
166  inline trackPart part() const;
167 
168  //- Access to the part of the tracking operation taking place
169  inline trackPart& part();
170  };
171 
172 
173 protected:
174 
175  // Protected data
176 
177  //- Force on particle due to collisions [N]
178  vector f_;
179 
180  //- Angular momentum of Parcel in global reference frame [kg m2/s]
182 
183  //- Torque on particle due to collisions in global
184  // reference frame [Nm]
185  vector torque_;
186 
187  //- Particle collision records
188  collisionRecordList collisionRecords_;
189 
190 
191 public:
192 
193  // Static Data Members
194 
195  //- String representation of properties
197  (
198  ParcelType,
199  " (fx fy fz)"
200  + " (angularMomentumx angularMomentumy angularMomentumz)"
201  + " (torquex torquey torquez)"
202  + " collisionRecordsPairAccessed"
203  + " collisionRecordsPairOrigProcOfOther"
204  + " collisionRecordsPairOrigIdOfOther"
205  + " (collisionRecordsPairData)"
206  + " collisionRecordsWallAccessed"
207  + " collisionRecordsWallPRel"
208  + " (collisionRecordsWallData)"
209  );
210 
211 
212  // Constructors
213 
214  //- Construct from mesh, coordinates and topology
215  // Other properties initialised as null
216  inline CollidingParcel
217  (
218  const polyMesh& mesh,
219  const barycentric& coordinates,
220  const label celli,
221  const label tetFacei,
222  const label tetPti
223  );
224 
225  //- Construct from a position and a cell, searching for the rest of the
226  // required topology. Other properties are initialised as null.
227  inline CollidingParcel
228  (
229  const polyMesh& mesh,
230  const vector& position,
231  const label celli
232  );
233 
234  //- Construct from Istream
236  (
237  const polyMesh& mesh,
238  Istream& is,
239  bool readFields = true
240  );
241 
242  //- Construct as a copy
244 
245  //- Construct as a copy
247 
248  //- Construct and return a (basic particle) clone
249  virtual autoPtr<particle> clone() const
250  {
251  return autoPtr<particle>(new CollidingParcel(*this));
252  }
253 
254  //- Construct and return a (basic particle) clone
255  virtual autoPtr<particle> clone(const polyMesh& mesh) const
256  {
257  return autoPtr<particle>(new CollidingParcel(*this, mesh));
258  }
259 
260  //- Factory class to read-construct particles used for
261  // parallel transfer
262  class iNew
263  {
264  const polyMesh& mesh_;
265 
266  public:
268  iNew(const polyMesh& mesh)
269  :
270  mesh_(mesh)
271  {}
273  autoPtr<CollidingParcel<ParcelType>> operator()(Istream& is) const
274  {
276  (
277  new CollidingParcel<ParcelType>(mesh_, is, true)
278  );
279  }
280  };
281 
282 
283  // Member Functions
284 
285  // Access
286 
287  //- Return const access to force
288  inline const vector& f() const;
289 
290  //- Return const access to angular momentum
291  inline const vector& angularMomentum() const;
292 
293  //- Return const access to torque
294  inline const vector& torque() const;
295 
296  //- Return const access to the collision records
297  inline const collisionRecordList& collisionRecords() const;
298 
299  //- Return access to force
300  inline vector& f();
301 
302  //- Return access to angular momentum
303  inline vector& angularMomentum();
304 
305  //- Return access to torque
306  inline vector& torque();
307 
308  //- Return access to collision records
309  inline collisionRecordList& collisionRecords();
310 
311  //- Particle angular velocity
312  inline vector omega() const;
313 
314 
315  // Tracking
316 
317  //- Move the parcel
318  template<class TrackCloudType>
319  bool move
320  (
321  TrackCloudType& cloud,
322  trackingData& td,
323  const scalar trackTime
324  );
325 
326  //- Transform the physical properties of the particle
327  // according to the given transformation tensor
328  virtual void transformProperties(const transformer&);
329 
330 
331  // I-O
332 
333  //- Read
334  template<class CloudType>
335  static void readFields(CloudType& c);
336 
337  //- Write
338  template<class CloudType>
339  static void writeFields(const CloudType& c);
340 
341 
342  // Ostream Operator
343 
344  friend Ostream& operator<< <ParcelType>
345  (
346  Ostream&,
348  );
349 };
350 
351 
352 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
353 
354 } // End namespace Foam
355 
356 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
357 
358 #include "CollidingParcelI.H"
360 
361 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
362 
363 #ifdef NoRepository
364  #include "CollidingParcel.C"
365 #endif
366 
367 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
368 
369 #endif
370 
371 // ************************************************************************* //
vector f_
Force on particle due to collisions [N].
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space...
Definition: transformer.H:83
CollisionRecordList< vector, vector > collisionRecordList
vectorFieldCompactIOField pairDataFieldCompactIOField
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Wrapper around parcel types to add collision modelling.
Class to hold thermo particle constant properties.
Class to hold temporary data during tracking.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
const vector & f() const
Return const access to force.
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:753
static void writeFields(const CloudType &c)
Write.
vectorFieldCompactIOField wallDataFieldCompactIOField
scalar poissonsRatio() const
Return const access to Poisson&#39;s ratio.
virtual void transformProperties(const transformer &)
Transform the physical properties of the particle.
scalar youngsModulus() const
Return const access to Young&#39;s Modulus.
collisionRecordList collisionRecords_
Particle collision records.
fvMesh & mesh
const dimensionedScalar c
Speed of light in a vacuum.
AddToPropertyList(ParcelType, " (fx fy fz)"+" (angularMomentumx angularMomentumy angularMomentumz)"+" (torquex torquey torquez)"+" collisionRecordsPairAccessed"+" collisionRecordsPairOrigProcOfOther"+" collisionRecordsPairOrigIdOfOther"+" (collisionRecordsPairData)"+" collisionRecordsWallAccessed"+" collisionRecordsWallPRel"+" (collisionRecordsWallData)")
String representation of properties.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
TemplateName(FvFaceCellWave)
bool move(TrackCloudType &cloud, trackingData &td, const scalar trackTime)
Move the parcel.
const vector & torque() const
Return const access to torque.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
const vector & angularMomentum() const
Return const access to angular momentum.
vector angularMomentum_
Angular momentum of Parcel in global reference frame [kg m2/s].
vector torque_
Torque on particle due to collisions in global.
vector omega() const
Particle angular velocity.
A Field of objects of type <T> with automated input and output using a compact storage. Behaves like IOField except when binary output in case it writes a CompactListList.
CollidingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
Factory class to read-construct particles used for.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
volScalarField & p
const collisionRecordList & collisionRecords() const
Return const access to the collision records.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
Namespace for OpenFOAM.
static void readFields(CloudType &c)
Read.