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-2020 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 )
36 :
37  ParcelType::trackingData(cloud),
38  volumeAverage_
39  (
41  (
42  IOobject
43  (
44  cloud.name() + ":volumeAverage",
45  cloud.db().time().timeName(),
46  cloud.mesh()
47  ),
48  cloud.solution().dict(),
49  cloud.mesh()
50  )
51  ),
52  radiusAverage_
53  (
55  (
56  IOobject
57  (
58  cloud.name() + ":radiusAverage",
59  cloud.db().time().timeName(),
60  cloud.mesh()
61  ),
62  cloud.solution().dict(),
63  cloud.mesh()
64  )
65  ),
66  rhoAverage_
67  (
69  (
70  IOobject
71  (
72  cloud.name() + ":rhoAverage",
73  cloud.db().time().timeName(),
74  cloud.mesh()
75  ),
76  cloud.solution().dict(),
77  cloud.mesh()
78  )
79  ),
80  uAverage_
81  (
83  (
84  IOobject
85  (
86  cloud.name() + ":uAverage",
87  cloud.db().time().timeName(),
88  cloud.mesh()
89  ),
90  cloud.solution().dict(),
91  cloud.mesh()
92  )
93  ),
94  uSqrAverage_
95  (
97  (
98  IOobject
99  (
100  cloud.name() + ":uSqrAverage",
101  cloud.db().time().timeName(),
102  cloud.mesh()
103  ),
104  cloud.solution().dict(),
105  cloud.mesh()
106  )
107  ),
108  frequencyAverage_
109  (
111  (
112  IOobject
113  (
114  cloud.name() + ":frequencyAverage",
115  cloud.db().time().timeName(),
116  cloud.mesh()
117  ),
118  cloud.solution().dict(),
119  cloud.mesh()
120  )
121  ),
122  massAverage_
123  (
125  (
126  IOobject
127  (
128  cloud.name() + ":massAverage",
129  cloud.db().time().timeName(),
130  cloud.mesh()
131  ),
132  cloud.solution().dict(),
133  cloud.mesh()
134  )
135  ),
136  part_(tpPredictTrack)
137 {}
138 
139 
140 template<class ParcelType>
141 template<class TrackCloudType>
143 (
144  const TrackCloudType& cloud
145 )
146 {
147  // zero the sums
148  volumeAverage_() = 0;
149  radiusAverage_() = 0;
150  rhoAverage_() = 0;
151  uAverage_() = Zero;
152  uSqrAverage_() = 0;
153  frequencyAverage_() = 0;
154  massAverage_() = 0;
155 
156  // temporary weights
157  autoPtr<AveragingMethod<scalar>> weightAveragePtr
158  (
160  (
161  IOobject
162  (
163  cloud.name() + ":weightAverage",
164  cloud.db().time().timeName(),
165  cloud.mesh()
166  ),
167  cloud.solution().dict(),
168  cloud.mesh()
169  )
170  );
171  AveragingMethod<scalar>& weightAverage = weightAveragePtr();
172 
173  // averaging sums
174  forAllConstIter(typename TrackCloudType, cloud, iter)
175  {
176  const typename TrackCloudType::parcelType& p = iter();
177  const tetIndices tetIs = p.currentTetIndices();
178 
179  const scalar m = p.nParticle()*p.mass();
180 
181  volumeAverage_->add(p.coordinates(), tetIs, p.nParticle()*p.volume());
182  rhoAverage_->add(p.coordinates(), tetIs, m*p.rho());
183  uAverage_->add(p.coordinates(), tetIs, m*p.U());
184  massAverage_->add(p.coordinates(), tetIs, m);
185  }
186  volumeAverage_->average();
187  massAverage_->average();
188  rhoAverage_->average(massAverage_);
189  uAverage_->average(massAverage_);
190 
191  // squared velocity deviation
192  forAllConstIter(typename TrackCloudType, cloud, iter)
193  {
194  const typename TrackCloudType::parcelType& p = iter();
195  const tetIndices tetIs = p.currentTetIndices();
196 
197  const vector u = uAverage_->interpolate(p.coordinates(), tetIs);
198 
199  uSqrAverage_->add
200  (
201  p.coordinates(),
202  tetIs,
203  p.nParticle()*p.mass()*magSqr(p.U() - u)
204  );
205  }
206  uSqrAverage_->average(massAverage_);
207 
208  // sauter mean radius
209  radiusAverage_() = volumeAverage_();
210  weightAverage = 0;
211  forAllConstIter(typename TrackCloudType, cloud, iter)
212  {
213  const typename TrackCloudType::parcelType& p = iter();
214  const tetIndices tetIs = p.currentTetIndices();
215 
216  weightAverage.add
217  (
218  p.coordinates(),
219  tetIs,
220  p.nParticle()*pow(p.volume(), 2.0/3.0)
221  );
222  }
223  weightAverage.average();
224  radiusAverage_->average(weightAverage);
225 
226  // collision frequency
227  weightAverage = 0;
228  forAllConstIter(typename TrackCloudType, cloud, iter)
229  {
230  const typename TrackCloudType::parcelType& p = iter();
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 
247 template<class ParcelType>
250 {
251  return part_;
252 }
253 
254 
255 template<class ParcelType>
258 {
259  return part_;
260 }
261 
262 
263 // ************************************************************************* //
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)
trackingData(const TrackCloudType &cloud)
Constructors.
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