MPPICParcelTrackingDataI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2018 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 TrackCloudType>
33 (
34  const TrackCloudType& cloud,
35  trackPart part
36 )
37 :
38  ParcelType::trackingData(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 TrackCloudType>
144 (
145  const TrackCloudType& cloud
146 )
147 {
148  // zero the sums
149  volumeAverage_() = 0;
150  radiusAverage_() = 0;
151  rhoAverage_() = 0;
152  uAverage_() = Zero;
153  uSqrAverage_() = 0;
154  frequencyAverage_() = 0;
155  massAverage_() = 0;
156 
157  // temporary weights
158  autoPtr<AveragingMethod<scalar>> weightAveragePtr
159  (
161  (
162  IOobject
163  (
164  cloud.name() + ":weightAverage",
165  cloud.db().time().timeName(),
166  cloud.mesh()
167  ),
168  cloud.solution().dict(),
169  cloud.mesh()
170  )
171  );
172  AveragingMethod<scalar>& weightAverage = weightAveragePtr();
173 
174  // averaging sums
175  forAllConstIter(typename TrackCloudType, cloud, iter)
176  {
177  const typename TrackCloudType::parcelType& p = iter();
178  const tetIndices tetIs = p.currentTetIndices();
179 
180  const scalar m = p.nParticle()*p.mass();
181 
182  volumeAverage_->add(p.coordinates(), tetIs, p.nParticle()*p.volume());
183  rhoAverage_->add(p.coordinates(), tetIs, m*p.rho());
184  uAverage_->add(p.coordinates(), tetIs, m*p.U());
185  massAverage_->add(p.coordinates(), tetIs, m);
186  }
187  volumeAverage_->average();
188  massAverage_->average();
189  rhoAverage_->average(massAverage_);
190  uAverage_->average(massAverage_);
191 
192  // squared velocity deviation
193  forAllConstIter(typename TrackCloudType, cloud, iter)
194  {
195  const typename TrackCloudType::parcelType& p = iter();
196  const tetIndices tetIs = p.currentTetIndices();
197 
198  const vector u = uAverage_->interpolate(p.coordinates(), tetIs);
199 
200  uSqrAverage_->add
201  (
202  p.coordinates(),
203  tetIs,
204  p.nParticle()*p.mass()*magSqr(p.U() - u)
205  );
206  }
207  uSqrAverage_->average(massAverage_);
208 
209  // sauter mean radius
210  radiusAverage_() = volumeAverage_();
211  weightAverage = 0;
212  forAllConstIter(typename TrackCloudType, cloud, iter)
213  {
214  const typename TrackCloudType::parcelType& p = iter();
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  forAllConstIter(typename TrackCloudType, cloud, iter)
230  {
231  const typename TrackCloudType::parcelType& p = iter();
232  const tetIndices tetIs = p.currentTetIndices();
233 
234  const scalar a = volumeAverage_->interpolate(p.coordinates(), tetIs);
235  const scalar r = radiusAverage_->interpolate(p.coordinates(), tetIs);
236  const vector u = uAverage_->interpolate(p.coordinates(), tetIs);
237 
238  const scalar f = 0.75*a/pow3(r)*sqr(0.5*p.d() + r)*mag(p.U() - u);
239 
240  frequencyAverage_->add(p.coordinates(), tetIs, p.nParticle()*f*f);
241 
242  weightAverage.add(p.coordinates(), tetIs, p.nParticle()*f);
243  }
244  frequencyAverage_->average(weightAverage);
245 }
246 
247 
248 template<class ParcelType>
251 {
252  return part_;
253 }
254 
255 
256 template<class ParcelType>
259 {
260  return part_;
261 }
262 
263 
264 // ************************************************************************* //
trackingData(const TrackCloudType &cloud, trackPart part=tpLinearTrack)
Constructors.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual void add(const barycentric &coordinates, const tetIndices &tetIs, const Type &value)=0
Member Functions.
Base class for lagrangian averaging methods.
Definition: MPPICParcel.H:55
static const zero Zero
Definition: zero.H:97
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
void updateAverages(const TrackCloudType &cloud)
Update the MPPIC averages.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
labelList f(nPoints)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92