All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-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 \*---------------------------------------------------------------------------*/
25 
26 #include "CollidingCloud.H"
27 #include "CollisionModel.H"
28 #include "NoCollision.H"
29 
30 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
31 
32 template<class CloudType>
34 {
35  collisionModel_.reset
36  (
38  (
39  this->subModelProperties(),
40  *this
41  ).ptr()
42  );
43 }
44 
45 
46 template<class CloudType>
47 template<class TrackCloudType>
49 (
50  TrackCloudType& cloud,
51  typename parcelType::trackingData& td,
52  const scalar deltaT
53 )
54 {
55  td.part() = parcelType::trackingData::tpVelocityHalfStep;
56  CloudType::move(cloud, td, deltaT);
57 
58  td.part() = parcelType::trackingData::tpLinearTrack;
59  CloudType::move(cloud, td, deltaT);
60 
61  // td.part() = parcelType::trackingData::tpRotationalTrack;
62  // CloudType::move(cloud, td, deltaT);
63 
64  this->updateCellOccupancy();
65 
66  this->collision().collide();
67 
68  td.part() = parcelType::trackingData::tpVelocityHalfStep;
69  CloudType::move(cloud, td, deltaT);
70 }
71 
72 
73 
74 template<class CloudType>
76 {
77  CloudType::cloudReset(c);
78 
79  collisionModel_.reset(c.collisionModel_.ptr());
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84 
85 template<class CloudType>
87 (
88  const word& cloudName,
89  const volScalarField& rho,
90  const volVectorField& U,
91  const volScalarField& mu,
92  const dimensionedVector& g,
93  bool readFields
94 )
95 :
96  CloudType(cloudName, rho, U, mu, g, false),
97  constProps_(this->particleProperties()),
98  collisionModel_(nullptr)
99 {
100  if (this->solution().active())
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 
124 template<class CloudType>
126 (
128  const word& name
129 )
130 :
131  CloudType(c, name),
132  collisionModel_(c.collisionModel_->clone())
133 {}
134 
135 
136 template<class CloudType>
138 (
139  const fvMesh& mesh,
140  const word& name,
142 )
143 :
144  CloudType(mesh, name, c),
145  collisionModel_(nullptr)
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
150 
151 template<class CloudType>
153 {}
154 
155 
156 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
157 
158 template<class CloudType>
160 {
161  cloudCopyPtr_.reset
162  (
163  static_cast<CollidingCloud<CloudType>*>
164  (
165  clone(this->name() + "Copy").ptr()
166  )
167  );
168 }
169 
170 
171 template<class CloudType>
173 {
174  cloudReset(cloudCopyPtr_());
175  cloudCopyPtr_.clear();
176 }
177 
178 
179 template<class CloudType>
181 {
182  if (this->solution().canEvolve())
183  {
184  typename parcelType::trackingData td(*this);
185 
186  this->solve(*this, td);
187  }
188 }
189 
190 
191 template<class CloudType>
192 template<class TrackCloudType>
194 (
195  TrackCloudType& cloud,
196  typename parcelType::trackingData& td
197 )
198 {
199  // Sympletic leapfrog integration of particle forces:
200  // + apply half deltaV with stored force
201  // + move positions with new velocity
202  // + calculate forces in new position
203  // + apply half deltaV with new force
204 
205  label nSubCycles = collision().nSubCycles();
206 
207  if (nSubCycles > 1)
208  {
209  Info<< " " << nSubCycles << " move-collide subCycles" << endl;
210 
211  subCycleTime moveCollideSubCycle
212  (
213  const_cast<Time&>(this->db().time()),
214  nSubCycles
215  );
216 
217  while(!(++moveCollideSubCycle).end())
218  {
219  moveCollide(cloud, td, this->db().time().deltaTValue());
220  }
221 
222  moveCollideSubCycle.endSubCycle();
223  }
224  else
225  {
226  moveCollide(cloud, td, this->db().time().deltaTValue());
227  }
228 }
229 
230 
231 template<class CloudType>
233 {
234  CloudType::info();
235 
236  scalar rotationalKineticEnergy = rotationalKineticEnergyOfSystem();
237  reduce(rotationalKineticEnergy, sumOp<scalar>());
238 
239  Info<< " Rotational kinetic energy = "
240  << rotationalKineticEnergy << nl;
241 }
242 
243 
244 // ************************************************************************* //
void cloudReset(CollidingCloud< CloudType > &c)
Reset state of cloud.
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
DSMCCloud< dsmcParcel > CloudType
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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:319
void evolve()
Evolve the cloud.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
CollidingCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const volScalarField &mu, const dimensionedVector &g, bool readFields=true)
Construct given carrier gas fields.
rhoEqn solve()
T clone(const T &t)
Definition: List.H:54
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
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.
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 kinematic 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.
void info()
Print cloud information.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
messageStream Info
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
void setModels()
Set cloud sub-models.
virtual ~CollidingCloud()
Destructor.
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.