MPPICCloud.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) 2013-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 "MPPICCloud.H"
27 #include "NoPacking.H"
28 #include "ParticleStressModel.H"
29 #include "NoDamping.H"
30 #include "NoIsotropy.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  const 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  setModels();
91 
92  if (readFields)
93  {
95  this->deleteLostParticles();
96  }
97 }
98 
99 
100 template<class CloudType>
102 (
103  const word& cloudName,
104  const volScalarField& rho,
105  const volVectorField& U,
106  const dimensionedVector& g,
107  const fluidThermo& carrierThermo,
108  const bool readFields
109 )
110 :
111  MPPICCloud(cloudName, rho, U, carrierThermo.mu(), g, readFields)
112 {}
113 
114 
115 template<class CloudType>
117 (
119  const word& name
120 )
121 :
122  CloudType(c, name),
123  packingModel_(c.packingModel_->clone()),
124  dampingModel_(c.dampingModel_->clone()),
125  isotropyModel_(c.isotropyModel_->clone())
126 {}
127 
128 
129 template<class CloudType>
131 (
132  const fvMesh& mesh,
133  const word& name,
134  const MPPICCloud<CloudType>& c
135 )
136 :
137  CloudType(mesh, name, c),
138  packingModel_(nullptr),
139  dampingModel_(nullptr),
140  isotropyModel_(nullptr)
141 {}
142 
143 
144 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
145 
146 template<class CloudType>
148 {}
149 
150 
151 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
152 
153 template<class CloudType>
155 {
156  cloudCopyPtr_.reset
157  (
158  static_cast<MPPICCloud<CloudType>*>
159  (
160  clone(this->name() + "Copy").ptr()
161  )
162  );
163 }
164 
165 
166 template<class CloudType>
168 {
169  this->cloudReset(cloudCopyPtr_());
170  cloudCopyPtr_.clear();
171 }
172 
173 
174 template<class CloudType>
176 {
177  if (this->solution().canEvolve())
178  {
179  typename parcelType::trackingData td(*this);
180 
181  this->solve(*this, td);
182  }
183 }
184 
185 
186 template<class CloudType>
187 template<class TrackCloudType>
189 (
190  TrackCloudType& cloud,
191  typename parcelType::trackingData& td
192 )
193 {
194  // Momentum
195  // ~~~~~~~~
196 
197  // force calculation and tracking
198  td.part() = parcelType::trackingData::tpPredictTrack;
199  CloudType::move(cloud, td, this->db().time().deltaTValue());
200 
201 
202  // Preliminary
203  // ~~~~~~~~~~~
204 
205  // switch forces off so they are not applied in corrector steps
206  this->forces().setCalcNonCoupled(false);
207  this->forces().setCalcCoupled(false);
208 
209 
210  // Damping
211  // ~~~~~~~
212 
213  if (!isType<DampingModels::NoDamping<CloudType>>(dampingModel_()))
214  {
215  if (this->mesh().moving())
216  {
218  << "MPPIC damping modelling does not support moving meshes."
219  << exit(FatalError);
220  }
221 
222  // update averages
223  td.updateAverages(cloud);
224 
225  // memory allocation and eulerian calculations
226  dampingModel_->cacheFields(true);
227 
228  // calculate the damping velocity corrections without moving the parcels
229  td.part() = parcelType::trackingData::tpDampingNoTrack;
230  CloudType::move(cloud, td, this->db().time().deltaTValue());
231 
232  // correct the parcel positions and velocities
233  td.part() = parcelType::trackingData::tpCorrectTrack;
234  CloudType::move(cloud, td, this->db().time().deltaTValue());
235 
236  // finalise and free memory
237  dampingModel_->cacheFields(false);
238  }
239 
240 
241  // Packing
242  // ~~~~~~~
243 
244  if (!isType<PackingModels::NoPacking<CloudType>>(packingModel_()))
245  {
246  if (this->mesh().moving())
247  {
249  << "MPPIC packing modelling does not support moving meshes."
250  << exit(FatalError);
251  }
252 
253  // same procedure as for damping
254  td.updateAverages(cloud);
255  packingModel_->cacheFields(true);
256  td.part() = parcelType::trackingData::tpPackingNoTrack;
257  CloudType::move(cloud, td, this->db().time().deltaTValue());
258  td.part() = parcelType::trackingData::tpCorrectTrack;
259  CloudType::move(cloud, td, this->db().time().deltaTValue());
260  packingModel_->cacheFields(false);
261  }
262 
263 
264  // Isotropy
265  // ~~~~~~~~
266 
267  if (!isType<IsotropyModels::NoIsotropy<CloudType>>(isotropyModel_()))
268  {
269  // update averages
270  td.updateAverages(cloud);
271 
272  // apply isotropy model
273  isotropyModel_->calculate();
274  }
275 
276 
277  // Final
278  // ~~~~~
279 
280  // update cell occupancy
281  this->updateCellOccupancy();
282 
283  // switch forces back on
284  this->forces().setCalcNonCoupled(true);
285  this->forces().setCalcCoupled(this->solution().coupled());
286 }
287 
288 
289 template<class CloudType>
291 {
292  CloudType::info();
293 
294  tmp<volScalarField> alpha = this->theta();
295 
296  const scalar alphaMin = gMin(alpha().primitiveField());
297  const scalar alphaMax = gMax(alpha().primitiveField());
298 
299  Info<< " Min cell volume fraction = " << alphaMin << endl;
300  Info<< " Max cell volume fraction = " << alphaMax << endl;
301 
302  if (alphaMax < small)
303  {
304  return;
305  }
306 
307  scalar nMin = great;
308 
309  forAll(this->mesh().cells(), celli)
310  {
311  const label n = this->cellOccupancy()[celli].size();
312 
313  if (n > 0)
314  {
315  const scalar nPack = n*alphaMax/alpha()[celli];
316 
317  if (nPack < nMin)
318  {
319  nMin = nPack;
320  }
321  }
322  }
323 
324  reduce(nMin, minOp<scalar>());
325 
326  Info<< " Min dense number of parcels = " << nMin << endl;
327 }
328 
329 
330 // ************************************************************************* //
virtual ~MPPICCloud()
Destructor.
Definition: MPPICCloud.C:147
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
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
Type gMin(const FieldField< Field, Type > &f)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void storeState()
Store the current cloud state.
Definition: MPPICCloud.C:154
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Base class for packing models.
Definition: MPPICCloud.H:53
autoPtr< DampingModel< MPPICCloud< CloudType > > > dampingModel_
Damping model.
Definition: MPPICCloud.H:113
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
Definition: MPPICCloud.C:189
fvMesh & mesh
void info()
I-O.
Definition: MPPICCloud.C:290
T clone(const T &t)
Definition: List.H:54
const cellShapeList & cells
autoPtr< IsotropyModel< MPPICCloud< CloudType > > > isotropyModel_
Exchange model.
Definition: MPPICCloud.H:117
A class for handling words, derived from string.
Definition: word.H:59
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:53
bool isType(const Type &t)
Check the typeid.
Definition: typeInfo.H:126
MPPICCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const volScalarField &mu, const dimensionedVector &g, const bool readFields=true)
Construct given carrier fields.
Definition: MPPICCloud.C:69
Type gMax(const FieldField< Field, Type > &f)
Adds MPPIC modelling to clouds.
Definition: MPPICCloud.H:74
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(viscosity->lookup("alphaMax"))
rhoEqn solve()
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:95
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:175
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedVector & g
Base class for collisional damping models.
Definition: MPPICCloud.H:56
autoPtr< PackingModel< MPPICCloud< CloudType > > > packingModel_
Packing model.
Definition: MPPICCloud.H:109
void restoreState()
Reset the current cloud to the previously stored state.
Definition: MPPICCloud.C:167