CollidingCloud.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-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::CollidingCloud
28
29Group
30 grpLagrangianIntermediateClouds
31
32Description
33 Adds coolisions to kinematic clouds
34
35SourceFiles
36 CollidingCloudI.H
37 CollidingCloud.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef CollidingCloud_H
42#define CollidingCloud_H
43
44#include "particle.H"
45#include "Cloud.H"
46#include "IOdictionary.H"
47#include "autoPtr.H"
48#include "fvMesh.H"
49#include "volFields.H"
50
51// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52
53namespace Foam
54{
55
56// Forward declaration of classes
57
58template<class CloudType>
59class CollisionModel;
60
61/*---------------------------------------------------------------------------*\
62 Class CollidingCloud Declaration
63\*---------------------------------------------------------------------------*/
64
65template<class CloudType>
67:
68 public CloudType
69{
70public:
71
72 // Public typedefs
73
74 //- Type of cloud this cloud was instantiated for
75 typedef CloudType cloudType;
76
77 //- Type of parcel the cloud was instantiated for
78 typedef typename CloudType::particleType parcelType;
79
80 //- Convenience typedef for this cloud type
82
83
84private:
85
86 // Private data
87
88 //- Cloud copy pointer
90
91
92 // Private Member Functions
93
94 //- No copy construct
95 CollidingCloud(const CollidingCloud&) = delete;
96
97 //- No copy assignment
98 void operator=(const CollidingCloud&) = delete;
99
100
101protected:
102
103 // Protected data
104
105 //- Thermo parcel constant properties
107
108
109 // References to the cloud sub-models
110
111 //- Collision model
114
115
116 // Initialisation
117
118 //- Set cloud sub-models
119 void setModels();
120
121
122 // Cloud evolution functions
123
124 //- Move-collide particles
125 template<class TrackCloudType>
126 void moveCollide
127 (
128 TrackCloudType& cloud,
129 typename parcelType::trackingData& td,
130 const scalar deltaT
131 );
132
133 //- Reset state of cloud
135
136
137public:
138
139 // Constructors
140
141 //- Construct given carrier gas fields
143 (
144 const word& cloudName,
145 const volScalarField& rho,
146 const volVectorField& U,
147 const volScalarField& mu,
148 const dimensionedVector& g,
149 bool readFields = true
150 );
151
152 //- Copy constructor with new name
154 (
156 const word& name
157 );
158
159 //- Copy constructor with new name - creates bare cloud
161 (
162 const fvMesh& mesh,
163 const word& name,
165 );
166
167 //- Construct and return clone based on (this) with new name
169 {
171 (
172 new CollidingCloud(*this, name)
173 );
174 }
175
176 //- Construct and return bare clone based on (this) with new name
177 virtual autoPtr<Cloud<parcelType>> cloneBare(const word& name) const
178 {
180 (
181 new CollidingCloud(this->mesh(), name, *this)
182 );
183 }
184
185
186 //- Destructor
187 virtual ~CollidingCloud();
188
189
190 // Member Functions
191
192 // Access
193
194 //- Return a reference to the cloud copy
195 inline const CollidingCloud& cloudCopy() const;
196
197 //- Return the constant properties
198 inline const typename parcelType::constantProperties&
199 constProps() const;
200
201
202 // Sub-models
203
204 //- Return const access to the collision model
206 collision() const;
207
208 //- Return reference to the collision model
210 collision();
211
212 // Check
213
214 //- Total rotational kinetic energy in the system
215 inline scalar rotationalKineticEnergyOfSystem() const;
216
217
218 // Cloud evolution functions
219
220 //- Store the current cloud state
221 void storeState();
222
223 //- Reset the current cloud to the previously stored state
224 void restoreState();
225
226 //- Evolve the cloud
227 void evolve();
228
229 //- Particle motion
230 template<class TrackCloudType>
231 void motion
232 (
233 TrackCloudType& cloud,
234 typename parcelType::trackingData& td
235 );
236
237
238 // I-O
239
240 //- Print cloud information
241 void info();
242};
243
244
245// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246
247} // End namespace Foam
248
249// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250
251#include "CollidingCloudI.H"
252
253// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254
255#ifdef NoRepository
256 #include "CollidingCloud.C"
257#endif
258
259// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260
261#endif
262
263// ************************************************************************* //
const uniformDimensionedVectorField & g
ParticleType particleType
Definition: Cloud.H:114
Adds coolisions to kinematic clouds.
const parcelType::constantProperties & constProps() const
Return the constant properties.
void setModels()
Set cloud sub-models.
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
void storeState()
Store the current cloud state.
autoPtr< CollisionModel< CollidingCloud< CloudType > > > collisionModel_
Collision model.
scalar rotationalKineticEnergyOfSystem() const
Total rotational kinetic energy in the system.
void cloudReset(CollidingCloud< CloudType > &c)
Reset state of cloud.
void moveCollide(TrackCloudType &cloud, typename parcelType::trackingData &td, const scalar deltaT)
Move-collide particles.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
const CollidingCloud & cloudCopy() const
Return a reference to the cloud copy.
CloudType cloudType
Type of cloud this cloud was instantiated for.
virtual autoPtr< Cloud< parcelType > > clone(const word &name)
Construct and return clone based on (this) with new name.
virtual autoPtr< Cloud< parcelType > > cloneBare(const word &name) const
Construct and return bare clone based on (this) with new name.
void evolve()
Evolve the cloud.
parcelType::constantProperties constProps_
Thermo parcel constant properties.
CollidingCloud< CloudType > collidingCloudType
Convenience typedef for this cloud type.
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
virtual ~CollidingCloud()
Destructor.
const CollisionModel< CollidingCloud< CloudType > > & collision() const
Return const access to the collision model.
Templated collision model class.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
const word & cloudName() const
Return the cloud type.
Definition: DSMCCloudI.H:37
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
Class to hold DSMC particle constant properties.
Definition: DSMCParcel.H:82
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:459
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 used to pass data into container.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
const volScalarField & mu
Namespace for OpenFOAM.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.