CollidingCloud.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) 2011-2018 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::CollidingCloud
26 
27 Description
28  Adds coolisions to kinematic clouds
29 
30 SourceFiles
31  CollidingCloudI.H
32  CollidingCloud.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef CollidingCloud_H
37 #define CollidingCloud_H
38 
39 #include "particle.H"
40 #include "Cloud.H"
41 #include "IOdictionary.H"
42 #include "autoPtr.H"
43 #include "fvMesh.H"
44 #include "volFields.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 // Forward declaration of classes
52 
53 template<class CloudType>
54 class CollisionModel;
55 
56 /*---------------------------------------------------------------------------*\
57  Class CollidingCloud Declaration
58 \*---------------------------------------------------------------------------*/
59 
60 template<class CloudType>
61 class CollidingCloud
62 :
63  public CloudType
64 {
65 public:
66 
67  // Public typedefs
68 
69  //- Type of cloud this cloud was instantiated for
70  typedef CloudType cloudType;
71 
72  //- Type of parcel the cloud was instantiated for
73  typedef typename CloudType::particleType parcelType;
74 
75  //- Convenience typedef for this cloud type
77 
78 
79 private:
80 
81  // Private data
82 
83  //- Cloud copy pointer
84  autoPtr<CollidingCloud<CloudType>> cloudCopyPtr_;
85 
86 
87  // Private Member Functions
88 
89  //- Disallow default bitwise copy construct
91 
92  //- Disallow default bitwise assignment
93  void operator=(const CollidingCloud&);
94 
95 
96 protected:
97 
98  // Protected data
99 
100  //- Thermo parcel constant properties
101  typename parcelType::constantProperties constProps_;
102 
103 
104  // References to the cloud sub-models
105 
106  //- Collision model
109 
110 
111  // Initialisation
112 
113  //- Set cloud sub-models
114  void setModels();
115 
116 
117  // Cloud evolution functions
118 
119  //- Move-collide particles
120  template<class TrackCloudType>
121  void moveCollide
122  (
123  TrackCloudType& cloud,
124  typename parcelType::trackingData& td,
125  const scalar deltaT
126  );
127 
128  //- Reset state of cloud
130 
131 
132 public:
133 
134  // Constructors
135 
136  //- Construct given carrier gas fields
138  (
139  const word& cloudName,
140  const volScalarField& rho,
141  const volVectorField& U,
142  const volScalarField& mu,
143  const dimensionedVector& g,
144  bool readFields = true
145  );
146 
147  //- Copy constructor with new name
149  (
151  const word& name
152  );
153 
154  //- Copy constructor with new name - creates bare cloud
156  (
157  const fvMesh& mesh,
158  const word& name,
160  );
161 
162  //- Construct and return clone based on (this) with new name
163  virtual autoPtr<Cloud<parcelType>> clone(const word& name)
164  {
166  (
167  new CollidingCloud(*this, name)
168  );
169  }
170 
171  //- Construct and return bare clone based on (this) with new name
172  virtual autoPtr<Cloud<parcelType>> cloneBare(const word& name) const
173  {
175  (
176  new CollidingCloud(this->mesh(), name, *this)
177  );
178  }
179 
180 
181  //- Destructor
182  virtual ~CollidingCloud();
183 
184 
185  // Member Functions
186 
187  // Access
188 
189  //- Return a reference to the cloud copy
190  inline const CollidingCloud& cloudCopy() const;
191 
192  //- Return the constant properties
193  inline const typename parcelType::constantProperties&
194  constProps() const;
195 
196 
197  // Sub-models
198 
199  //- Return const access to the collision model
201  collision() const;
202 
203  //- Return reference to the collision model
205  collision();
206 
207  // Check
208 
209  //- Total rotational kinetic energy in the system
210  inline scalar rotationalKineticEnergyOfSystem() const;
211 
212 
213  // Cloud evolution functions
214 
215  //- Store the current cloud state
216  void storeState();
217 
218  //- Reset the current cloud to the previously stored state
219  void restoreState();
220 
221  //- Evolve the cloud
222  void evolve();
223 
224  //- Particle motion
225  template<class TrackCloudType>
226  void motion
227  (
228  TrackCloudType& cloud,
229  typename parcelType::trackingData& td
230  );
231 
232 
233  // I-O
234 
235  //- Print cloud information
236  void info();
237 };
238 
239 
240 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242 } // End namespace Foam
243 
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 
246 #include "CollidingCloudI.H"
247 
248 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249 
250 #ifdef NoRepository
251  #include "CollidingCloud.C"
252 #endif
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 #endif
257 
258 // ************************************************************************* //
ParcelType particleType
Definition: Cloud.H:105
void cloudReset(CollidingCloud< CloudType > &c)
Reset state of cloud.
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:270
const word & name() const
Return name.
Definition: IOobject.H:297
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
void evolve()
Evolve the cloud.
const word & cloudName() const
Return the cloud type.
Definition: DSMCCloudI.H:34
const CollisionModel< CollidingCloud< CloudType > > & collision() const
Return const access to the collision model.
CloudType cloudType
Type of cloud this cloud was instantiated for.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
CollidingCloud< CloudType > collidingCloudType
Convenience typedef for this cloud type.
A class for handling words, derived from string.
Definition: word.H:59
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
void moveCollide(TrackCloudType &cloud, typename parcelType::trackingData &td, const scalar deltaT)
Move-collide particles.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
Templated collision model class.
Adds coolisions to kinematic clouds.
const parcelType::constantProperties & constProps() const
Return the constant properties.
const dimensionedScalar mu
Atomic mass unit.
void restoreState()
Reset the current cloud to the previously stored state.
const fvMesh & mesh() const
Return references to the mesh.
Definition: DSMCCloudI.H:41
void info()
Print cloud information.
U
Definition: pEqn.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionedScalar c
Speed of light in a vacuum.
parcelType::constantProperties constProps_
Thermo parcel constant properties.
virtual autoPtr< Cloud< parcelType > > cloneBare(const word &name) const
Construct and return bare clone based on (this) with new name.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
void setModels()
Set cloud sub-models.
virtual ~CollidingCloud()
Destructor.
const dimensionedVector & g
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
const CollidingCloud & cloudCopy() const
Return a reference to the cloud copy.
autoPtr< CollisionModel< CollidingCloud< CloudType > > > collisionModel_
Collision model.
void storeState()
Store the current cloud state.
scalar rotationalKineticEnergyOfSystem() const
Total rotational kinetic energy in the system.
Namespace for OpenFOAM.