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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-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 "MPPICCloud.H"
29 #include "PackingModel.H"
30 #include "ParticleStressModel.H"
31 #include "DampingModel.H"
32 #include "IsotropyModel.H"
33 #include "TimeScaleModel.H"
34 
35 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
36 
37 template<class CloudType>
39 {
40  packingModel_.reset
41  (
43  (
44  this->subModelProperties(),
45  *this
46  ).ptr()
47  );
48  dampingModel_.reset
49  (
51  (
52  this->subModelProperties(),
53  *this
54  ).ptr()
55  );
56  isotropyModel_.reset
57  (
59  (
60  this->subModelProperties(),
61  *this
62  ).ptr()
63  );
64 }
65 
66 
67 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68 
69 template<class CloudType>
71 (
72  const word& cloudName,
73  const volScalarField& rho,
74  const volVectorField& U,
75  const volScalarField& mu,
76  const dimensionedVector& g,
77  bool readFields
78 )
79 :
80  CloudType(cloudName, rho, U, mu, g, false),
81  packingModel_(nullptr),
82  dampingModel_(nullptr),
83  isotropyModel_(nullptr)
84 {
85  if (this->solution().active())
86  {
87  if (this->solution().steadyState())
88  {
90  << "MPPIC modelling not available for steady state calculations"
91  << exit(FatalError);
92  }
93 
94  setModels();
95 
96  if (readFields)
97  {
99  this->deleteLostParticles();
100  }
101  }
102 }
103 
104 
105 template<class CloudType>
107 (
108  MPPICCloud<CloudType>& c,
109  const word& name
110 )
111 :
112  CloudType(c, name),
113  packingModel_(c.packingModel_->clone()),
114  dampingModel_(c.dampingModel_->clone()),
115  isotropyModel_(c.isotropyModel_->clone())
116 {}
117 
118 
119 template<class CloudType>
121 (
122  const fvMesh& mesh,
123  const word& name,
124  const MPPICCloud<CloudType>& c
125 )
126 :
127  CloudType(mesh, name, c),
128  packingModel_(nullptr),
129  dampingModel_(nullptr),
130  isotropyModel_(nullptr)
131 {}
132 
133 
134 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
135 
136 template<class CloudType>
138 {}
139 
140 
141 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
142 
143 template<class CloudType>
145 {
146  cloudCopyPtr_.reset
147  (
148  static_cast<MPPICCloud<CloudType>*>
149  (
150  clone(this->name() + "Copy").ptr()
151  )
152  );
153 }
154 
155 
156 template<class CloudType>
158 {
159  this->cloudReset(cloudCopyPtr_());
160  cloudCopyPtr_.clear();
161 }
162 
163 
164 template<class CloudType>
166 {
167  if (this->solution().canEvolve())
168  {
169  typename parcelType::trackingData td(*this);
170 
171  this->solve(*this, td);
172  }
173 }
174 
175 
176 template<class CloudType>
177 template<class TrackCloudType>
179 (
180  TrackCloudType& cloud,
181  typename parcelType::trackingData& td
182 )
183 {
184  // Kinematic
185  // ~~~~~~~~~
186 
187  // force calculation and tracking
188  td.part() = parcelType::trackingData::tpLinearTrack;
189  CloudType::move(cloud, td, this->db().time().deltaTValue());
190 
191 
192  // Preliminary
193  // ~~~~~~~~~~~
194 
195  // switch forces off so they are not applied in corrector steps
196  this->forces().setCalcNonCoupled(false);
197  this->forces().setCalcCoupled(false);
198 
199 
200  // Damping
201  // ~~~~~~~
202 
203  if (dampingModel_->active())
204  {
205  if (this->mesh().moving())
206  {
208  << "MPPIC damping modelling does not support moving meshes."
209  << exit(FatalError);
210  }
211 
212  // update averages
213  td.updateAverages(cloud);
214 
215  // memory allocation and eulerian calculations
216  dampingModel_->cacheFields(true);
217 
218  // calculate the damping velocity corrections without moving the parcels
219  td.part() = parcelType::trackingData::tpDampingNoTrack;
220  CloudType::move(cloud, td, this->db().time().deltaTValue());
221 
222  // correct the parcel positions and velocities
223  td.part() = parcelType::trackingData::tpCorrectTrack;
224  CloudType::move(cloud, td, this->db().time().deltaTValue());
225 
226  // finalise and free memory
227  dampingModel_->cacheFields(false);
228  }
229 
230 
231  // Packing
232  // ~~~~~~~
233 
234  if (packingModel_->active())
235  {
236  if (this->mesh().moving())
237  {
239  << "MPPIC packing modelling does not support moving meshes."
240  << exit(FatalError);
241  }
242 
243  // same procedure as for damping
244  td.updateAverages(cloud);
245  packingModel_->cacheFields(true);
246  td.part() = parcelType::trackingData::tpPackingNoTrack;
247  CloudType::move(cloud, td, this->db().time().deltaTValue());
248  td.part() = parcelType::trackingData::tpCorrectTrack;
249  CloudType::move(cloud, td, this->db().time().deltaTValue());
250  packingModel_->cacheFields(false);
251  }
252 
253 
254  // Isotropy
255  // ~~~~~~~~
256 
257  if (isotropyModel_->active())
258  {
259  // update averages
260  td.updateAverages(cloud);
261 
262  // apply isotropy model
263  isotropyModel_->calculate();
264  }
265 
266 
267  // Final
268  // ~~~~~
269 
270  // update cell occupancy
271  this->updateCellOccupancy();
272 
273  // switch forces back on
274  this->forces().setCalcNonCoupled(true);
275  this->forces().setCalcCoupled(this->solution().coupled());
276 }
277 
278 
279 template<class CloudType>
281 {
282  CloudType::info();
283 
284  tmp<volScalarField> alpha = this->theta();
285 
286  const scalar alphaMin = gMin(alpha().primitiveField());
287  const scalar alphaMax = gMax(alpha().primitiveField());
288 
289  Log_ << " Min cell volume fraction = " << alphaMin << nl
290  << " Max cell volume fraction = " << alphaMax << endl;
291 
292  if (alphaMax < SMALL)
293  {
294  return;
295  }
296 
297  scalar nMin = GREAT;
298 
299  forAll(this->mesh().cells(), celli)
300  {
301  const label n = this->cellOccupancy()[celli].size();
302 
303  if (n > 0)
304  {
305  const scalar nPack = n*alphaMax/alpha()[celli];
306 
307  if (nPack < nMin)
308  {
309  nMin = nPack;
310  }
311  }
312  }
313 
314  reduce(nMin, minOp<scalar>());
315 
316  Log_<< " Min dense number of parcels = " << nMin << endl;
317 }
318 
319 
320 // ************************************************************************* //
virtual ~MPPICCloud()
Destructor.
Definition: MPPICCloud.C:130
DSMCCloud< dsmcParcel > CloudType
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Type gMin(const FieldField< Field, Type > &f)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
void storeState()
Store the current cloud state.
Definition: MPPICCloud.C:137
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Base class for packing models.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
Definition: MPPICCloud.C:172
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
CEqn solve()
void info()
I-O.
Definition: MPPICCloud.C:273
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
#define Log_
Report write to Foam::Info if the class log switch is true.
const cellShapeList & cells
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
const word cloudName(propsDict.get< word >("cloud"))
A class for handling words, derived from Foam::string.
Definition: word.H:63
const uniformDimensionedVectorField & g
Type gMax(const FieldField< Field, Type > &f)
Adds MPPIC modelling to kinematic clouds.
Definition: MPPICCloud.H:67
const dimensionedScalar mu
Atomic mass unit.
Base class for collisional return-to-isotropy models.
const List< DynamicList< molecule * > > & cellOccupancy
U
Definition: pEqn.H:72
void setModels()
Set cloud sub-models.
Definition: MPPICCloud.C:31
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionedScalar c
Speed of light in a vacuum.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
label n
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:92
void evolve()
Evolve the cloud.
Definition: MPPICCloud.C:158
A class for managing temporary objects.
Definition: HashPtrTable.H:50
bool coupled
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
Base class for collisional damping models.
void restoreState()
Reset the current cloud to the previously stored state.
Definition: MPPICCloud.C:150