KinematicParcelTrackingDataI.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 template<class ParcelType>
34 template<class TrackCloudType>
36 (
37  const TrackCloudType& cloud,
38  trackPart part
39 )
40 :
41  ParcelType::trackingData(cloud),
42  rhoInterp_
43  (
44  interpolation<scalar>::New
45  (
46  cloud.solution().interpolationSchemes(),
47  cloud.rho()
48  )
49  ),
50  UInterp_
51  (
53  (
54  cloud.solution().interpolationSchemes(),
55  cloud.U()
56  )
57  ),
58  muInterp_
59  (
60  interpolation<scalar>::New
61  (
62  cloud.solution().interpolationSchemes(),
63  cloud.mu()
64  )
65  ),
66  rhoc_(Zero),
67  Uc_(Zero),
68  muc_(Zero),
69 
70  volumeAverage_
71  (
72  AveragingMethod<scalar>::New
73  (
74  IOobject
75  (
76  cloud.name() + ":volumeAverage",
77  cloud.db().time().timeName(),
78  cloud.mesh()
79  ),
80  cloud.solution().dict(),
81  cloud.mesh()
82  )
83  ),
84  radiusAverage_
85  (
86  AveragingMethod<scalar>::New
87  (
88  IOobject
89  (
90  cloud.name() + ":radiusAverage",
91  cloud.db().time().timeName(),
92  cloud.mesh()
93  ),
94  cloud.solution().dict(),
95  cloud.mesh()
96  )
97  ),
98  rhoAverage_
99  (
100  AveragingMethod<scalar>::New
101  (
102  IOobject
103  (
104  cloud.name() + ":rhoAverage",
105  cloud.db().time().timeName(),
106  cloud.mesh()
107  ),
108  cloud.solution().dict(),
109  cloud.mesh()
110  )
111  ),
112  uAverage_
113  (
115  (
116  IOobject
117  (
118  cloud.name() + ":uAverage",
119  cloud.db().time().timeName(),
120  cloud.mesh()
121  ),
122  cloud.solution().dict(),
123  cloud.mesh()
124  )
125  ),
126  uSqrAverage_
127  (
128  AveragingMethod<scalar>::New
129  (
130  IOobject
131  (
132  cloud.name() + ":uSqrAverage",
133  cloud.db().time().timeName(),
134  cloud.mesh()
135  ),
136  cloud.solution().dict(),
137  cloud.mesh()
138  )
139  ),
140  frequencyAverage_
141  (
142  AveragingMethod<scalar>::New
143  (
144  IOobject
145  (
146  cloud.name() + ":frequencyAverage",
147  cloud.db().time().timeName(),
148  cloud.mesh()
149  ),
150  cloud.solution().dict(),
151  cloud.mesh()
152  )
153  ),
154  massAverage_
155  (
156  AveragingMethod<scalar>::New
157  (
158  IOobject
159  (
160  cloud.name() + ":massAverage",
161  cloud.db().time().timeName(),
162  cloud.mesh()
163  ),
164  cloud.solution().dict(),
165  cloud.mesh()
166  )
167  ),
168 
169  g_(cloud.g().value()),
170  part_(part)
171 {}
172 
173 
174 template<class ParcelType>
177 {
178  return *rhoInterp_;
179 }
180 
181 
182 template<class ParcelType>
185 {
186  return *UInterp_;
187 }
188 
189 
190 template<class ParcelType>
193 {
194  return *muInterp_;
195 }
196 
197 
198 template<class ParcelType>
199 inline const Foam::vector&
201 {
202  return g_;
203 }
204 
205 
206 template<class ParcelType>
207 inline Foam::scalar
209 {
210  return rhoc_;
211 }
212 
213 
214 template<class ParcelType>
216 {
217  return rhoc_;
218 }
219 
220 
221 template<class ParcelType>
222 inline const Foam::vector&
224 {
225  return Uc_;
226 }
227 
228 
229 template<class ParcelType>
231 {
232  return Uc_;
233 }
234 
235 
236 template<class ParcelType>
238 {
239  return muc_;
240 }
241 
242 
243 template<class ParcelType>
245 {
246  return muc_;
247 }
248 
249 
250 template<class ParcelType>
253 {
254  return part_;
255 }
256 
257 
258 template<class ParcelType>
261 {
262  return part_;
263 }
264 
265 
266 template<class ParcelType>
267 template<class TrackCloudType>
270 (
271  const TrackCloudType& cloud
272 )
273 {
274  // zero the sums
275  volumeAverage_() = 0;
276  radiusAverage_() = 0;
277  rhoAverage_() = 0;
278  uAverage_() = Zero;
279  uSqrAverage_() = 0;
280  frequencyAverage_() = 0;
281  massAverage_() = 0;
282 
283  // temporary weights
284  autoPtr<AveragingMethod<scalar>> weightAveragePtr
285  (
287  (
288  IOobject
289  (
290  cloud.name() + ":weightAverage",
291  cloud.db().time().timeName(),
292  cloud.mesh()
293  ),
294  cloud.solution().dict(),
295  cloud.mesh()
296  )
297  );
298  AveragingMethod<scalar>& weightAverage = weightAveragePtr();
299 
300  // averaging sums
301  for (const typename TrackCloudType::parcelType& p : cloud)
302  {
303  const tetIndices tetIs = p.currentTetIndices();
304 
305  const scalar m = p.nParticle()*p.mass();
306 
307  volumeAverage_->add(p.coordinates(), tetIs, p.nParticle()*p.volume());
308  rhoAverage_->add(p.coordinates(), tetIs, m*p.rho());
309  uAverage_->add(p.coordinates(), tetIs, m*p.U());
310  massAverage_->add(p.coordinates(), tetIs, m);
311  }
312  volumeAverage_->average();
313  massAverage_->average();
314  rhoAverage_->average(*massAverage_);
315  uAverage_->average(*massAverage_);
316 
317  // squared velocity deviation
318  for (const typename TrackCloudType::parcelType& p : cloud)
319  {
320  const tetIndices tetIs = p.currentTetIndices();
321 
322  const vector u = uAverage_->interpolate(p.coordinates(), tetIs);
323 
324  uSqrAverage_->add
325  (
326  p.coordinates(),
327  tetIs,
328  p.nParticle()*p.mass()*magSqr(p.U() - u)
329  );
330  }
331  uSqrAverage_->average(*massAverage_);
332 
333  // sauter mean radius
334  radiusAverage_() = volumeAverage_();
335  weightAverage = 0;
336  for (const typename TrackCloudType::parcelType& p : cloud)
337  {
338  const tetIndices tetIs = p.currentTetIndices();
339 
340  weightAverage.add
341  (
342  p.coordinates(),
343  tetIs,
344  p.nParticle()*pow(p.volume(), 2.0/3.0)
345  );
346  }
347  weightAverage.average();
348  radiusAverage_->average(weightAverage);
349 
350  // collision frequency
351  weightAverage = 0;
352  for (const typename TrackCloudType::parcelType& p : cloud)
353  {
354  const tetIndices tetIs = p.currentTetIndices();
355 
356  const scalar a = volumeAverage_->interpolate(p.coordinates(), tetIs);
357  const scalar r = radiusAverage_->interpolate(p.coordinates(), tetIs);
358  const vector u = uAverage_->interpolate(p.coordinates(), tetIs);
359 
360  const scalar f = 0.75*a/pow3(r)*sqr(0.5*p.d() + r)*mag(p.U() - u);
361 
362  frequencyAverage_->add(p.coordinates(), tetIs, p.nParticle()*f*f);
363 
364  weightAverage.add(p.coordinates(), tetIs, p.nParticle()*f);
365  }
366  frequencyAverage_->average(weightAverage);
367 }
368 
369 // ************************************************************************* //
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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.
const interpolation< vector > & UInterp() const
Return const access to the interpolator for continuous phase velocity field.
virtual void add(const barycentric &coordinates, const tetIndices &tetIs, const Type &value)=0
Member Functions.
word timeName
Definition: getTimeIndex.H:3
trackPart part() const
Return the part of the tracking operation taking place.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const Time & time() const noexcept
Return time registry.
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
Base class for lagrangian averaging methods.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:450
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:78
Vector< scalar > vector
Definition: vector.H:57
const uniformDimensionedVectorField & g
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
scalar rhoc() const
Return the continuous phase density.
scalar muc() const
Return the continuous phase viscosity.
labelList f(nPoints)
trackingData(const TrackCloudType &cloud, trackPart part=tpLinearTrack)
Construct from components.
const dimensionedScalar mu
Atomic mass unit.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
virtual void average()
Calculate the average.
const interpolation< scalar > & rhoInterp() const
Return const access to the interpolator for continuous phase density field.
Abstract base class for volume field interpolation.
const interpolation< scalar > & muInterp() const
Return const access to the interpolator for continuous phase dynamic viscosity field.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:92
volScalarField & p
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
const vector & Uc() const
Return the continuous phase velocity.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
void updateAverages(const TrackCloudType &cloud)
Update the MPPIC averages.