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-2014 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_(NULL),
80  dampingModel_(NULL),
81  isotropyModel_(NULL)
82 {
83  if (this->solution().steadyState())
84  {
86  (
87  "Foam::MPPICCloud<CloudType>::MPPICCloud"
88  "("
89  "const word&, "
90  "const volScalarField&, "
91  "const volVectorField&, "
92  "const volScalarField&, "
93  "const dimensionedVector&, "
94  "bool"
95  ")"
96  ) << "MPPIC modelling not available for steady state calculations"
97  << exit(FatalError);
98  }
99 
100  if (this->solution().active())
101  {
102  setModels();
103 
104  if (readFields)
105  {
106  parcelType::readFields(*this);
107  }
108  }
109 }
110 
111 
112 template<class CloudType>
114 (
116  const word& name
117 )
118 :
119  CloudType(c, name),
120  packingModel_(c.packingModel_->clone()),
121  dampingModel_(c.dampingModel_->clone()),
122  isotropyModel_(c.isotropyModel_->clone())
123 {}
124 
125 
126 template<class CloudType>
128 (
129  const fvMesh& mesh,
130  const word& name,
131  const MPPICCloud<CloudType>& c
132 )
133 :
134  CloudType(mesh, name, c),
135  packingModel_(NULL),
136  dampingModel_(NULL),
137  isotropyModel_(NULL)
138 {}
139 
140 
141 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
142 
143 template<class CloudType>
145 {}
146 
147 
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
149 
150 template<class CloudType>
152 {
153  cloudCopyPtr_.reset
154  (
155  static_cast<MPPICCloud<CloudType>*>
156  (
157  clone(this->name() + "Copy").ptr()
158  )
159  );
160 }
161 
162 
163 template<class CloudType>
165 {
166  this->cloudReset(cloudCopyPtr_());
167  cloudCopyPtr_.clear();
168 }
169 
170 
171 template<class CloudType>
173 {
174  if (this->solution().canEvolve())
175  {
176  typename parcelType::template
177  TrackingData<MPPICCloud<CloudType> > td(*this);
178 
179  this->solve(td);
180  }
181 }
182 
183 
184 template<class CloudType>
185 template<class TrackData>
187 {
188  // Kinematic
189  // ~~~~~~~~~
190 
191  // force calculation and tracking
192  td.part() = TrackData::tpLinearTrack;
193  CloudType::move(td, this->db().time().deltaTValue());
194 
195 
196  // Preliminary
197  // ~~~~~~~~~~~
198 
199  // switch forces off so they are not applied in corrector steps
200  this->forces().setCalcNonCoupled(false);
201  this->forces().setCalcCoupled(false);
202 
203 
204  // Damping
205  // ~~~~~~~
206 
207  if (dampingModel_->active())
208  {
209  // update averages
210  td.updateAverages(*this);
211 
212  // memory allocation and eulerian calculations
213  dampingModel_->cacheFields(true);
214 
215  // calculate the damping velocity corrections without moving the parcels
216  td.part() = TrackData::tpDampingNoTrack;
217  CloudType::move(td, this->db().time().deltaTValue());
218 
219  // correct the parcel positions and velocities
220  td.part() = TrackData::tpCorrectTrack;
221  CloudType::move(td, this->db().time().deltaTValue());
222 
223  // finalise and free memory
224  dampingModel_->cacheFields(false);
225  }
226 
227 
228  // Packing
229  // ~~~~~~~
230 
231  if (packingModel_->active())
232  {
233  // same procedure as for damping
234  td.updateAverages(*this);
235  packingModel_->cacheFields(true);
236  td.part() = TrackData::tpPackingNoTrack;
237  CloudType::move(td, this->db().time().deltaTValue());
238  td.part() = TrackData::tpCorrectTrack;
239  CloudType::move(td, this->db().time().deltaTValue());
240  packingModel_->cacheFields(false);
241  }
242 
243 
244  // Isotropy
245  // ~~~~~~~~
246 
247  if (isotropyModel_->active())
248  {
249  // update averages
250  td.updateAverages(*this);
251 
252  // apply isotropy model
253  isotropyModel_->calculate();
254  }
255 
256 
257  // Final
258  // ~~~~~
259 
260  // update cell occupancy
261  this->updateCellOccupancy();
262 
263  // switch forces back on
264  this->forces().setCalcNonCoupled(true);
265  this->forces().setCalcCoupled(this->solution().coupled());
266 }
267 
268 
269 template<class CloudType>
271 {
272  CloudType::info();
273 
274  tmp<volScalarField> alpha = this->theta();
275 
276  const scalar alphaMin = gMin(alpha().internalField());
277  const scalar alphaMax = gMax(alpha().internalField());
278 
279  Info<< " Min cell volume fraction = " << alphaMin << endl;
280  Info<< " Max cell volume fraction = " << alphaMax << endl;
281 
282  if (alphaMax < SMALL)
283  {
284  return;
285  }
286 
287  scalar nMin = GREAT;
288 
289  forAll(this->mesh().cells(), cellI)
290  {
291  const label n = this->cellOccupancy()[cellI].size();
292 
293  if (n > 0)
294  {
295  const scalar nPack = n*alphaMax/alpha()[cellI];
296 
297  if (nPack < nMin)
298  {
299  nMin = nPack;
300  }
301  }
302  }
303 
304  reduce(nMin, minOp<scalar>());
305 
306  Info<< " Min dense number of parcels = " << nMin << endl;
307 }
308 
309 
310 // ************************************************************************* //
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Adds MPPIC modelling to kinematic clouds.
Definition: MPPICCloud.H:66
virtual ~MPPICCloud()
Destructor.
Definition: MPPICCloud.C:144
Base class for collisional damping models.
Definition: MPPICCloud.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void evolve()
Evolve the cloud.
Definition: MPPICCloud.C:172
rhoEqn solve()
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void motion(TrackData &td)
Particle motion.
Definition: MPPICCloud.C:186
This function object calculates the forces and moments by integrating the pressure and skin-friction ...
Definition: forces.H:206
void storeState()
Store the current cloud state.
Definition: MPPICCloud.C:151
A class for handling words, derived from string.
Definition: word.H:59
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
void info()
I-O.
Definition: MPPICCloud.C:270
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
DSMCCloud< dsmcParcel > CloudType
messageStream Info
Base class for packing models.
Definition: MPPICCloud.H:53
dynamicFvMesh & mesh
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
void setModels()
Set cloud sub-models.
Definition: MPPICCloud.C:36
label n
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
void readFields(const Mesh &mesh, const IOobjectList &objects, PtrList< GeoField > &fields)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
autoPtr< IsotropyModel< MPPICCloud< CloudType > > > isotropyModel_
Exchange model.
Definition: MPPICCloud.H:116
#define forAll(list, i)
Definition: UList.H:421
autoPtr< DampingModel< MPPICCloud< CloudType > > > dampingModel_
Damping model.
Definition: MPPICCloud.H:112
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
Type gMin(const FieldField< Field, Type > &f)
void restoreState()
Reset the current cloud to the previously stored state.
Definition: MPPICCloud.C:164
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
const cellShapeList & cells
error FatalError
const List< DynamicList< molecule * > > & cellOccupancy
autoPtr< PackingModel< MPPICCloud< CloudType > > > packingModel_
Packing model.
Definition: MPPICCloud.H:108
Type gMax(const FieldField< Field, Type > &f)
A class for managing temporary objects.
Definition: PtrList.H:118
Base class for collisional return-to-isotropy models.
Definition: MPPICCloud.H:59
conserve internalField()+