CollidingParcel.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 | Copyright (C) 2011-2016 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 kinematic 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 "CollisionRecordList.H"
42 #include "labelFieldIOField.H"
43 #include "vectorFieldIOField.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
53 
54 template<class ParcelType>
55 class CollidingParcel;
56 
57 // Forward declaration of friend functions
58 
59 template<class ParcelType>
60 Ostream& operator<<
61 (
62  Ostream&,
64 );
65 
66 /*---------------------------------------------------------------------------*\
67  Class CollidingParcel Declaration
68 \*---------------------------------------------------------------------------*/
69 
70 template<class ParcelType>
71 class CollidingParcel
72 :
73  public ParcelType
74 {
75  // Private member data
76 
77  //- Size in bytes of the fields
78  static const std::size_t sizeofFields_;
79 
80 
81 public:
82 
83  //- Class to hold thermo particle constant properties
84  class constantProperties
85  :
86  public ParcelType::constantProperties
87  {
88 
89  // Private data
90 
91  //- Young's modulus [N/m2]
92  demandDrivenEntry<scalar> youngsModulus_;
93 
94  //- Poisson's ratio
95  demandDrivenEntry<scalar> poissonsRatio_;
96 
97 
98  public:
99 
100  // Constructors
101 
102  //- Null constructor
104 
105  //- Copy constructor
107 
108  //- Construct from dictionary
109  constantProperties(const dictionary& parentDict);
110 
111 
112  // Member functions
113 
114  //- Return const access to Young's Modulus
115  inline scalar youngsModulus() const;
116 
117  //- Return const access to Poisson's ratio
118  inline scalar poissonsRatio() const;
119  };
120 
121 
122 protected:
123 
124  // Protected data
125 
126  //- Force on particle due to collisions [N]
127  vector f_;
128 
129  //- Angular momentum of Parcel in global reference frame [kg m2/s]
131 
132  //- Torque on particle due to collisions in global
133  // reference frame [Nm]
134  vector torque_;
135 
136  //- Particle collision records
137  collisionRecordList collisionRecords_;
138 
139 
140 public:
141 
142  // Static data members
143 
144  //- Runtime type information
145  TypeName("CollidingParcel");
146 
147  //- String representation of properties
149  (
150  ParcelType,
151  " (fx fy fz)"
152  + " (angularMomentumx angularMomentumy angularMomentumz)"
153  + " (torquex torquey torquez)"
154  + " collisionRecordsPairAccessed"
155  + " collisionRecordsPairOrigProcOfOther"
156  + " collisionRecordsPairOrigIdOfOther"
157  + " (collisionRecordsPairData)"
158  + " collisionRecordsWallAccessed"
159  + " collisionRecordsWallPRel"
160  + " (collisionRecordsWallData)"
161  );
162 
163 
164  // Constructors
165 
166  //- Construct from owner, position, and cloud owner
167  // Other properties initialised as null
168  inline CollidingParcel
169  (
170  const polyMesh& mesh,
171  const vector& position,
172  const label celli,
173  const label tetFacei,
174  const label tetPtI
175  );
176 
177  //- Construct from components
178  inline CollidingParcel
179  (
180  const polyMesh& mesh,
181  const vector& position,
182  const label celli,
183  const label tetFacei,
184  const label tetPtI,
185  const label typeId,
186  const scalar nParticle0,
187  const scalar d0,
188  const scalar dTarget0,
189  const vector& U0,
190  const vector& f0,
191  const vector& angularMomentum0,
192  const vector& torque0,
193  const typename ParcelType::constantProperties& constProps
194  );
195 
196  //- Construct from Istream
198  (
199  const polyMesh& mesh,
200  Istream& is,
201  bool readFields = true
202  );
203 
204  //- Construct as a copy
206 
207  //- Construct as a copy
209 
210  //- Construct and return a (basic particle) clone
211  virtual autoPtr<particle> clone() const
212  {
213  return autoPtr<particle>(new CollidingParcel(*this));
214  }
215 
216  //- Construct and return a (basic particle) clone
217  virtual autoPtr<particle> clone(const polyMesh& mesh) const
218  {
219  return autoPtr<particle>(new CollidingParcel(*this, mesh));
220  }
221 
222  //- Factory class to read-construct particles used for
223  // parallel transfer
224  class iNew
225  {
226  const polyMesh& mesh_;
227 
228  public:
230  iNew(const polyMesh& mesh)
231  :
232  mesh_(mesh)
233  {}
235  autoPtr<CollidingParcel<ParcelType>> operator()(Istream& is) const
236  {
238  (
239  new CollidingParcel<ParcelType>(mesh_, is, true)
240  );
241  }
242  };
243 
244 
245  // Member Functions
246 
247  // Access
248 
249  //- Return const access to force
250  inline const vector& f() const;
251 
252  //- Return const access to angular momentum
253  inline const vector& angularMomentum() const;
254 
255  //- Return const access to torque
256  inline const vector& torque() const;
257 
258  //- Return const access to the collision records
259  inline const collisionRecordList& collisionRecords() const;
260 
261  //- Return access to force
262  inline vector& f();
263 
264  //- Return access to angular momentum
265  inline vector& angularMomentum();
266 
267  //- Return access to torque
268  inline vector& torque();
269 
270  //- Return access to collision records
271  inline collisionRecordList& collisionRecords();
272 
273  //- Particle angular velocity
274  inline vector omega() const;
275 
276 
277  // Tracking
278 
279  //- Move the parcel
280  template<class TrackData>
281  bool move(TrackData& td, const scalar trackTime);
282 
283  //- Transform the physical properties of the particle
284  // according to the given transformation tensor
285  virtual void transformProperties(const tensor& T);
286 
287  //- Transform the physical properties of the particle
288  // according to the given separation vector
289  virtual void transformProperties(const vector& separation);
290 
291 
292  // I-O
293 
294  //- Read
295  template<class CloudType>
296  static void readFields(CloudType& c);
297 
298  //- Write
299  template<class CloudType>
300  static void writeFields(const CloudType& c);
301 
302 
303  // Ostream Operator
304 
305  friend Ostream& operator<< <ParcelType>
306  (
307  Ostream&,
309  );
310 };
311 
312 
313 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
314 
315 } // End namespace Foam
316 
317 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
318 
319 #include "CollidingParcelI.H"
320 
321 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
322 
323 #ifdef NoRepository
324  #include "CollidingParcel.C"
325 #endif
326 
327 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
328 
329 #endif
330 
331 // ************************************************************************* //
const vector & torque() const
Return const access to torque.
AddToPropertyList(ParcelType," (fx fy fz)"+" (angularMomentumx angularMomentumy angularMomentumz)"+" (torquex torquey torquez)"+" collisionRecordsPairAccessed"+" collisionRecordsPairOrigProcOfOther"+" collisionRecordsPairOrigIdOfOther"+" (collisionRecordsPairData)"+" collisionRecordsWallAccessed"+" collisionRecordsWallPRel"+" (collisionRecordsWallData)")
String representation of properties.
vector f_
Force on particle due to collisions [N].
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
CollisionRecordList< vector, vector > collisionRecordList
vectorFieldCompactIOField pairDataFieldCompactIOField
scalar poissonsRatio() const
Return const access to Poisson&#39;s ratio.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Wrapper around kinematic parcel types to add collision modelling.
Class to hold thermo particle constant properties.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
static void writeFields(const CloudType &c)
Write.
vectorFieldCompactIOField wallDataFieldCompactIOField
bool cp(const fileName &src, const fileName &dst)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:620
const vector & f() const
Return const access to force.
collisionRecordList collisionRecords_
Particle collision records.
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
const collisionRecordList & collisionRecords() const
Return const access to the collision records.
dynamicFvMesh & mesh
CollidingParcel(const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPtI)
Construct from owner, position, and cloud owner.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
vector omega() const
Particle angular velocity.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
vector angularMomentum_
Angular momentum of Parcel in global reference frame [kg m2/s].
vector torque_
Torque on particle due to collisions in global.
const vector & angularMomentum() const
Return const access to angular momentum.
scalar youngsModulus() const
Return const access to Young&#39;s Modulus.
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.
const dimensionedScalar c
Speed of light in a vacuum.
bool move(TrackData &td, const scalar trackTime)
Move the parcel.
TypeName("CollidingParcel")
Runtime type information.
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:53
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Namespace for OpenFOAM.
static void readFields(CloudType &c)
Read.