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-2022 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 template<class CloudType>
67 {
68  CloudType::cloudReset(c);
69 
70  packingModel_.reset(c.packingModel_.ptr());
71  dampingModel_.reset(c.dampingModel_.ptr());
72  isotropyModel_.reset(c.isotropyModel_.ptr());
73 }
74 
75 
76 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
77 
78 template<class CloudType>
80 (
81  const word& cloudName,
82  const volScalarField& rho,
83  const volVectorField& U,
84  const volScalarField& mu,
85  const dimensionedVector& g,
86  const bool readFields
87 )
88 :
89  CloudType(cloudName, rho, U, mu, g, false),
90  packingModel_(nullptr),
91  dampingModel_(nullptr),
92  isotropyModel_(nullptr)
93 {
94  if (this->solution().steadyState())
95  {
97  << "MPPIC modelling not available for steady state calculations"
98  << exit(FatalError);
99  }
100 
101  setModels();
102 
103  if (readFields)
104  {
105  parcelType::readFields(*this);
106  this->deleteLostParticles();
107  }
108 }
109 
110 
111 template<class CloudType>
113 (
114  const word& cloudName,
115  const volScalarField& rho,
116  const volVectorField& U,
117  const dimensionedVector& g,
118  const fluidThermo& carrierThermo,
119  const bool readFields
120 )
121 :
122  MPPICCloud(cloudName, rho, U, carrierThermo.mu(), g, readFields)
123 {}
124 
125 
126 template<class CloudType>
128 (
130  const word& name
131 )
132 :
133  CloudType(c, name),
134  packingModel_(c.packingModel_->clone()),
135  dampingModel_(c.dampingModel_->clone()),
136  isotropyModel_(c.isotropyModel_->clone())
137 {}
138 
139 
140 template<class CloudType>
142 (
143  const fvMesh& mesh,
144  const word& name,
145  const MPPICCloud<CloudType>& c
146 )
147 :
148  CloudType(mesh, name, c),
149  packingModel_(nullptr),
150  dampingModel_(nullptr),
151  isotropyModel_(nullptr)
152 {}
153 
154 
155 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
156 
157 template<class CloudType>
159 {}
160 
161 
162 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
163 
164 template<class CloudType>
166 {
167  cloudCopyPtr_.reset
168  (
169  static_cast<MPPICCloud<CloudType>*>
170  (
171  clone(this->name() + "Copy").ptr()
172  )
173  );
174 }
175 
176 
177 template<class CloudType>
179 {
180  cloudReset(cloudCopyPtr_());
181 
182  cloudCopyPtr_.clear();
183 }
184 
185 
186 template<class CloudType>
188 {
189  if (this->solution().canEvolve())
190  {
191  typename parcelType::trackingData td(*this);
192 
193  this->solve(*this, td);
194  }
195 }
196 
197 
198 template<class CloudType>
199 template<class TrackCloudType>
201 (
202  TrackCloudType& cloud,
203  typename parcelType::trackingData& td
204 )
205 {
206  // Assign parcel ID-s
207  label i = 0;
208  forAllIter(typename MPPICCloud<CloudType>, *this, iter)
209  {
210  iter().id() = labelPair(Pstream::myProcNo(), i ++);
211  }
212 
213  // Create a copy of all parcels and sources to use as a predictor
214  autoPtr<MPPICCloud<CloudType>> predictorCloudPtr
215  (
216  static_cast<MPPICCloud<CloudType>*>
217  (
218  clone(this->name() + "Predictor").ptr()
219  )
220  );
221  MPPICCloud<CloudType>& predictorCloud = predictorCloudPtr();
222 
223  // Predictor move
224  predictorCloud.CloudType::move(predictorCloud, td);
225 
226  // Calculate correction velocities
227  const scalar trackTime = td.trackTime();
228  td.updateAverages(predictorCloud);
229  predictorCloud.dampingModel().cacheFields(true);
230  predictorCloud.packingModel().cacheFields(true);
231  vectorField UCorr(this->size(), Zero);
234  forAllIter(typename MPPICCloud<CloudType>, predictorCloud, iter)
235  {
236  const labelPair& id = iter().id();
237 
238  const vector dU =
239  predictorCloud.packingModel().velocityCorrection(iter(), trackTime)
240  + predictorCloud.dampingModel().velocityCorrection(iter(), trackTime);
241 
242  if (id.first() == Pstream::myProcNo())
243  {
244  UCorr[id.second()] = dU;
245  }
246  else
247  {
248  UCorrProc[id.first()].append(dU);
249  IDProc[id.first()].append(id.second());
250  }
251  }
252  predictorCloud.dampingModel().cacheFields(false);
253  predictorCloud.packingModel().cacheFields(false);
254 
255  // Distribute the correction velocities
257  if (Pstream::parRun())
258  {
259  forAll(UCorrProc, proci)
260  {
261  if (proci == Pstream::myProcNo()) continue;
262 
263  UOPstream os(proci, pBufs);
264 
265  os << UCorrProc[proci] << IDProc[proci];
266  }
267 
268  pBufs.finishedSends();
269 
270  forAll(UCorrProc, proci)
271  {
272  if (proci == Pstream::myProcNo()) continue;
273 
274  UIPstream is(proci, pBufs);
275 
276  is >> UCorrProc[proci] >> IDProc[proci];
277  }
278  }
279  forAll(UCorrProc, proci)
280  {
281  if (proci == Pstream::myProcNo()) continue;
282 
283  forAll(UCorrProc[proci], i)
284  {
285  UCorr[IDProc[proci][i]] = UCorrProc[proci][i];
286  }
287  }
288 
289  // Apply the correction velocities to the parcels
290  forAllIter(typename MPPICCloud<CloudType>, *this, iter)
291  {
292  iter().U() += UCorr[iter().id().second()];
293  }
294 
295  // Corrector
296  CloudType::move(cloud, td);
297 
298  // Apply isotropy model
299  td.updateAverages(cloud);
300  isotropyModel_->calculate();
301 
302  // Update cell occupancy
303  this->updateCellOccupancy();
304 }
305 
306 
307 template<class CloudType>
309 {
310  CloudType::info();
311 
312  tmp<volScalarField> alpha = this->theta();
313 
314  const scalar alphaMin = gMin(alpha().primitiveField());
315  const scalar alphaMax = gMax(alpha().primitiveField());
316 
317  Info<< " Min cell volume fraction = " << alphaMin << endl;
318  Info<< " Max cell volume fraction = " << alphaMax << endl;
319 
320  if (alphaMax < small)
321  {
322  return;
323  }
324 
325  scalar nMin = great;
326 
327  forAll(this->mesh().cells(), celli)
328  {
329  const label n = this->cellOccupancy()[celli].size();
330 
331  if (n > 0)
332  {
333  const scalar nPack = n*alphaMax/alpha()[celli];
334 
335  if (nPack < nMin)
336  {
337  nMin = nPack;
338  }
339  }
340  }
341 
342  reduce(nMin, minOp<scalar>());
343 
344  Info<< " Min dense number of parcels = " << nMin << endl;
345 }
346 
347 
348 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
const List< DynamicList< molecule * > > & cellOccupancy
void deleteLostParticles()
Remove lost particles from cloud and delete.
Definition: Cloud.C:232
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td)
Move the particles.
Definition: Cloud.C:275
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:79
void info() const
Print cloud information.
Definition: DSMCCloud.C:970
Base class for collisional damping models.
Definition: DampingModel.H:60
Generic GeometricField class.
Base class for collisional return-to-isotropy models.
Definition: IsotropyModel.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
Adds MPPIC modelling to clouds.
Definition: MPPICCloud.H:78
const DampingModel< MPPICCloud< CloudType > > & dampingModel() const
Return condt access to the damping model.
Definition: MPPICCloudI.H:54
void setModels()
Set cloud sub-models.
Definition: MPPICCloud.C:36
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
Definition: MPPICCloud.C:201
void storeState()
Store the current cloud state.
Definition: MPPICCloud.C:165
virtual ~MPPICCloud()
Destructor.
Definition: MPPICCloud.C:158
void evolve()
Evolve the cloud.
Definition: MPPICCloud.C:187
void info()
I-O.
Definition: MPPICCloud.C:308
void restoreState()
Reset the current cloud to the previously stored state.
Definition: MPPICCloud.C:178
void cloudReset(MPPICCloud< CloudType > &c)
Reset state of cloud.
Definition: MPPICCloud.C:66
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:80
const PackingModel< MPPICCloud< CloudType > > & packingModel() const
Return const access to the packing model.
Definition: MPPICCloudI.H:38
Base class for packing models.
Definition: PackingModel.H:65
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void finishedSends(const bool block=true)
Mark all sends as having been done. This will start receives.
Input inter-processor communications stream operating on external buffer.
Definition: UIPstream.H:57
Output inter-processor communications stream operating on external buffer.
Definition: UOPstream.H:58
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:51
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const cellShapeList & cells
U
Definition: pEqn.H:72
dimensionedScalar alphaMax(viscosity->lookup("alphaMax"))
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar mu
Atomic mass unit.
const dimensionedScalar c
Speed of light in a vacuum.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
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
messageStream Info
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T clone(const T &t)
Definition: List.H:55
error FatalError
Type gMin(const FieldField< Field, Type > &f)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Type gMax(const FieldField< Field, Type > &f)
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
const word cloudName(propsDict.lookup("cloudName"))