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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2024 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::CollidingParcel
29 
30 Group
31  grpLagrangianIntermediateParcels
32 
33 Description
34  Wrapper around kinematic parcel types to add collision modelling
35 
36 SourceFiles
37  CollidingParcelI.H
38  CollidingParcel.C
39  CollidingParcelIO.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef CollidingParcel_H
44 #define CollidingParcel_H
45 
46 #include "particle.H"
47 #include "CollisionRecordList.H"
48 #include "labelFieldIOField.H"
49 #include "vectorFieldIOField.H"
50 #include "demandDrivenEntry.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 
60 
61 template<class ParcelType>
62 class CollidingParcel;
63 
64 // Forward declaration of friend functions
65 
66 template<class ParcelType>
67 Ostream& operator<<
68 (
69  Ostream&,
71 );
72 
73 /*---------------------------------------------------------------------------*\
74  Class CollidingParcel Declaration
75 \*---------------------------------------------------------------------------*/
76 
77 template<class ParcelType>
78 class CollidingParcel
79 :
80  public ParcelType
81 {
82 public:
83 
84  //- Size in bytes of the fields
85  static const std::size_t sizeofFields;
86 
87 
88  //- Class to hold thermo particle constant properties
89  class constantProperties
90  :
91  public ParcelType::constantProperties
92  {
93 
94  // Private data
95 
96  //- Young's modulus [N/m2]
97  demandDrivenEntry<scalar> youngsModulus_;
98 
99  //- Poisson's ratio
100  demandDrivenEntry<scalar> poissonsRatio_;
101 
102 
103  public:
104 
105  // Constructors
106 
107  //- Null constructor
109 
110  //- Copy constructor
112 
113  //- Construct from dictionary
114  constantProperties(const dictionary& parentDict);
115 
116 
117  // Member functions
118 
119  //- Return const access to Young's Modulus
120  inline scalar youngsModulus() const;
121 
122  //- Return const access to Poisson's ratio
123  inline scalar poissonsRatio() const;
124  };
125 
126 
127  //- Use base tracking data
128  typedef typename ParcelType::trackingData trackingData;
129 
130 
131 protected:
132 
133  // Protected data
134 
135  //- Force on particle due to collisions [N]
136  vector f_;
137 
138  //- Angular momentum of Parcel in global reference frame [kg m2/s]
140 
141  //- Torque on particle due to collisions in global
142  // reference frame [Nm]
144 
145  //- Particle collision records
147 
148 
149 public:
150 
151  // Static data members
152 
153  //- Runtime type information
154  TypeName("CollidingParcel");
155 
156  //- String representation of properties
158  (
159  ParcelType,
160  " (fx fy fz)"
161  + " (angularMomentumx angularMomentumy angularMomentumz)"
162  + " (torquex torquey torquez)"
163  + " collisionRecordsPairAccessed"
164  + " collisionRecordsPairOrigProcOfOther"
165  + " collisionRecordsPairOrigIdOfOther"
166  + " (collisionRecordsPairData)"
167  + " collisionRecordsWallAccessed"
168  + " collisionRecordsWallPRel"
169  + " (collisionRecordsWallData)"
170  );
171 
172 
173  // Constructors
174 
175  //- Construct from mesh, coordinates and topology
176  // Other properties initialised as null
177  inline CollidingParcel
178  (
179  const polyMesh& mesh,
180  const barycentric& coordinates,
181  const label celli,
182  const label tetFacei,
183  const label tetPti
184  );
185 
186  //- Construct from a position and a cell, searching for the rest of the
187  // required topology. Other properties are initialised as null.
188  inline CollidingParcel
189  (
190  const polyMesh& mesh,
191  const vector& position,
192  const label celli
193  );
194 
195  //- Construct from components
196  inline CollidingParcel
197  (
198  const polyMesh& mesh,
199  const barycentric& coordinates,
200  const label celli,
201  const label tetFacei,
202  const label tetPti,
203  const label typeId,
204  const scalar nParticle0,
205  const scalar d0,
206  const scalar dTarget0,
207  const vector& U0,
208  const vector& f0,
209  const vector& angularMomentum0,
210  const vector& torque0,
211  const typename ParcelType::constantProperties& constProps
212  );
213 
214  //- Construct from Istream
216  (
217  const polyMesh& mesh,
218  Istream& is,
219  bool readFields = true,
220  bool newFormat = true
221  );
222 
223  //- Construct as a copy
225 
226  //- Construct as a copy
228 
229  //- Return a (basic particle) clone
230  virtual autoPtr<particle> clone() const
231  {
232  return particle::Clone(*this);
233  }
234 
235  //- Return a (basic particle) clone
236  virtual autoPtr<particle> clone(const polyMesh& mesh) const
237  {
238  return particle::Clone(*this, mesh);
239  }
240 
241  //- Factory class to read-construct particles (for parallel transfer)
242  class iNew
243  {
244  const polyMesh& mesh_;
245 
246  public:
247 
248  iNew(const polyMesh& mesh)
249  :
250  mesh_(mesh)
251  {}
252 
253  autoPtr<CollidingParcel<ParcelType>> operator()(Istream& is) const
254  {
255  return autoPtr<CollidingParcel<ParcelType>>
256  (
257  new CollidingParcel<ParcelType>(mesh_, is, true)
258  );
259  }
260  };
261 
262 
263  // Member Functions
264 
265  // Access
266 
267  //- Return const access to force
268  inline const vector& f() const;
269 
270  //- Return const access to angular momentum
271  inline const vector& angularMomentum() const;
272 
273  //- Return const access to torque
274  inline const vector& torque() const;
275 
276  //- Return const access to the collision records
277  inline const collisionRecordList& collisionRecords() const;
278 
279  //- Return access to force
280  inline vector& f();
281 
282  //- Return access to angular momentum
283  inline vector& angularMomentum();
284 
285  //- Return access to torque
286  inline vector& torque();
287 
288  //- Return access to collision records
291  //- Particle angular velocity
292  inline vector omega() const;
293 
294 
295  // Tracking
297  //- Move the parcel
298  template<class TrackCloudType>
299  bool move
300  (
301  TrackCloudType& cloud,
302  trackingData& td,
303  const scalar trackTime
304  );
305 
306  //- Transform the physical properties of the particle
307  // according to the given transformation tensor
308  virtual void transformProperties(const tensor& T);
309 
310  //- Transform the physical properties of the particle
311  // according to the given separation vector
312  virtual void transformProperties(const vector& separation);
313 
314 
315  // I-O
316 
317  //- Read
318  template<class CloudType>
319  static void readFields(CloudType& c);
320 
321  //- Write
322  template<class CloudType>
323  static void writeFields(const CloudType& c);
324 
325  //- Write individual parcel properties to stream
326  void writeProperties
327  (
328  Ostream& os,
329  const wordRes& filters,
330  const word& delim,
331  const bool namesOnly
332  ) const;
333 
334  //- Read particle fields as objects from the obr registry
335  template<class CloudType>
336  static void readObjects(CloudType& c, const objectRegistry& obr);
337 
338  //- Write particle fields as objects into the obr registry
339  template<class CloudType>
340  static void writeObjects(const CloudType& c, objectRegistry& obr);
341 
342 
343  // Ostream Operator
344 
345  friend Ostream& operator<< <ParcelType>
346  (
347  Ostream&,
349  );
350 };
351 
352 
353 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
354 
355 } // End namespace Foam
356 
357 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
358 
359 #include "CollidingParcelI.H"
360 
361 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
362 
363 #ifdef NoRepository
364  #include "CollidingParcel.C"
365 #endif
366 
367 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
368 
369 #endif
370 
371 // ************************************************************************* //
ParcelType::trackingData trackingData
Use base tracking data.
vector f_
Force on particle due to collisions [N].
iNew(const polyMesh &mesh)
CollisionRecordList< vector, vector > collisionRecordList
vectorFieldCompactIOField pairDataFieldCompactIOField
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Wrapper around kinematic parcel types to add collision modelling.
Class to hold thermo particle constant properties.
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
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 the source to the destination (recursively if necessary).
Definition: POSIX.C:1063
static const std::size_t sizeofFields
Size in bytes of the fields.
static void writeFields(const CloudType &c)
Write.
vectorFieldCompactIOField wallDataFieldCompactIOField
scalar poissonsRatio() const
Return const access to Poisson&#39;s ratio.
scalar youngsModulus() const
Return const access to Young&#39;s Modulus.
collisionRecordList collisionRecords_
Particle collision records.
autoPtr< CollidingParcel< ParcelType > > operator()(Istream &is) const
dynamicFvMesh & mesh
AddToPropertyList(ParcelType, " (fx fy fz)"+" (angularMomentumx angularMomentumy angularMomentumz)"+" (torquex torquey torquez)"+" collisionRecordsPairAccessed"+" collisionRecordsPairOrigProcOfOther"+" collisionRecordsPairOrigIdOfOther"+" (collisionRecordsPairData)"+" collisionRecordsWallAccessed"+" collisionRecordsWallPRel"+" (collisionRecordsWallData)")
String representation of properties.
A class for handling words, derived from Foam::string.
Definition: word.H:63
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
Vector< scalar > vector
Definition: vector.H:57
bool move(TrackCloudType &cloud, trackingData &td, const scalar trackTime)
Move the parcel.
const vector & torque() const
Return const access to torque.
static autoPtr< particle > Clone(const Derived &p)
Clone a particle.
Definition: particle.H:552
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
const vector & angularMomentum() const
Return const access to angular momentum.
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.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
vector omega() const
Particle angular velocity.
PtrList< coordinateSystem > coordinates(solidRegions.size())
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.
CollidingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
TypeName("CollidingParcel")
Runtime type information.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
volScalarField & p
Registry of regIOobjects.
Tensor of scalars, i.e. Tensor<scalar>.
const collisionRecordList & collisionRecords() const
Return const access to the collision records.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
virtual autoPtr< particle > clone() const
Return a (basic particle) clone.
Namespace for OpenFOAM.
static void readFields(CloudType &c)
Read.