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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2017 OpenFOAM Foundation
9 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
29#include "AveragingMethod.H"
30
31// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32
33template<class ParcelType>
34template<class TrackCloudType>
36(
37 const TrackCloudType& cloud,
38 trackPart part
39)
40:
41 ParcelType::trackingData(cloud),
42 volumeAverage_
43 (
44 AveragingMethod<scalar>::New
45 (
47 (
48 cloud.name() + ":volumeAverage",
49 cloud.db().time().timeName(),
50 cloud.mesh()
51 ),
53 cloud.mesh()
54 )
55 ),
56 radiusAverage_
57 (
58 AveragingMethod<scalar>::New
59 (
61 (
62 cloud.name() + ":radiusAverage",
63 cloud.db().time().timeName(),
64 cloud.mesh()
65 ),
67 cloud.mesh()
68 )
69 ),
70 rhoAverage_
71 (
72 AveragingMethod<scalar>::New
73 (
75 (
76 cloud.name() + ":rhoAverage",
77 cloud.db().time().timeName(),
78 cloud.mesh()
79 ),
81 cloud.mesh()
82 )
83 ),
84 uAverage_
85 (
87 (
89 (
90 cloud.name() + ":uAverage",
91 cloud.db().time().timeName(),
92 cloud.mesh()
93 ),
95 cloud.mesh()
96 )
97 ),
98 uSqrAverage_
99 (
100 AveragingMethod<scalar>::New
101 (
103 (
104 cloud.name() + ":uSqrAverage",
105 cloud.db().time().timeName(),
106 cloud.mesh()
107 ),
108 cloud.solution().dict(),
109 cloud.mesh()
110 )
111 ),
112 frequencyAverage_
113 (
114 AveragingMethod<scalar>::New
115 (
117 (
118 cloud.name() + ":frequencyAverage",
119 cloud.db().time().timeName(),
120 cloud.mesh()
121 ),
122 cloud.solution().dict(),
123 cloud.mesh()
124 )
125 ),
126 massAverage_
127 (
128 AveragingMethod<scalar>::New
129 (
131 (
132 cloud.name() + ":massAverage",
133 cloud.db().time().timeName(),
134 cloud.mesh()
135 ),
136 cloud.solution().dict(),
137 cloud.mesh()
138 )
139 ),
140 part_(part)
141{}
142
143
144template<class ParcelType>
145template<class TrackCloudType>
147(
148 const TrackCloudType& cloud
149)
150{
151 // zero the sums
152 volumeAverage_() = 0;
153 radiusAverage_() = 0;
154 rhoAverage_() = 0;
155 uAverage_() = Zero;
156 uSqrAverage_() = 0;
157 frequencyAverage_() = 0;
158 massAverage_() = 0;
159
160 // temporary weights
161 autoPtr<AveragingMethod<scalar>> weightAveragePtr
162 (
164 (
166 (
167 cloud.name() + ":weightAverage",
168 cloud.db().time().timeName(),
169 cloud.mesh()
170 ),
171 cloud.solution().dict(),
172 cloud.mesh()
173 )
174 );
175 AveragingMethod<scalar>& weightAverage = weightAveragePtr();
176
177 // averaging sums
178 for (const typename TrackCloudType::parcelType& p : cloud)
179 {
180 const tetIndices tetIs = p.currentTetIndices();
181
182 const scalar m = p.nParticle()*p.mass();
183
184 volumeAverage_->add(p.coordinates(), tetIs, p.nParticle()*p.volume());
185 rhoAverage_->add(p.coordinates(), tetIs, m*p.rho());
186 uAverage_->add(p.coordinates(), tetIs, m*p.U());
187 massAverage_->add(p.coordinates(), tetIs, m);
188 }
189 volumeAverage_->average();
190 massAverage_->average();
191 rhoAverage_->average(*massAverage_);
192 uAverage_->average(*massAverage_);
193
194 // squared velocity deviation
195 for (const typename TrackCloudType::parcelType& p : cloud)
196 {
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 for (const typename TrackCloudType::parcelType& p : cloud)
214 {
215 const tetIndices tetIs = p.currentTetIndices();
216
217 weightAverage.add
218 (
219 p.coordinates(),
220 tetIs,
221 p.nParticle()*pow(p.volume(), 2.0/3.0)
222 );
223 }
224 weightAverage.average();
225 radiusAverage_->average(weightAverage);
226
227 // collision frequency
228 weightAverage = 0;
229 for (const typename TrackCloudType::parcelType& p : cloud)
230 {
231 const tetIndices tetIs = p.currentTetIndices();
232
233 const scalar a = volumeAverage_->interpolate(p.coordinates(), tetIs);
234 const scalar r = radiusAverage_->interpolate(p.coordinates(), tetIs);
235 const vector u = uAverage_->interpolate(p.coordinates(), tetIs);
236
237 const scalar f = 0.75*a/pow3(r)*sqr(0.5*p.d() + r)*mag(p.U() - u);
238
239 frequencyAverage_->add(p.coordinates(), tetIs, p.nParticle()*f*f);
240
241 weightAverage.add(p.coordinates(), tetIs, p.nParticle()*f);
242 }
243 frequencyAverage_->average(weightAverage);
244}
245
246
247template<class ParcelType>
250{
251 return part_;
252}
253
254
255template<class ParcelType>
258{
259 return part_;
260}
261
262
263// ************************************************************************* //
Base class for lagrangian averaging methods.
virtual void add(const barycentric &coordinates, const tetIndices &tetIs, const Type &value)=0
Member Functions.
virtual void average()
Calculate the average.
dimensioned< Type > average() const
Calculate and return arithmetic average.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:500
void updateAverages(const TrackCloudType &cloud)
Update the MPPIC averages.
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
Class used to pass data into container.
const Time & time() const noexcept
Return time registry.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:66
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:84
volScalarField & p
dynamicFvMesh & mesh
word timeName
Definition: getTimeIndex.H:3
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
labelList f(nPoints)
dictionary dict