MPPICParcelTrackingDataI.H
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-2017 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 "AveragingMethod.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 template<class ParcelType>
31 template<class CloudType>
33 (
34  CloudType& cloud,
35  trackPart part
36 )
37 :
38  ParcelType::template TrackingData<CloudType>(cloud),
39  volumeAverage_
40  (
42  (
43  IOobject
44  (
45  cloud.name() + ":volumeAverage",
46  cloud.db().time().timeName(),
47  cloud.mesh()
48  ),
49  cloud.solution().dict(),
50  cloud.mesh()
51  )
52  ),
53  radiusAverage_
54  (
56  (
57  IOobject
58  (
59  cloud.name() + ":radiusAverage",
60  cloud.db().time().timeName(),
61  cloud.mesh()
62  ),
63  cloud.solution().dict(),
64  cloud.mesh()
65  )
66  ),
67  rhoAverage_
68  (
70  (
71  IOobject
72  (
73  cloud.name() + ":rhoAverage",
74  cloud.db().time().timeName(),
75  cloud.mesh()
76  ),
77  cloud.solution().dict(),
78  cloud.mesh()
79  )
80  ),
81  uAverage_
82  (
84  (
85  IOobject
86  (
87  cloud.name() + ":uAverage",
88  cloud.db().time().timeName(),
89  cloud.mesh()
90  ),
91  cloud.solution().dict(),
92  cloud.mesh()
93  )
94  ),
95  uSqrAverage_
96  (
98  (
99  IOobject
100  (
101  cloud.name() + ":uSqrAverage",
102  cloud.db().time().timeName(),
103  cloud.mesh()
104  ),
105  cloud.solution().dict(),
106  cloud.mesh()
107  )
108  ),
109  frequencyAverage_
110  (
112  (
113  IOobject
114  (
115  cloud.name() + ":frequencyAverage",
116  cloud.db().time().timeName(),
117  cloud.mesh()
118  ),
119  cloud.solution().dict(),
120  cloud.mesh()
121  )
122  ),
123  massAverage_
124  (
126  (
127  IOobject
128  (
129  cloud.name() + ":massAverage",
130  cloud.db().time().timeName(),
131  cloud.mesh()
132  ),
133  cloud.solution().dict(),
134  cloud.mesh()
135  )
136  ),
137  part_(part)
138 {}
139 
140 
141 template<class ParcelType>
142 template<class CloudType>
143 inline void
145 (
146  CloudType& cloud
147 )
148 {
149  // zero the sums
150  volumeAverage_() = 0;
151  radiusAverage_() = 0;
152  rhoAverage_() = 0;
153  uAverage_() = Zero;
154  uSqrAverage_() = 0;
155  frequencyAverage_() = 0;
156  massAverage_() = 0;
157 
158  // temporary weights
159  autoPtr<AveragingMethod<scalar>> weightAveragePtr
160  (
162  (
163  IOobject
164  (
165  cloud.name() + ":weightAverage",
166  cloud.db().time().timeName(),
167  cloud.mesh()
168  ),
169  cloud.solution().dict(),
170  cloud.mesh()
171  )
172  );
173  AveragingMethod<scalar>& weightAverage = weightAveragePtr();
174 
175  // averaging sums
176  forAllConstIter(typename CloudType, cloud, iter)
177  {
178  const typename CloudType::parcelType& p = iter();
179  const tetIndices tetIs = p.currentTetIndices();
180 
181  const scalar m = p.nParticle()*p.mass();
182 
183  volumeAverage_->add(p.coordinates(), tetIs, p.nParticle()*p.volume());
184  rhoAverage_->add(p.coordinates(), tetIs, m*p.rho());
185  uAverage_->add(p.coordinates(), tetIs, m*p.U());
186  massAverage_->add(p.coordinates(), tetIs, m);
187  }
188  volumeAverage_->average();
189  massAverage_->average();
190  rhoAverage_->average(massAverage_);
191  uAverage_->average(massAverage_);
192 
193  // squared velocity deviation
194  forAllConstIter(typename CloudType, cloud, iter)
195  {
196  const typename CloudType::parcelType& p = iter();
197  const tetIndices tetIs = p.currentTetIndices();
198 
199  const vector u = uAverage_->interpolate(p.coordinates(), tetIs);
200 
201  uSqrAverage_->add
202  (
203  p.coordinates(),
204  tetIs,
205  p.nParticle()*p.mass()*magSqr(p.U() - u)
206  );
207  }
208  uSqrAverage_->average(massAverage_);
209 
210  // sauter mean radius
211  radiusAverage_() = volumeAverage_();
212  weightAverage = 0;
213  forAllConstIter(typename CloudType, cloud, iter)
214  {
215  const typename CloudType::parcelType& p = iter();
216  const tetIndices tetIs = p.currentTetIndices();
217 
218  weightAverage.add
219  (
220  p.coordinates(),
221  tetIs,
222  p.nParticle()*pow(p.volume(), 2.0/3.0)
223  );
224  }
225  weightAverage.average();
226  radiusAverage_->average(weightAverage);
227 
228  // collision frequency
229  weightAverage = 0;
230  forAllConstIter(typename CloudType, cloud, iter)
231  {
232  const typename CloudType::parcelType& p = iter();
233  const tetIndices tetIs = p.currentTetIndices();
234 
235  const scalar a = volumeAverage_->interpolate(p.coordinates(), tetIs);
236  const scalar r = radiusAverage_->interpolate(p.coordinates(), tetIs);
237  const vector u = uAverage_->interpolate(p.coordinates(), tetIs);
238 
239  const scalar f = 0.75*a/pow3(r)*sqr(0.5*p.d() + r)*mag(p.U() - u);
240 
241  frequencyAverage_->add(p.coordinates(), tetIs, p.nParticle()*f*f);
242 
243  weightAverage.add(p.coordinates(), tetIs, p.nParticle()*f);
244  }
245  frequencyAverage_->average(weightAverage);
246 }
247 
248 
249 template<class ParcelType>
250 template<class CloudType>
254 {
255  return part_;
256 }
257 
258 
259 template<class ParcelType>
260 template<class CloudType>
264 {
265  return part_;
266 }
267 
268 
269 // ************************************************************************* //
const word & name() const
Return name.
Definition: IOobject.H:291
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:644
virtual void add(const barycentric &coordinates, const tetIndices &tetIs, const Type &value)=0
Member Functions.
Wrapper around kinematic parcel types to add MPPIC modelling.
Definition: MPPICParcel.H:52
Base class for lagrangian averaging methods.
Definition: MPPICParcel.H:55
static const zero Zero
Definition: zero.H:91
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:218
dimensioned< scalar > magSqr(const dimensioned< Type > &)
const Time & time() const
Return time.
labelList f(nPoints)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const fvMesh & mesh() const
Return refernce to the mesh.
Definition: DSMCCloudI.H:41
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual void average()
Calculate the average.
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
volScalarField & p
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:331
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69