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-2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::CollidingParcel
29
30Group
31 grpLagrangianIntermediateParcels
32
33Description
34 Wrapper around kinematic parcel types to add collision modelling
35
36SourceFiles
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
54namespace Foam
55{
60
61template<class ParcelType>
62class CollidingParcel;
63
64// Forward declaration of friend functions
65
66template<class ParcelType>
67Ostream& operator<<
68(
69 Ostream&,
71);
72
73/*---------------------------------------------------------------------------*\
74 Class CollidingParcel Declaration
75\*---------------------------------------------------------------------------*/
76
77template<class ParcelType>
79:
80 public ParcelType
81{
82public:
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
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
131protected:
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
149public:
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,
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,
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 //- Construct and return a (basic particle) clone
230 virtual autoPtr<particle> clone() const
231 {
232 return autoPtr<particle>(new CollidingParcel(*this));
233 }
234
235 //- Construct and return a (basic particle) clone
236 virtual autoPtr<particle> clone(const polyMesh& mesh) const
237 {
238 return autoPtr<particle>(new CollidingParcel(*this, mesh));
239 }
240
241 //- Factory class to read-construct particles used for
242 // parallel transfer
243 class iNew
244 {
245 const polyMesh& mesh_;
246
247 public:
249 iNew(const polyMesh& mesh)
250 :
251 mesh_(mesh)
252 {}
255 {
257 (
258 new CollidingParcel<ParcelType>(mesh_, is, true)
259 );
260 }
261 };
262
263
264 // Member Functions
265
266 // Access
267
268 //- Return const access to force
269 inline const vector& f() const;
270
271 //- Return const access to angular momentum
272 inline const vector& angularMomentum() const;
273
274 //- Return const access to torque
275 inline const vector& torque() const;
276
277 //- Return const access to the collision records
278 inline const collisionRecordList& collisionRecords() const;
279
280 //- Return access to force
281 inline vector& f();
282
283 //- Return access to angular momentum
284 inline vector& angularMomentum();
285
286 //- Return access to torque
287 inline vector& torque();
288
289 //- Return access to collision records
291
292 //- Particle angular velocity
293 inline vector omega() const;
294
295
296 // Tracking
297
298 //- Move the parcel
299 template<class TrackCloudType>
300 bool move
301 (
302 TrackCloudType& cloud,
303 trackingData& td,
304 const scalar trackTime
305 );
306
307 //- Transform the physical properties of the particle
308 // according to the given transformation tensor
309 virtual void transformProperties(const tensor& T);
310
311 //- Transform the physical properties of the particle
312 // according to the given separation vector
313 virtual void transformProperties(const vector& separation);
314
315
316 // I-O
317
318 //- Read
319 template<class CloudType>
320 static void readFields(CloudType& c);
321
322 //- Write
323 template<class CloudType>
324 static void writeFields(const CloudType& c);
325
326 //- Write individual parcel properties to stream
327 void writeProperties
328 (
329 Ostream& os,
330 const wordRes& filters,
331 const word& delim,
332 const bool namesOnly
333 ) const;
334
335 //- Read particle fields as objects from the obr registry
336 template<class CloudType>
337 static void readObjects(CloudType& c, const objectRegistry& obr);
338
339 //- Write particle fields as objects into the obr registry
340 template<class CloudType>
341 static void writeObjects(const CloudType& c, objectRegistry& obr);
342
343
344 // Ostream Operator
346 friend Ostream& operator<< <ParcelType>
347 (
348 Ostream&,
350 );
351};
352
353
354// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
355
356} // End namespace Foam
357
358// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
359
360#include "CollidingParcelI.H"
361
362// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
363
364#ifdef NoRepository
365 #include "CollidingParcel.C"
366#endif
367
368// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
369
370#endif
371
372// ************************************************************************* //
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.
Factory class to read-construct particles used for.
iNew(const polyMesh &mesh)
autoPtr< CollidingParcel< ParcelType > > operator()(Istream &is) const
Wrapper around kinematic parcel types to add collision modelling.
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
virtual autoPtr< particle > clone(const polyMesh &mesh) const
Construct and return a (basic particle) clone.
CollidingParcel(const CollidingParcel &p)
Construct as a copy.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
const vector & angularMomentum() const
Return const access to angular momentum.
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
static const std::size_t sizeofFields
Size in bytes of the fields.
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.
ParcelType::trackingData trackingData
Use base tracking data.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
vector omega() const
Particle angular velocity.
const vector & f() const
Return const access to force.
collisionRecordList collisionRecords_
Particle collision records.
static void writeFields(const CloudType &c)
Write.
const vector & torque() const
Return const access to torque.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
CollidingParcel(const CollidingParcel &p, const polyMesh &mesh)
Construct as a copy.
bool move(TrackCloudType &cloud, trackingData &td, const scalar trackTime)
Move the parcel.
TypeName("CollidingParcel")
Runtime type information.
vector angularMomentum_
Angular momentum of Parcel in global reference frame [kg m2/s].
const collisionRecordList & collisionRecords() const
Return const access to the collision records.
static void readFields(CloudType &c)
Read.
vector torque_
Torque on particle due to collisions in global.
A Field of objects of type <T> with automated input and output using a compact storage....
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
Class for demand-driven dictionary entries.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Class used to pass data into container.
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
const volScalarField & T
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
vectorFieldCompactIOField pairDataFieldCompactIOField
CollisionRecordList< vector, vector > collisionRecordList
vectorFieldCompactIOField wallDataFieldCompactIOField
#define AddToPropertyList(ParcelType, str)
Add to existing static 'propertyList' for particle properties.
const volScalarField & cp
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73