MPPICParcel.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) 2013-2019 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::MPPICParcel
26 
27 Description
28  Wrapper around kinematic parcel types to add MPPIC modelling
29 
30 SourceFiles
31  MPPICParcelI.H
32  MPPICParcelTrackingDataI.H
33  MPPICParcel.C
34  MPPICParcelIO.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef MPPICParcel_H
39 #define MPPICParcel_H
40 
41 #include "particle.H"
42 #include "labelFieldIOField.H"
43 #include "vectorFieldIOField.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 // Forward declaration of clases
51 
52 template<class ParcelType>
53 class MPPICParcel;
54 
55 template<class Type>
56 class AveragingMethod;
57 
58 // Forward declaration of friend functions
59 
60 template<class ParcelType>
61 Ostream& operator<<
62 (
63  Ostream&,
65 );
66 
67 /*---------------------------------------------------------------------------*\
68  Class MPPICParcel Declaration
69 \*---------------------------------------------------------------------------*/
70 
71 template<class ParcelType>
72 class MPPICParcel
73 :
74  public ParcelType
75 {
76  // Private Data
77 
78  //- Size in bytes of the fields
79  static const std::size_t sizeofFields_;
80 
81 
82 public:
83 
84  class trackingData
85  :
86  public ParcelType::trackingData
87  {
88 
89  public:
90 
91  enum trackPart
92  {
97  };
98 
99 
100  private:
101 
102  // Private Data
103 
104  // MPPIC Averages
105 
106  //- Volume average
107  autoPtr<AveragingMethod<scalar>> volumeAverage_;
108 
109  //- Radius average [ volume^(1/3) ]
110  autoPtr<AveragingMethod<scalar>> radiusAverage_;
111 
112  //- Density average
113  autoPtr<AveragingMethod<scalar>> rhoAverage_;
114 
115  //- Velocity average
117 
118  //- Magnitude velocity sqyuared average
119  autoPtr<AveragingMethod<scalar>> uSqrAverage_;
120 
121  //- Frequency average
122  autoPtr<AveragingMethod<scalar>> frequencyAverage_;
123 
124  //- Mass average
125  autoPtr<AveragingMethod<scalar>> massAverage_;
126 
127 
128  //- Label specifying the current part of the tracking process
129  trackPart part_;
130 
131 
132  public:
133 
134  //- Constructors
135 
136  //- Construct from components
137  template<class TrackCloudType>
138  inline trackingData
139  (
140  const TrackCloudType& cloud,
142  );
143 
144 
145  //- Update the MPPIC averages
146  template<class TrackCloudType>
147  inline void updateAverages(const TrackCloudType& cloud);
148 
149 
150  //- Access
151 
152  //- Const access to the tracking part label
153  inline trackPart part() const;
154 
155  //- Non const access to the tracking part label
156  inline trackPart& part();
157  };
158 
159 
160 protected:
161 
162  // Protected data
163 
164  //- Velocity correction due to collisions [m/s]
166 
167 
168 public:
169 
170  // Static Data Members
171 
172  //- Runtime type information
173  TypeName("MPPICParcel");
174 
175  //- String representation of properties
177  (
178  ParcelType,
179  "(UCorrectx UCorrecty UCorrectz)"
180  );
181 
182 
183  // Constructors
184 
185  //- Construct from mesh, coordinates and topology
186  // Other properties initialised as null
187  inline MPPICParcel
188  (
189  const polyMesh& mesh,
190  const barycentric& coordinates,
191  const label celli,
192  const label tetFacei,
193  const label tetPti
194  );
195 
196  //- Construct from a position and a cell, searching for the rest of the
197  // required topology. Other properties are initialised as null.
198  inline MPPICParcel
199  (
200  const polyMesh& mesh,
201  const vector& position,
202  const label celli
203  );
204 
205  //- Construct from components
206  inline MPPICParcel
207  (
208  const polyMesh& mesh,
209  const barycentric& coordinates,
210  const label celli,
211  const label tetFacei,
212  const label tetPti,
213  const label typeId,
214  const scalar nParticle0,
215  const scalar d0,
216  const scalar dTarget0,
217  const vector& U0,
218  const vector& UCorrect0,
219  const typename ParcelType::constantProperties& constProps
220  );
221 
222  //- Construct from Istream
224  (
225  const polyMesh& mesh,
226  Istream& is,
227  bool readFields = true
228  );
229 
230  //- Construct as a copy
231  MPPICParcel(const MPPICParcel& p);
232 
233  //- Construct as a copy
234  MPPICParcel(const MPPICParcel& p, const polyMesh& mesh);
235 
236  //- Construct and return a (basic particle) clone
237  virtual autoPtr<particle> clone() const
238  {
239  return autoPtr<particle>(new MPPICParcel(*this));
240  }
241 
242  //- Construct and return a (basic particle) clone
243  virtual autoPtr<particle> clone(const polyMesh& mesh) const
244  {
245  return autoPtr<particle>(new MPPICParcel(*this, mesh));
246  }
247 
248  //- Factory class to read-construct particles used for
249  // parallel transfer
250  class iNew
251  {
252  const polyMesh& mesh_;
253 
254  public:
256  iNew(const polyMesh& mesh)
257  :
258  mesh_(mesh)
259  {}
261  autoPtr<MPPICParcel<ParcelType>> operator()(Istream& is) const
262  {
264  (
265  new MPPICParcel<ParcelType>(mesh_, is, true)
266  );
267  }
268  };
269 
270 
271  // Member Functions
272 
273  // Access
274 
275  //- Return const access to correction velocity
276  inline const vector& UCorrect() const;
277 
278  //- Return access to correction velocity
279  inline vector& UCorrect();
280 
281 
282  // Tracking
283 
284  //- Move the parcel
285  template<class TrackCloudType>
286  bool move
287  (
288  TrackCloudType& cloud,
289  trackingData& td,
290  const scalar trackTime
291  );
292 
293 
294  // Friend Functions
295 
296  // I-O
297 
298  //- Read
299  template<class CloudType>
300  static void readFields(CloudType& c);
301 
302  //- Write
303  template<class CloudType>
304  static void writeFields(const CloudType& c);
305 
306 
307  // Ostream operator
308 
309  friend Ostream& operator<< <ParcelType>
310  (
311  Ostream&,
313  );
314 };
315 
316 
317 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
318 
319 } // End namespace Foam
320 
321 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
322 
323 #include "MPPICParcelI.H"
325 
326 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
327 
328 #ifdef NoRepository
329  #include "MPPICParcel.C"
330 #endif
331 
332 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
333 
334 #endif
335 
336 // ************************************************************************* //
bool move(TrackCloudType &cloud, trackingData &td, const scalar trackTime)
Move the parcel.
Definition: MPPICParcel.C:58
trackingData(const TrackCloudType &cloud, trackPart part=tpLinearTrack)
Constructors.
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
Definition: MPPICParcel.H:236
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
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
dynamicFvMesh & mesh
Wrapper around kinematic parcel types to add MPPIC modelling.
Definition: MPPICParcel.H:52
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
Base class for lagrangian averaging methods.
Definition: MPPICParcel.H:55
void updateAverages(const TrackCloudType &cloud)
Update the MPPIC averages.
const vector & UCorrect() const
Return const access to correction velocity.
Definition: MPPICParcelI.H:94
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static void writeFields(const CloudType &c)
Write.
PtrList< coordinateSystem > coordinates(solidRegions.size())
vector UCorrect_
Velocity correction due to collisions [m/s].
Definition: MPPICParcel.H:164
MPPICParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
Definition: MPPICParcelI.H:30
const dimensionedScalar c
Speed of light in a vacuum.
static void readFields(CloudType &c)
Read.
Definition: MPPICParcelIO.C:78
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:74
volScalarField & p
Factory class to read-construct particles used for.
Definition: MPPICParcel.H:249
TypeName("MPPICParcel")
Runtime type information.
AddToPropertyList(ParcelType, "(UCorrectx UCorrecty UCorrectz)")
String representation of properties.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
Namespace for OpenFOAM.