CollidingCloud.C
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 -------------------------------------------------------------------------------
10 License
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 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "CollidingCloud.H"
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  bool readFields
96 )
97 :
98  CloudType(cloudName, rho, U, mu, g, false),
99  constProps_(this->particleProperties()),
100  collisionModel_(nullptr)
101 {
102  if (this->solution().active())
103  {
104  setModels();
105 
106  if (readFields)
107  {
108  parcelType::readFields(*this);
109  this->deleteLostParticles();
110  }
111 
112  if
113  (
114  this->solution().steadyState()
115  && !isType<NoCollision<CollidingCloud<CloudType>>>(collisionModel_())
116  )
117  {
119  << "Collision modelling not currently available "
120  << "for steady state calculations" << exit(FatalError);
121  }
122  }
123 }
124 
125 
126 template<class CloudType>
128 (
130  const word& name
131 )
132 :
133  CloudType(c, name),
134  collisionModel_(c.collisionModel_->clone())
135 {}
136 
137 
138 template<class CloudType>
140 (
141  const fvMesh& mesh,
142  const word& name,
144 )
145 :
146  CloudType(mesh, name, c),
147  collisionModel_(nullptr)
148 {}
149 
150 
151 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
152 
153 template<class CloudType>
155 {}
156 
157 
158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
159 
160 template<class CloudType>
162 {
163  cloudCopyPtr_.reset
164  (
165  static_cast<CollidingCloud<CloudType>*>
166  (
167  clone(this->name() + "Copy").ptr()
168  )
169  );
170 }
171 
172 
173 template<class CloudType>
175 {
176  cloudReset(cloudCopyPtr_());
177  cloudCopyPtr_.clear();
178 }
179 
180 
181 template<class CloudType>
183 {
184  if (this->solution().canEvolve())
185  {
186  typename parcelType::trackingData td(*this);
187 
188  this->solve(*this, td);
189  }
190 }
191 
192 
193 template<class CloudType>
194 template<class TrackCloudType>
196 (
197  TrackCloudType& cloud,
198  typename parcelType::trackingData& td
199 )
200 {
201  // Sympletic leapfrog integration of particle forces:
202  // + apply half deltaV with stored force
203  // + move positions with new velocity
204  // + calculate forces in new position
205  // + apply half deltaV with new force
206 
207  label nSubCycles = collision().nSubCycles();
208 
209  if (nSubCycles > 1)
210  {
211  Info<< " " << nSubCycles << " move-collide subCycles" << endl;
212 
213  subCycleTime moveCollideSubCycle
214  (
215  const_cast<Time&>(this->db().time()),
216  nSubCycles
217  );
218 
219  while(!(++moveCollideSubCycle).end())
220  {
221  moveCollide(cloud, td, this->db().time().deltaTValue());
222  }
223 
224  moveCollideSubCycle.endSubCycle();
225  }
226  else
227  {
228  moveCollide(cloud, td, this->db().time().deltaTValue());
229  }
230 }
231 
232 
233 template<class CloudType>
235 {
236  CloudType::info();
237 
238  scalar rotationalKineticEnergy = rotationalKineticEnergyOfSystem();
239  reduce(rotationalKineticEnergy, sumOp<scalar>());
240 
241  Info<< " Rotational kinetic energy = "
242  << rotationalKineticEnergy << nl;
243 }
244 
245 
246 // ************************************************************************* //
Foam::CollidingCloud
Adds coolisions to kinematic clouds.
Definition: CollidingCloud.H:65
Foam::CollidingCloud::cloudReset
void cloudReset(CollidingCloud< CloudType > &c)
Reset state of cloud.
Definition: CollidingCloud.C:77
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::solution
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:51
cloudName
const word cloudName(propsDict.get< word >("cloud"))
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
CollidingCloud.H
Foam::CollidingCloud::evolve
void evolve()
Evolve the cloud.
Definition: CollidingCloud.C:182
Foam::CollidingCloud::info
void info()
Print cloud information.
Definition: CollidingCloud.C:234
Foam::CollidingCloud::moveCollide
void moveCollide(TrackCloudType &cloud, typename parcelType::trackingData &td, const scalar deltaT)
Move-collide particles.
Definition: CollidingCloud.C:51
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
CollisionModel.H
Foam::subCycleTime
A class for managing sub-cycling times.
Definition: subCycleTime.H:50
rho
rho
Definition: readInitialConditions.H:88
solve
CEqn solve()
Foam::sumOp
Definition: ops.H:213
Foam::CollidingCloud::restoreState
void restoreState()
Reset the current cloud to the previously stored state.
Definition: CollidingCloud.C:174
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::subCycleTime::endSubCycle
void endSubCycle()
End the sub-cycling and reset the time-state.
Definition: subCycleTime.C:68
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::CollidingCloud::motion
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
Definition: CollidingCloud.C:196
Foam::CloudType
DSMCCloud< dsmcParcel > CloudType
Definition: makeDSMCParcelBinaryCollisionModels.C:38
Foam::CollidingCloud::setModels
void setModels()
Set cloud sub-models.
Definition: CollidingCloud.C:35
Foam::FatalError
error FatalError
reduce
reduce(hasMovingMesh, orOp< bool >())
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned< vector >
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:24
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::readFields
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Definition: ReadFieldsTemplates.C:312
Foam::isType
bool isType(const Type &t)
Check is typeid is identical to the TargetType.
Definition: typeInfo.H:218
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::CollidingCloud::storeState
void storeState()
Store the current cloud state.
Definition: CollidingCloud.C:161
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:57
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
NoCollision.H
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::CollidingCloud::~CollidingCloud
virtual ~CollidingCloud()
Destructor.
Definition: CollidingCloud.C:154
Foam::CollidingParcel::trackingData
ParcelType::trackingData trackingData
Use base tracking data.
Definition: CollidingParcel.H:127
Foam::NoCollision
Place holder for 'none' option.
Definition: NoCollision.H:54
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::CollisionModel
Templated collision model class.
Definition: CollidingCloud.H:58