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-2020 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 classes
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  //- Which part of the integration algorithm is taking place
129  trackPart part_;
130 
131 
132  public:
133 
134  //- Constructors
135 
136  //- Construct from components
137  template<class TrackCloudType>
138  inline trackingData(const TrackCloudType& cloud);
139 
140 
141  //- Update the MPPIC averages
142  template<class TrackCloudType>
143  inline void updateAverages(const TrackCloudType& cloud);
144 
145 
146  //- Access
147 
148  //- Const access to the tracking part
149  inline trackPart part() const;
150 
151  //- Non const access to the tracking part
152  inline trackPart& part();
153  };
154 
155 
156 protected:
157 
158  // Protected data
159 
160  //- Velocity correction due to collisions [m/s]
162 
163 
164 public:
165 
166  // Static Data Members
167 
168  //- Runtime type information
169  TypeName("MPPICParcel");
170 
171  //- String representation of properties
173  (
174  ParcelType,
175  "(UCorrectx UCorrecty UCorrectz)"
176  );
177 
178 
179  // Constructors
180 
181  //- Construct from mesh, coordinates and topology
182  // Other properties initialised as null
183  inline MPPICParcel
184  (
185  const polyMesh& mesh,
186  const barycentric& coordinates,
187  const label celli,
188  const label tetFacei,
189  const label tetPti
190  );
191 
192  //- Construct from a position and a cell, searching for the rest of the
193  // required topology. Other properties are initialised as null.
194  inline MPPICParcel
195  (
196  const polyMesh& mesh,
197  const vector& position,
198  const label celli
199  );
200 
201  //- Construct from components
202  inline MPPICParcel
203  (
204  const polyMesh& mesh,
205  const barycentric& coordinates,
206  const label celli,
207  const label tetFacei,
208  const label tetPti,
209  const label typeId,
210  const scalar nParticle0,
211  const scalar d0,
212  const scalar dTarget0,
213  const vector& U0,
214  const vector& UCorrect0,
215  const typename ParcelType::constantProperties& constProps
216  );
217 
218  //- Construct from Istream
220  (
221  const polyMesh& mesh,
222  Istream& is,
223  bool readFields = true
224  );
225 
226  //- Construct as a copy
227  MPPICParcel(const MPPICParcel& p);
228 
229  //- Construct as a copy
230  MPPICParcel(const MPPICParcel& p, const polyMesh& mesh);
231 
232  //- Construct and return a (basic particle) clone
233  virtual autoPtr<particle> clone() const
234  {
235  return autoPtr<particle>(new MPPICParcel(*this));
236  }
237 
238  //- Construct and return a (basic particle) clone
239  virtual autoPtr<particle> clone(const polyMesh& mesh) const
240  {
241  return autoPtr<particle>(new MPPICParcel(*this, mesh));
242  }
243 
244  //- Factory class to read-construct particles used for
245  // parallel transfer
246  class iNew
247  {
248  const polyMesh& mesh_;
249 
250  public:
252  iNew(const polyMesh& mesh)
253  :
254  mesh_(mesh)
255  {}
257  autoPtr<MPPICParcel<ParcelType>> operator()(Istream& is) const
258  {
260  (
261  new MPPICParcel<ParcelType>(mesh_, is, true)
262  );
263  }
264  };
265 
266 
267  // Member Functions
268 
269  // Access
270 
271  //- Return const access to correction velocity
272  inline const vector& UCorrect() const;
273 
274  //- Return access to correction velocity
275  inline vector& UCorrect();
276 
277 
278  // Tracking
279 
280  //- Move the parcel
281  template<class TrackCloudType>
282  bool move
283  (
284  TrackCloudType& cloud,
285  trackingData& td,
286  const scalar trackTime
287  );
288 
289 
290  // Friend Functions
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 "MPPICParcelI.H"
321 
322 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
323 
324 #ifdef NoRepository
325  #include "MPPICParcel.C"
326 #endif
327 
328 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
329 
330 #endif
331 
332 // ************************************************************************* //
bool move(TrackCloudType &cloud, trackingData &td, const scalar trackTime)
Move the parcel.
Definition: MPPICParcel.C:58
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
Definition: MPPICParcel.H:232
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
const dimensionedScalar & c
Speed of light in a vacuum.
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:54
static void writeFields(const CloudType &c)
Write.
trackingData(const TrackCloudType &cloud)
Constructors.
PtrList< coordinateSystem > coordinates(solidRegions.size())
vector UCorrect_
Velocity correction due to collisions [m/s].
Definition: MPPICParcel.H:160
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
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:245
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.