MPPICCloud.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 | Copyright (C) 2013-2016 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 "MPPICCloud.H"
27 #include "PackingModel.H"
28 #include "ParticleStressModel.H"
29 #include "DampingModel.H"
30 #include "IsotropyModel.H"
31 #include "TimeScaleModel.H"
32 
33 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
34 
35 template<class CloudType>
37 {
38  packingModel_.reset
39  (
41  (
42  this->subModelProperties(),
43  *this
44  ).ptr()
45  );
46  dampingModel_.reset
47  (
49  (
50  this->subModelProperties(),
51  *this
52  ).ptr()
53  );
54  isotropyModel_.reset
55  (
57  (
58  this->subModelProperties(),
59  *this
60  ).ptr()
61  );
62 }
63 
64 
65 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
66 
67 template<class CloudType>
69 (
70  const word& cloudName,
71  const volScalarField& rho,
72  const volVectorField& U,
73  const volScalarField& mu,
74  const dimensionedVector& g,
75  bool readFields
76 )
77 :
78  CloudType(cloudName, rho, U, mu, g, false),
79  packingModel_(nullptr),
80  dampingModel_(nullptr),
81  isotropyModel_(nullptr)
82 {
83  if (this->solution().steadyState())
84  {
86  << "MPPIC modelling not available for steady state calculations"
87  << exit(FatalError);
88  }
89 
90  if (this->solution().active())
91  {
92  setModels();
93 
94  if (readFields)
95  {
97  this->deleteLostParticles();
98  }
99  }
100 }
101 
102 
103 template<class CloudType>
105 (
107  const word& name
108 )
109 :
110  CloudType(c, name),
111  packingModel_(c.packingModel_->clone()),
112  dampingModel_(c.dampingModel_->clone()),
113  isotropyModel_(c.isotropyModel_->clone())
114 {}
115 
116 
117 template<class CloudType>
119 (
120  const fvMesh& mesh,
121  const word& name,
122  const MPPICCloud<CloudType>& c
123 )
124 :
125  CloudType(mesh, name, c),
126  packingModel_(nullptr),
127  dampingModel_(nullptr),
128  isotropyModel_(nullptr)
129 {}
130 
131 
132 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
133 
134 template<class CloudType>
136 {}
137 
138 
139 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140 
141 template<class CloudType>
143 {
144  cloudCopyPtr_.reset
145  (
146  static_cast<MPPICCloud<CloudType>*>
147  (
148  clone(this->name() + "Copy").ptr()
149  )
150  );
151 }
152 
153 
154 template<class CloudType>
156 {
157  this->cloudReset(cloudCopyPtr_());
158  cloudCopyPtr_.clear();
159 }
160 
161 
162 template<class CloudType>
164 {
165  if (this->solution().canEvolve())
166  {
167  typename parcelType::template
168  TrackingData<MPPICCloud<CloudType>> td(*this);
169 
170  this->solve(td);
171  }
172 }
173 
174 
175 template<class CloudType>
176 template<class TrackData>
178 {
179  // Kinematic
180  // ~~~~~~~~~
181 
182  // force calculation and tracking
183  td.part() = TrackData::tpLinearTrack;
184  CloudType::move(td, this->db().time().deltaTValue());
185 
186 
187  // Preliminary
188  // ~~~~~~~~~~~
189 
190  // switch forces off so they are not applied in corrector steps
191  this->forces().setCalcNonCoupled(false);
192  this->forces().setCalcCoupled(false);
193 
194 
195  // Damping
196  // ~~~~~~~
197 
198  if (dampingModel_->active())
199  {
200  // update averages
201  td.updateAverages(*this);
202 
203  // memory allocation and eulerian calculations
204  dampingModel_->cacheFields(true);
205 
206  // calculate the damping velocity corrections without moving the parcels
207  td.part() = TrackData::tpDampingNoTrack;
208  CloudType::move(td, this->db().time().deltaTValue());
209 
210  // correct the parcel positions and velocities
211  td.part() = TrackData::tpCorrectTrack;
212  CloudType::move(td, this->db().time().deltaTValue());
213 
214  // finalise and free memory
215  dampingModel_->cacheFields(false);
216  }
217 
218 
219  // Packing
220  // ~~~~~~~
221 
222  if (packingModel_->active())
223  {
224  // same procedure as for damping
225  td.updateAverages(*this);
226  packingModel_->cacheFields(true);
227  td.part() = TrackData::tpPackingNoTrack;
228  CloudType::move(td, this->db().time().deltaTValue());
229  td.part() = TrackData::tpCorrectTrack;
230  CloudType::move(td, this->db().time().deltaTValue());
231  packingModel_->cacheFields(false);
232  }
233 
234 
235  // Isotropy
236  // ~~~~~~~~
237 
238  if (isotropyModel_->active())
239  {
240  // update averages
241  td.updateAverages(*this);
242 
243  // apply isotropy model
244  isotropyModel_->calculate();
245  }
246 
247 
248  // Final
249  // ~~~~~
250 
251  // update cell occupancy
252  this->updateCellOccupancy();
253 
254  // switch forces back on
255  this->forces().setCalcNonCoupled(true);
256  this->forces().setCalcCoupled(this->solution().coupled());
257 }
258 
259 
260 template<class CloudType>
262 {
263  CloudType::info();
264 
265  tmp<volScalarField> alpha = this->theta();
266 
267  const scalar alphaMin = gMin(alpha().primitiveField());
268  const scalar alphaMax = gMax(alpha().primitiveField());
269 
270  Info<< " Min cell volume fraction = " << alphaMin << endl;
271  Info<< " Max cell volume fraction = " << alphaMax << endl;
272 
273  if (alphaMax < SMALL)
274  {
275  return;
276  }
277 
278  scalar nMin = GREAT;
279 
280  forAll(this->mesh().cells(), celli)
281  {
282  const label n = this->cellOccupancy()[celli].size();
283 
284  if (n > 0)
285  {
286  const scalar nPack = n*alphaMax/alpha()[celli];
287 
288  if (nPack < nMin)
289  {
290  nMin = nPack;
291  }
292  }
293  }
294 
295  reduce(nMin, minOp<scalar>());
296 
297  Info<< " Min dense number of parcels = " << nMin << endl;
298 }
299 
300 
301 // ************************************************************************* //
virtual ~MPPICCloud()
Destructor.
Definition: MPPICCloud.C:135
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
Type gMin(const FieldField< Field, Type > &f)
void storeState()
Store the current cloud state.
Definition: MPPICCloud.C:142
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Base class for packing models.
Definition: MPPICCloud.H:53
autoPtr< DampingModel< MPPICCloud< CloudType > > > dampingModel_
Damping model.
Definition: MPPICCloud.H:112
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
void info()
I-O.
Definition: MPPICCloud.C:261
dynamicFvMesh & mesh
const cellShapeList & cells
autoPtr< IsotropyModel< MPPICCloud< CloudType > > > isotropyModel_
Exchange model.
Definition: MPPICCloud.H:116
A class for handling words, derived from string.
Definition: word.H:59
rhoEqn solve()
Type gMax(const FieldField< Field, Type > &f)
Adds MPPIC modelling to kinematic clouds.
Definition: MPPICCloud.H:66
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
Base class for collisional return-to-isotropy models.
Definition: MPPICCloud.H:59
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
const List< DynamicList< molecule * > > & cellOccupancy
void setModels()
Set cloud sub-models.
Definition: MPPICCloud.C:36
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void motion(TrackData &td)
Particle motion.
Definition: MPPICCloud.C:177
messageStream Info
label n
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
void evolve()
Evolve the cloud.
Definition: MPPICCloud.C:163
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Base class for collisional damping models.
Definition: MPPICCloud.H:56
autoPtr< PackingModel< MPPICCloud< CloudType > > > packingModel_
Packing model.
Definition: MPPICCloud.H:108
void restoreState()
Reset the current cloud to the previously stored state.
Definition: MPPICCloud.C:155