CollidingCloud.C
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-2021 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 \*---------------------------------------------------------------------------*/
25 
26 #include "CollidingCloud.H"
27 #include "subCycleTime.H"
28 
29 #include "CollisionModel.H"
30 #include "NoCollision.H"
31 
32 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
33 
34 template<class CloudType>
36 {
37  collisionModel_.reset
38  (
40  (
41  this->subModelProperties(),
42  *this
43  ).ptr()
44  );
45 }
46 
47 
48 template<class CloudType>
49 template<class TrackCloudType>
51 (
52  TrackCloudType& cloud,
53  typename parcelType::trackingData& td,
54  const scalar deltaT
55 )
56 {
57  td.part() = parcelType::trackingData::tpVelocityHalfStep;
58  CloudType::move(cloud, td, deltaT);
59 
60  td.part() = parcelType::trackingData::tpLinearTrack;
61  CloudType::move(cloud, td, deltaT);
62 
63  // td.part() = parcelType::trackingData::tpRotationalTrack;
64  // CloudType::move(cloud, td, deltaT);
65 
66  this->updateCellOccupancy();
67 
68  this->collision().collide();
69 
70  td.part() = parcelType::trackingData::tpVelocityHalfStep;
71  CloudType::move(cloud, td, deltaT);
72 }
73 
74 
75 
76 template<class CloudType>
78 {
79  CloudType::cloudReset(c);
80 
81  collisionModel_.reset(c.collisionModel_.ptr());
82 }
83 
84 
85 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
86 
87 template<class CloudType>
89 (
90  const word& cloudName,
91  const volScalarField& rho,
92  const volVectorField& U,
93  const volScalarField& mu,
94  const dimensionedVector& g,
95  const bool readFields
96 )
97 :
98  CloudType(cloudName, rho, U, mu, g, false),
99  constProps_(this->particleProperties()),
100  collisionModel_(nullptr)
101 {
102  setModels();
103 
104  if (readFields)
105  {
106  parcelType::readFields(*this);
107  this->deleteLostParticles();
108  }
109 
110  if
111  (
112  this->solution().steadyState()
113  && !isType<NoCollision<CollidingCloud<CloudType>>>(collisionModel_())
114  )
115  {
117  << "Collision modelling not currently available "
118  << "for steady state calculations" << exit(FatalError);
119  }
120 }
121 
122 
123 template<class CloudType>
125 (
126  const word& cloudName,
127  const volScalarField& rho,
128  const volVectorField& U,
129  const dimensionedVector& g,
130  const fluidThermo& carrierThermo,
131  const bool readFields
132 )
133 :
134  CollidingCloud(cloudName, rho, U, carrierThermo.mu(), g, readFields)
135 {}
136 
137 
138 template<class CloudType>
140 (
142  const word& name
143 )
144 :
145  CloudType(c, name),
146  collisionModel_(c.collisionModel_->clone())
147 {}
148 
149 
150 template<class CloudType>
152 (
153  const fvMesh& mesh,
154  const word& name,
156 )
157 :
158  CloudType(mesh, name, c),
159  collisionModel_(nullptr)
160 {}
161 
162 
163 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
164 
165 template<class CloudType>
167 {}
168 
169 
170 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
171 
172 template<class CloudType>
174 {
175  cloudCopyPtr_.reset
176  (
177  static_cast<CollidingCloud<CloudType>*>
178  (
179  clone(this->name() + "Copy").ptr()
180  )
181  );
182 }
183 
184 
185 template<class CloudType>
187 {
188  cloudReset(cloudCopyPtr_());
189  cloudCopyPtr_.clear();
190 }
191 
192 
193 template<class CloudType>
195 {
196  if (this->solution().canEvolve())
197  {
198  typename parcelType::trackingData td(*this);
199 
200  this->solve(*this, td);
201  }
202 }
203 
204 
205 template<class CloudType>
206 template<class TrackCloudType>
208 (
209  TrackCloudType& cloud,
210  typename parcelType::trackingData& td
211 )
212 {
213  // Sympletic leapfrog integration of particle forces:
214  // + apply half deltaV with stored force
215  // + move positions with new velocity
216  // + calculate forces in new position
217  // + apply half deltaV with new force
218 
219  label nSubCycles = collision().nSubCycles();
220 
221  if (nSubCycles > 1)
222  {
223  Info<< " " << nSubCycles << " move-collide subCycles" << endl;
224 
225  subCycleTime moveCollideSubCycle
226  (
227  const_cast<Time&>(this->db().time()),
228  nSubCycles
229  );
230 
231  while(!(++moveCollideSubCycle).end())
232  {
233  moveCollide(cloud, td, this->db().time().deltaTValue());
234  }
235 
236  moveCollideSubCycle.endSubCycle();
237  }
238  else
239  {
240  moveCollide(cloud, td, this->db().time().deltaTValue());
241  }
242 }
243 
244 
245 template<class CloudType>
247 {
248  CloudType::info();
249 
250  scalar rotationalKineticEnergy = rotationalKineticEnergyOfSystem();
251  reduce(rotationalKineticEnergy, sumOp<scalar>());
252 
253  Info<< " Rotational kinetic energy = "
254  << rotationalKineticEnergy << nl;
255 }
256 
257 
258 // ************************************************************************* //
void cloudReset(CollidingCloud< CloudType > &c)
Reset state of cloud.
virtual tmp< volScalarField > mu() const =0
Dynamic viscosity of mixture [kg/m/s].
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
DSMCCloud< dsmcParcel > CloudType
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
CollidingCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const volScalarField &mu, const dimensionedVector &g, const bool readFields=true)
Construct given carrier fields.
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
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
void evolve()
Evolve the cloud.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
T clone(const T &t)
Definition: List.H:54
A class for handling words, derived from string.
Definition: word.H:59
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:53
void moveCollide(TrackCloudType &cloud, typename parcelType::trackingData &td, const scalar deltaT)
Move-collide particles.
bool isType(const Type &t)
Check the typeid.
Definition: typeInfo.H:126
Templated collision model class.
void endSubCycle()
End the sub-cycling and reset the time-state.
Definition: subCycleTime.C:59
Adds collisions to clouds.
static const char nl
Definition: Ostream.H:260
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Place holder for &#39;none&#39; option.
Definition: NoCollision.H:49
void restoreState()
Reset the current cloud to the previously stored state.
rhoEqn solve()
void info()
Print cloud information.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
messageStream Info
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
void setModels()
Set cloud sub-models.
virtual ~CollidingCloud()
Destructor.
const dimensionedVector & g
A class for managing sub-cycling times.
Definition: subCycleTime.H:48
autoPtr< CollisionModel< CollidingCloud< CloudType > > > collisionModel_
Collision model.
void storeState()
Store the current cloud state.