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-2022 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:
137 
138  enum trackPart
139  {
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
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  const label facei
224  );
225 
226  //- Construct from a position and a cell, searching for the rest of the
227  // required topology. Other properties are initialised as null.
228  inline CollidingParcel
229  (
230  const polyMesh& mesh,
231  const vector& position,
232  const label celli
233  );
234 
235  //- Construct from Istream
236  CollidingParcel(Istream& is, bool readFields = true);
237 
238  //- Construct as a copy
240 
241  //- Construct and return a clone
242  virtual autoPtr<particle> clone() const
243  {
244  return autoPtr<particle>(new CollidingParcel(*this));
245  }
246 
247  //- Construct from Istream and return
249  {
251  }
252 
253 
254  // Member Functions
255 
256  // Access
257 
258  //- Return const access to force
259  inline const vector& f() const;
260 
261  //- Return const access to angular momentum
262  inline const vector& angularMomentum() const;
263 
264  //- Return const access to torque
265  inline const vector& torque() const;
266 
267  //- Return const access to the collision records
268  inline const collisionRecordList& collisionRecords() const;
269 
270  //- Return access to force
271  inline vector& f();
272 
273  //- Return access to angular momentum
274  inline vector& angularMomentum();
275 
276  //- Return access to torque
277  inline vector& torque();
278 
279  //- Return access to collision records
281 
282  //- Particle angular velocity
283  inline vector omega() const;
284 
285 
286  // Tracking
287 
288  //- Move the parcel
289  template<class TrackCloudType>
290  bool move(TrackCloudType& cloud, trackingData& td);
291 
292  //- Transform the physical properties of the particle
293  // according to the given transformation tensor
294  virtual void transformProperties(const transformer&);
295 
296 
297  // I-O
298 
299  //- Read
300  template<class CloudType>
301  static void readFields(CloudType& c);
302 
303  //- Write
304  template<class CloudType>
305  static void writeFields(const CloudType& c);
306 
307 
308  // Ostream Operator
309 
310  friend Ostream& operator<< <ParcelType>
311  (
312  Ostream&,
314  );
315 };
316 
317 
318 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319 
320 } // End namespace Foam
321 
322 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
323 
324 #include "CollidingParcelI.H"
326 
327 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
328 
329 #ifdef NoRepository
330  #include "CollidingParcel.C"
331 #endif
332 
333 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
334 
335 #endif
336 
337 // ************************************************************************* //
Class to hold thermo particle constant properties.
scalar youngsModulus() const
Return const access to Young's Modulus.
scalar poissonsRatio() const
Return const access to Poisson's ratio.
Class to hold temporary data during tracking.
trackingData(const TrackCloudType &cloud)
Construct from components.
trackPart part() const
Return the part of the tracking operation taking place.
Wrapper around parcel types to add collision modelling.
bool move(TrackCloudType &cloud, trackingData &td)
Move the parcel.
const vector & angularMomentum() const
Return const access to angular momentum.
static autoPtr< CollidingParcel > New(Istream &is)
Construct from Istream and return.
vector f_
Force on particle due to collisions [N].
AddToPropertyList(ParcelType, " (fx fy fz)"+" (angularMomentumx angularMomentumy angularMomentumz)"+" (torquex torquey torquez)"+" collisionRecordsPairAccessed"+" collisionRecordsPairOrigProcOfOther"+" collisionRecordsPairOrigIdOfOther"+" (collisionRecordsPairData)"+" collisionRecordsWallAccessed"+" collisionRecordsWallPRel"+" (collisionRecordsWallData)")
String representation of properties.
CollidingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label facei)
Construct from mesh, coordinates and topology.
vector omega() const
Particle angular velocity.
const vector & f() const
Return const access to force.
virtual autoPtr< particle > clone() const
Construct and return a clone.
collisionRecordList collisionRecords_
Particle collision records.
static void writeFields(const CloudType &c)
Write.
const vector & torque() const
Return const access to torque.
vector angularMomentum_
Angular momentum of Parcel in global reference frame [kg m2/s].
const collisionRecordList & collisionRecords() const
Return const access to the collision records.
virtual void transformProperties(const transformer &)
Transform the physical properties of the particle.
static void readFields(CloudType &c)
Read.
vector torque_
Torque on particle due to collisions in global.
A Field of objects of type <Type> with automated input and output using a compact storage....
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:79
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:57
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space.
Definition: transformer.H:84
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
vectorFieldCompactIOField pairDataFieldCompactIOField
CollisionRecordList< vector, vector > collisionRecordList
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
TemplateName(FvFaceCellWave)
vectorFieldCompactIOField wallDataFieldCompactIOField
volScalarField & p