ReactingParcel.C
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) 2011-2022 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 "ReactingParcel.H"
27 #include "specie.H"
28 #include "CompositionModel.H"
29 #include "PhaseChangeModel.H"
30 #include "mathematicalConstants.H"
31 
32 using namespace Foam::constant::mathematical;
33 
34 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
35 
36 template<class ParcelType>
37 template<class TrackCloudType>
39 (
40  TrackCloudType& cloud,
41  trackingData& td,
42  const scalar dt,
43  const scalar Re,
44  const scalar Pr,
45  const scalar Ts,
46  const scalar nus,
47  const scalar d,
48  const scalar T,
49  const scalar mass,
50  const label idPhase,
51  const scalar YPhase,
52  const scalarField& Y,
53  scalarField& dMassPC,
54  scalar& Sh,
55  scalar& N,
56  scalar& NCpW,
57  scalarField& Cs
58 )
59 {
60  typedef typename TrackCloudType::thermoCloudType thermoCloudType;
62  cloud.composition();
63 
64  typedef typename TrackCloudType::reactingCloudType reactingCloudType;
65  PhaseChangeModel<reactingCloudType>& phaseChange = cloud.phaseChange();
66 
67  if (YPhase < small)
68  {
69  return;
70  }
71 
72  scalarField X(composition.liquids().X(Y));
73 
74  scalar Tvap = phaseChange.Tvap(X);
75 
76  if (T < Tvap)
77  {
78  return;
79  }
80 
81  const scalar TMax = phaseChange.TMax(td.pc(), X);
82  const scalar Tdash = min(T, TMax);
83  const scalar Tsdash = min(Ts, TMax);
84 
85  const typename TrackCloudType::parcelType& p =
86  static_cast<const typename TrackCloudType::parcelType&>(*this);
87  typename TrackCloudType::parcelType::trackingData& ttd =
88  static_cast<typename TrackCloudType::parcelType::trackingData&>(td);
89 
90  // Calculate mass transfer due to phase change
91  phaseChange.calculate
92  (
93  p,
94  ttd,
95  dt,
96  Re,
97  Pr,
98  d,
99  nus,
100  Tdash,
101  Tsdash,
102  td.pc(),
103  td.Tc(),
104  X,
105  dMassPC
106  );
107 
108  // Limit phase change mass by availability of each specie
109  dMassPC = min(mass*YPhase*Y, dMassPC);
110 
111  const scalar dMassTot = sum(dMassPC);
112 
113  // Add to cumulative phase change mass
114  phaseChange.addToPhaseChangeMass(this->nParticle_*dMassTot);
115 
116  forAll(dMassPC, i)
117  {
118  const label cid = composition.localToCarrierId(idPhase, i);
119 
120  const scalar dh = phaseChange.dh(cid, i, td.pc(), Tdash);
121  Sh -= dMassPC[i]*dh/dt;
122  }
123 
124 
125  // Update molar emissions
126  if (cloud.heatTransfer().BirdCorrection())
127  {
128  // Average molecular weight of carrier mix - assumes perfect gas
129  const scalar Wc = td.rhoc()*RR*td.Tc()/td.pc();
130 
131  forAll(dMassPC, i)
132  {
133  const label cid = composition.localToCarrierId(idPhase, i);
134 
135  const scalar Cp = composition.carrier().Cp(cid, td.pc(), Tsdash);
136  const scalar W = composition.carrier().Wi(cid);
137  const scalar Ni = dMassPC[i]/(this->areaS(d)*dt*W);
138 
139  const scalar Dab =
140  composition.liquids().properties()[i].D(td.pc(), Tsdash, Wc);
141 
142  // Molar flux of species coming from the particle (kmol/m^2/s)
143  N += Ni;
144 
145  // Sum of Ni*Cpi*Wi of emission species
146  NCpW += Ni*Cp*W;
147 
148  // Concentrations of emission species
149  Cs[cid] += Ni*d/(2.0*Dab);
150  }
151  }
152 }
153 
154 
155 template<class ParcelType>
157 (
158  const scalar mass0,
159  const scalarField& dMass,
160  scalarField& Y
161 ) const
162 {
163  scalar mass1 = mass0 - sum(dMass);
164 
165  // only update the mass fractions if the new particle mass is finite
166  if (mass1 > rootVSmall)
167  {
168  forAll(Y, i)
169  {
170  Y[i] = (Y[i]*mass0 - dMass[i])/mass1;
171  }
172  }
173 
174  return mass1;
175 }
176 
177 
178 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
179 
180 template<class ParcelType>
182 (
184 )
185 :
186  ParcelType(p),
187  mass0_(p.mass0_),
188  Y_(p.Y_)
189 {}
190 
191 
192 template<class ParcelType>
194 (
196  const polyMesh& mesh
197 )
198 :
199  ParcelType(p, mesh),
200  mass0_(p.mass0_),
201  Y_(p.Y_)
202 {}
203 
204 
205 // * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * * //
206 
207 template<class ParcelType>
208 template<class TrackCloudType>
210 (
211  TrackCloudType& cloud,
212  trackingData& td
213 )
214 {
215  ParcelType::setCellValues(cloud, td);
216 
217  td.pc() = td.pInterp().interpolate
218  (
219  this->coordinates(),
220  this->currentTetIndices()
221  );
222 
223  if (td.pc() < cloud.constProps().pMin())
224  {
225  if (debug)
226  {
228  << "Limiting observed pressure in cell " << this->cell()
229  << " to " << cloud.constProps().pMin() << nl << endl;
230  }
231 
232  td.pc() = cloud.constProps().pMin();
233  }
234 }
235 
236 
237 template<class ParcelType>
238 template<class TrackCloudType>
240 (
241  TrackCloudType& cloud,
242  trackingData& td,
243  const scalar dt
244 )
245 {
246  scalar addedMass = 0.0;
247  scalar maxMassI = 0.0;
248  forAll(cloud.rhoTrans(), i)
249  {
250  scalar dm = cloud.rhoTrans(i)[this->cell()];
251  maxMassI = max(maxMassI, mag(dm));
252  addedMass += dm;
253  }
254 
255  if (maxMassI < rootVSmall)
256  {
257  return;
258  }
259 
260  const scalar massCell = this->massCell(td);
261 
262  td.rhoc() += addedMass/cloud.pMesh().cellVolumes()[this->cell()];
263 
264  const scalar massCellNew = massCell + addedMass;
265  td.Uc() = (td.Uc()*massCell + cloud.UTransRef()[this->cell()])/massCellNew;
266 
267  scalar CpEff = 0.0;
268  forAll(cloud.rhoTrans(), i)
269  {
270  scalar Y = cloud.rhoTrans(i)[this->cell()]/addedMass;
271  CpEff += Y*cloud.composition().carrier().Cp(i, td.pc(), td.Tc());
272  }
273 
274  const scalar Cpc = td.CpInterp().psi()[this->cell()];
275  td.Cpc() = (massCell*Cpc + addedMass*CpEff)/massCellNew;
276 
277  td.Tc() += cloud.hsTransRef()[this->cell()]/(td.Cpc()*massCellNew);
278 
279  if (td.Tc() < cloud.constProps().TMin())
280  {
281  if (debug)
282  {
284  << "Limiting observed temperature in cell " << this->cell()
285  << " to " << cloud.constProps().TMin() << nl << endl;
286  }
287 
288  td.Tc() = cloud.constProps().TMin();
289  }
290 }
291 
292 
293 template<class ParcelType>
294 template<class TrackCloudType>
296 (
297  TrackCloudType& cloud,
298  trackingData& td,
299  const scalar T,
300  const scalarField& Cs,
301  scalar& rhos,
302  scalar& mus,
303  scalar& Prs,
304  scalar& kappas
305 )
306 {
307  // No correction if total concentration of emitted species is small
308  if (!cloud.heatTransfer().BirdCorrection() || (sum(Cs) < small))
309  {
310  return;
311  }
312 
313  const basicSpecieMixture& carrier = cloud.composition().carrier();
314 
315  // Far field carrier molar fractions
316  scalarField Xinf(carrier.species().size());
317 
318  forAll(Xinf, i)
319  {
320  Xinf[i] = carrier.Y(i)[this->cell()]/carrier.Wi(i);
321  }
322  Xinf /= sum(Xinf);
323 
324  // Molar fraction of far field species at particle surface
325  const scalar Xsff = 1.0 - min(sum(Cs)*RR*this->T_/td.pc(), 1.0);
326 
327  // Surface carrier total molar concentration
328  const scalar CsTot = td.pc()/(RR*this->T_);
329 
330  // Surface carrier composition (molar fraction)
331  scalarField Xs(Xinf.size());
332 
333  // Surface carrier composition (mass fraction)
334  scalarField Ys(Xinf.size());
335 
336  forAll(Xs, i)
337  {
338  // Molar concentration of species at particle surface
339  const scalar Csi = Cs[i] + Xsff*Xinf[i]*CsTot;
340 
341  Xs[i] = (2.0*Csi + Xinf[i]*CsTot)/3.0;
342  Ys[i] = Xs[i]*carrier.Wi(i);
343  }
344  Xs /= sum(Xs);
345  Ys /= sum(Ys);
346 
347 
348  rhos = 0;
349  mus = 0;
350  kappas = 0;
351  scalar Cps = 0;
352  scalar sumYiSqrtW = 0;
353  scalar sumYiCbrtW = 0;
354 
355  forAll(Ys, i)
356  {
357  const scalar W = carrier.Wi(i);
358  const scalar sqrtW = sqrt(W);
359  const scalar cbrtW = cbrt(W);
360 
361  rhos += Xs[i]*W;
362  mus += Ys[i]*sqrtW*carrier.mu(i, td.pc(), T);
363  kappas += Ys[i]*cbrtW*carrier.kappa(i, td.pc(), T);
364  Cps += Xs[i]*carrier.Cp(i, td.pc(), T);
365 
366  sumYiSqrtW += Ys[i]*sqrtW;
367  sumYiCbrtW += Ys[i]*cbrtW;
368  }
369 
370  Cps = max(Cps, rootVSmall);
371 
372  rhos *= td.pc()/(RR*T);
373  rhos = max(rhos, rootVSmall);
374 
375  mus /= sumYiSqrtW;
376  mus = max(mus, rootVSmall);
377 
378  kappas /= sumYiCbrtW;
379  kappas = max(kappas, rootVSmall);
380 
381  Prs = Cps*mus/kappas;
382 }
383 
384 
385 template<class ParcelType>
386 template<class TrackCloudType>
388 (
389  TrackCloudType& cloud,
390  trackingData& td,
391  const scalar dt
392 )
393 {
394  typedef typename TrackCloudType::thermoCloudType thermoCloudType;
396  cloud.composition();
397 
398 
399  // Define local properties at beginning of time step
400  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
401 
402  const scalar np0 = this->nParticle_;
403  const scalar d0 = this->d_;
404  const vector& U0 = this->U_;
405  const scalar T0 = this->T_;
406  const scalar mass0 = this->mass();
407 
408 
409  // Calc surface values
410  scalar Ts, rhos, mus, Prs, kappas;
411  this->calcSurfaceValues(cloud, td, T0, Ts, rhos, mus, Prs, kappas);
412  scalar Res = this->Re(rhos, U0, td.Uc(), d0, mus);
413 
414 
415  // Sources
416  // ~~~~~~~
417 
418  // Explicit momentum source for particle
419  vector Su = Zero;
420 
421  // Linearised momentum source coefficient
422  scalar Spu = 0.0;
423 
424  // Momentum transfer from the particle to the carrier phase
425  vector dUTrans = Zero;
426 
427  // Explicit enthalpy source for particle
428  scalar Sh = 0.0;
429 
430  // Linearised enthalpy source coefficient
431  scalar Sph = 0.0;
432 
433  // Sensible enthalpy transfer from the particle to the carrier phase
434  scalar dhsTrans = 0.0;
435 
436 
437  // 1. Compute models that contribute to mass transfer - U, T held constant
438  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
439 
440  // Phase change
441  // ~~~~~~~~~~~~
442 
443  // Mass transfer due to phase change
444  scalarField dMassPC(Y_.size(), 0.0);
445 
446  // Molar flux of species emitted from the particle (kmol/m^2/s)
447  scalar Ne = 0.0;
448 
449  // Sum Ni*Cpi*Wi of emission species
450  scalar NCpW = 0.0;
451 
452  // Surface concentrations of emitted species
453  scalarField Cs(composition.carrier().species().size(), 0.0);
454 
455  // Calc mass and enthalpy transfer due to phase change
456  calcPhaseChange
457  (
458  cloud,
459  td,
460  dt,
461  Res,
462  Prs,
463  Ts,
464  mus/rhos,
465  d0,
466  T0,
467  mass0,
468  0,
469  1.0,
470  Y_,
471  dMassPC,
472  Sh,
473  Ne,
474  NCpW,
475  Cs
476  );
477 
478 
479  // 2. Update the parcel properties due to change in mass
480  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
481 
482  scalarField dMass(dMassPC);
483  scalar mass1 = updateMassFraction(mass0, dMass, Y_);
484 
485  this->Cp_ = composition.Cp(0, Y_, td.pc(), T0);
486 
487  // Update particle density or diameter
488  if (cloud.constProps().constantVolume())
489  {
490  this->rho_ = mass1/this->volume();
491  }
492  else
493  {
494  this->d_ = cbrt(mass1/this->rho_*6.0/pi);
495  }
496 
497  // Remove the particle when mass falls below minimum threshold
498  if (np0*mass1 < cloud.constProps().minParcelMass())
499  {
500  td.keepParticle = false;
501 
502  if (cloud.solution().coupled())
503  {
504  scalar dm = np0*mass0;
505 
506  // Absorb parcel into carrier phase
507  forAll(Y_, i)
508  {
509  scalar dmi = dm*Y_[i];
510  label gid = composition.localToCarrierId(0, i);
511  scalar hs = composition.carrier().Hs(gid, td.pc(), T0);
512 
513  cloud.rhoTrans(gid)[this->cell()] += dmi;
514  cloud.hsTransRef()[this->cell()] += dmi*hs;
515  }
516  cloud.UTransRef()[this->cell()] += dm*U0;
517 
518  cloud.phaseChange().addToPhaseChangeMass(np0*mass1);
519  }
520 
521  return;
522  }
523 
524  // Correct surface values due to emitted species
525  correctSurfaceValues(cloud, td, Ts, Cs, rhos, mus, Prs, kappas);
526  Res = this->Re(rhos, U0, td.Uc(), this->d(), mus);
527 
528 
529  // 3. Compute heat- and momentum transfers
530  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
531 
532  // Heat transfer
533  // ~~~~~~~~~~~~~
534 
535  // Calculate new particle temperature
536  this->T_ =
537  this->calcHeatTransfer
538  (
539  cloud,
540  td,
541  dt,
542  Res,
543  Prs,
544  kappas,
545  NCpW,
546  Sh,
547  dhsTrans,
548  Sph
549  );
550 
551  this->Cp_ = composition.Cp(0, Y_, td.pc(), T0);
552 
553 
554  // Motion
555  // ~~~~~~
556 
557  // Calculate new particle velocity
558  this->U_ =
559  this->calcVelocity(cloud, td, dt, Res, mus, mass1, Su, dUTrans, Spu);
560 
561 
562  // 4. Accumulate carrier phase source terms
563  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
564 
565  if (cloud.solution().coupled())
566  {
567  // Transfer mass lost to carrier mass, momentum and enthalpy sources
568  forAll(dMass, i)
569  {
570  scalar dm = np0*dMass[i];
571  label gid = composition.localToCarrierId(0, i);
572  scalar hs = composition.carrier().Hs(gid, td.pc(), T0);
573 
574  cloud.rhoTrans(gid)[this->cell()] += dm;
575  cloud.UTransRef()[this->cell()] += dm*U0;
576  cloud.hsTransRef()[this->cell()] += dm*hs;
577  }
578 
579  // Update momentum transfer
580  cloud.UTransRef()[this->cell()] += np0*dUTrans;
581  cloud.UCoeffRef()[this->cell()] += np0*Spu;
582 
583  // Update sensible enthalpy transfer
584  cloud.hsTransRef()[this->cell()] += np0*dhsTrans;
585  cloud.hsCoeffRef()[this->cell()] += np0*Sph;
586 
587  // Update radiation fields
588  if (cloud.radiation())
589  {
590  const scalar ap = this->areaP();
591  const scalar T4 = pow4(T0);
592  cloud.radAreaP()[this->cell()] += dt*np0*ap;
593  cloud.radT4()[this->cell()] += dt*np0*T4;
594  cloud.radAreaPT4()[this->cell()] += dt*np0*ap*T4;
595  }
596  }
597 }
598 
599 
600 // * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * * * //
601 
602 #include "ReactingParcelIO.C"
603 
604 // ************************************************************************* //
virtual scalar TMax(const scalar p, const scalarField &X) const
Return maximum/limiting temperature.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Templated phase change model class.
Definition: ReactingCloud.H:57
virtual scalar kappa(const label speciei, const scalar p, const scalar T) const =0
Thermal conductivity [W/m/K].
void calcPhaseChange(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Re, const scalar Pr, const scalar Ts, const scalar nus, const scalar d, const scalar T, const scalar mass, const label idPhase, const scalar YPhase, const scalarField &YComponents, scalarField &dMassPC, scalar &Sh, scalar &N, scalar &NCpW, scalarField &Cs)
Calculate Phase change.
basicSpecieMixture & composition
scalar updateMassFraction(const scalar mass0, const scalarField &dMass, scalarField &Y) const
Update mass fraction.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
scalarField Y_
Mass fractions of mixture [].
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
ReactingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
twoPhaseChangeModel & phaseChange
virtual scalar Wi(const label speciei) const =0
Molecular weight of the given specie [kg/kmol].
label localToCarrierId(const label phaseI, const label id, const bool allowNotFound=false) const
Return carrier id of component given local id.
void addToPhaseChangeMass(const scalar dMass)
Add to phase change mass.
virtual scalar Cp(const label phaseI, const scalarField &Y, const scalar p, const scalar T) const
Return specific heat capacity for the phase phaseI.
Specialisation of basicMixture for a mixture consisting of a number for molecular species...
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
ParcelType::trackingData trackingData
Use base tracking data.
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
virtual scalar Cp(const label speciei, const scalar p, const scalar T) const =0
Heat capacity at constant pressure [J/kg/K].
virtual scalar mu(const label speciei, const scalar p, const scalar T) const =0
Dynamic viscosity [kg/m/s].
mathematical constants.
dimensionedScalar cbrt(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void correctSurfaceValues(TrackCloudType &cloud, trackingData &td, const scalar T, const scalarField &Cs, scalar &rhos, scalar &mus, scalar &Prs, scalar &kappas)
Correct surface values due to emitted species.
static const zero Zero
Definition: zero.H:97
virtual void calculate(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar dt, const scalar Re, const scalar Pr, const scalar d, const scalar nu, const scalar T, const scalar Ts, const scalar pc, const scalar Tc, const scalarField &X, scalarField &dMassPC) const =0
Update model.
static const char nl
Definition: Ostream.H:260
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual scalar dh(const label idc, const label idl, const scalar p, const scalar T) const
Return the enthalpy per unit mass.
#define WarningInFunction
Report a warning using Foam::Warning.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
const scalarList W(::W(thermo))
dimensionedScalar pow4(const dimensionedScalar &ds)
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
Reacting parcel class with one/two-way coupling with the continuous phase.
scalar mass0_
Initial mass [kg].
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
Templated reacting parcel composition model class Consists of carrier species (via thermo package)...
Definition: ThermoCloud.H:61
virtual scalar Tvap(const scalarField &X) const
Return vapourisation temperature.
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
const basicSpecieMixture & carrier() const
Return the carrier components (wrapper function)
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
const speciesTable & species() const
Return the table of species.
const liquidMixtureProperties & liquids() const
Return the global (additional) liquids.
scalar T0
Definition: createFields.H:22