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 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 MPPICParcelName Declaration
69 \*---------------------------------------------------------------------------*/
70 
72 
73 
74 /*---------------------------------------------------------------------------*\
75  Class MPPICParcel Declaration
76 \*---------------------------------------------------------------------------*/
77 
78 template<class ParcelType>
79 class MPPICParcel
80 :
81  public ParcelType,
82  public MPPICParcelName
83 {
84  // Private Data
85 
86  //- Size in bytes of the fields
87  static const std::size_t sizeofFields_;
88 
89 
90 public:
91 
92  class trackingData
93  :
94  public ParcelType::trackingData
95  {
96 
97  public:
98 
99  enum trackPart
100  {
105  };
106 
107 
108  private:
109 
110  // Private Data
111 
112  // MPPIC Averages
113 
114  //- Volume average
115  autoPtr<AveragingMethod<scalar>> volumeAverage_;
116 
117  //- Radius average [ volume^(1/3) ]
118  autoPtr<AveragingMethod<scalar>> radiusAverage_;
119 
120  //- Density average
121  autoPtr<AveragingMethod<scalar>> rhoAverage_;
122 
123  //- Velocity average
125 
126  //- Magnitude velocity sqyuared average
127  autoPtr<AveragingMethod<scalar>> uSqrAverage_;
128 
129  //- Frequency average
130  autoPtr<AveragingMethod<scalar>> frequencyAverage_;
131 
132  //- Mass average
133  autoPtr<AveragingMethod<scalar>> massAverage_;
134 
135 
136  //- Which part of the integration algorithm is taking place
137  trackPart part_;
138 
139 
140  public:
141 
142  //- Constructors
143 
144  //- Construct from components
145  template<class TrackCloudType>
146  inline trackingData(const TrackCloudType& cloud);
147 
148 
149  //- Update the MPPIC averages
150  template<class TrackCloudType>
151  inline void updateAverages(const TrackCloudType& cloud);
152 
153 
154  //- Access
155 
156  //- Const access to the tracking part
157  inline trackPart part() const;
158 
159  //- Non const access to the tracking part
160  inline trackPart& part();
161  };
162 
163 
164 protected:
165 
166  // Protected data
167 
168  //- Velocity correction due to collisions [m/s]
170 
171 
172 public:
173 
174  // Static Data Members
175 
176  //- String representation of properties
178  (
179  ParcelType,
180  "(UCorrectx UCorrecty UCorrectz)"
181  );
182 
183 
184  // Constructors
185 
186  //- Construct from mesh, coordinates and topology
187  // Other properties initialised as null
188  inline MPPICParcel
189  (
190  const polyMesh& mesh,
191  const barycentric& coordinates,
192  const label celli,
193  const label tetFacei,
194  const label tetPti
195  );
196 
197  //- Construct from a position and a cell, searching for the rest of the
198  // required topology. Other properties are initialised as null.
199  inline MPPICParcel
200  (
201  const polyMesh& mesh,
202  const vector& position,
203  const label celli
204  );
205 
206  //- Construct from Istream
208  (
209  const polyMesh& mesh,
210  Istream& is,
211  bool readFields = true
212  );
213 
214  //- Construct as a copy
215  MPPICParcel(const MPPICParcel& p);
216 
217  //- Construct as a copy
218  MPPICParcel(const MPPICParcel& p, const polyMesh& mesh);
219 
220  //- Construct and return a (basic particle) clone
221  virtual autoPtr<particle> clone() const
222  {
223  return autoPtr<particle>(new MPPICParcel(*this));
224  }
225 
226  //- Construct and return a (basic particle) clone
227  virtual autoPtr<particle> clone(const polyMesh& mesh) const
228  {
229  return autoPtr<particle>(new MPPICParcel(*this, mesh));
230  }
231 
232  //- Factory class to read-construct particles used for
233  // parallel transfer
234  class iNew
235  {
236  const polyMesh& mesh_;
237 
238  public:
240  iNew(const polyMesh& mesh)
241  :
242  mesh_(mesh)
243  {}
245  autoPtr<MPPICParcel<ParcelType>> operator()(Istream& is) const
246  {
248  (
249  new MPPICParcel<ParcelType>(mesh_, is, true)
250  );
251  }
252  };
253 
254 
255  // Member Functions
256 
257  // Access
258 
259  //- Return const access to correction velocity
260  inline const vector& UCorrect() const;
261 
262  //- Return access to correction velocity
263  inline vector& UCorrect();
264 
265 
266  // Tracking
267 
268  //- Move the parcel
269  template<class TrackCloudType>
270  bool move
271  (
272  TrackCloudType& cloud,
273  trackingData& td,
274  const scalar trackTime
275  );
276 
277 
278  // Friend Functions
279 
280  // I-O
281 
282  //- Read
283  template<class CloudType>
284  static void readFields(CloudType& c);
285 
286  //- Write
287  template<class CloudType>
288  static void writeFields(const CloudType& c);
289 
290 
291  // Ostream operator
292 
293  friend Ostream& operator<< <ParcelType>
294  (
295  Ostream&,
297  );
298 };
299 
300 
301 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
302 
303 } // End namespace Foam
304 
305 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306 
307 #include "MPPICParcelI.H"
309 
310 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
311 
312 #ifdef NoRepository
313  #include "MPPICParcel.C"
314 #endif
315 
316 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
317 
318 #endif
319 
320 // ************************************************************************* //
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:220
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
fvMesh & mesh
const dimensionedScalar c
Speed of light in a vacuum.
Wrapper around parcel types to add MPPIC modelling.
Definition: MPPICParcel.H:52
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
TemplateName(FvFaceCellWave)
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:59
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.
vector UCorrect_
Velocity correction due to collisions [m/s].
Definition: MPPICParcel.H:168
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:76
volScalarField & p
Factory class to read-construct particles used for.
Definition: MPPICParcel.H:233
AddToPropertyList(ParcelType, "(UCorrectx UCorrecty UCorrectz)")
String representation of properties.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
Namespace for OpenFOAM.