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-2023 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;
61  const CompositionModel<thermoCloudType>& composition =
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().Cpi(cid, td.pc(), Tsdash);
136  const scalar W = composition.carrier().WiValue(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  Y_(p.Y_)
188 {}
189 
190 
191 // * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * * //
192 
193 template<class ParcelType>
194 template<class TrackCloudType>
196 (
197  TrackCloudType& cloud,
198  trackingData& td
199 )
200 {
201  ParcelType::setCellValues(cloud, td);
202 
203  td.pc() = td.pInterp().interpolate
204  (
205  this->coordinates(),
206  this->currentTetIndices(td.mesh)
207  );
208 
209  if (td.pc() < cloud.constProps().pMin())
210  {
211  if (debug)
212  {
214  << "Limiting observed pressure in cell " << this->cell()
215  << " to " << cloud.constProps().pMin() << nl << endl;
216  }
217 
218  td.pc() = cloud.constProps().pMin();
219  }
220 }
221 
222 
223 template<class ParcelType>
224 template<class TrackCloudType>
226 (
227  TrackCloudType& cloud,
228  trackingData& td,
229  const scalar dt
230 )
231 {
232  scalar addedMass = 0.0;
233  scalar maxMassI = 0.0;
234  forAll(cloud.rhoTrans(), i)
235  {
236  scalar dm = cloud.rhoTrans(i)[this->cell()];
237  maxMassI = max(maxMassI, mag(dm));
238  addedMass += dm;
239  }
240 
241  if (maxMassI < rootVSmall)
242  {
243  return;
244  }
245 
246  const scalar massCell = this->massCell(td);
247 
248  td.rhoc() += addedMass/cloud.pMesh().cellVolumes()[this->cell()];
249 
250  const scalar massCellNew = massCell + addedMass;
251  td.Uc() = (td.Uc()*massCell + cloud.UTransRef()[this->cell()])/massCellNew;
252 
253  scalar CpEff = 0.0;
254  forAll(cloud.rhoTrans(), i)
255  {
256  scalar Y = cloud.rhoTrans(i)[this->cell()]/addedMass;
257  CpEff += Y*cloud.composition().carrier().Cpi(i, td.pc(), td.Tc());
258  }
259 
260  const scalar Cpc = td.CpInterp().psi()[this->cell()];
261  td.Cpc() = (massCell*Cpc + addedMass*CpEff)/massCellNew;
262 
263  td.Tc() += cloud.hsTransRef()[this->cell()]/(td.Cpc()*massCellNew);
264 
265  if (td.Tc() < cloud.constProps().TMin())
266  {
267  if (debug)
268  {
270  << "Limiting observed temperature in cell " << this->cell()
271  << " to " << cloud.constProps().TMin() << nl << endl;
272  }
273 
274  td.Tc() = cloud.constProps().TMin();
275  }
276 }
277 
278 
279 template<class ParcelType>
280 template<class TrackCloudType>
282 (
283  TrackCloudType& cloud,
284  trackingData& td,
285  const scalar T,
286  const scalarField& Cs,
287  scalar& rhos,
288  scalar& mus,
289  scalar& Prs,
290  scalar& kappas
291 )
292 {
293  // No correction if total concentration of emitted species is small
294  if (!cloud.heatTransfer().BirdCorrection() || (sum(Cs) < small))
295  {
296  return;
297  }
298 
299  const fluidMulticomponentThermo& carrierThermo =
300  cloud.composition().carrier();
301 
302  // Far field carrier molar fractions
303  scalarField Xinf(carrierThermo.species().size());
304 
305  forAll(Xinf, i)
306  {
307  Xinf[i] = carrierThermo.Y(i)[this->cell()]/carrierThermo.WiValue(i);
308  }
309  Xinf /= sum(Xinf);
310 
311  // Molar fraction of far field species at particle surface
312  const scalar Xsff = 1.0 - min(sum(Cs)*RR*this->T_/td.pc(), 1.0);
313 
314  // Surface carrier total molar concentration
315  const scalar CsTot = td.pc()/(RR*this->T_);
316 
317  // Surface carrier composition (molar fraction)
318  scalarField Xs(Xinf.size());
319 
320  // Surface carrier composition (mass fraction)
321  scalarField Ys(Xinf.size());
322 
323  forAll(Xs, i)
324  {
325  // Molar concentration of species at particle surface
326  const scalar Csi = Cs[i] + Xsff*Xinf[i]*CsTot;
327 
328  Xs[i] = (2.0*Csi + Xinf[i]*CsTot)/3.0;
329  Ys[i] = Xs[i]*carrierThermo.WiValue(i);
330  }
331  Xs /= sum(Xs);
332  Ys /= sum(Ys);
333 
334 
335  rhos = 0;
336  mus = 0;
337  kappas = 0;
338  scalar Cps = 0;
339  scalar sumYiSqrtW = 0;
340  scalar sumYiCbrtW = 0;
341 
342  forAll(Ys, i)
343  {
344  const scalar W = carrierThermo.WiValue(i);
345  const scalar sqrtW = sqrt(W);
346  const scalar cbrtW = cbrt(W);
347 
348  rhos += Xs[i]*W;
349  mus += Ys[i]*sqrtW*carrierThermo.mui(i, td.pc(), T);
350  kappas += Ys[i]*cbrtW*carrierThermo.kappai(i, td.pc(), T);
351  Cps += Xs[i]*carrierThermo.Cpi(i, td.pc(), T);
352 
353  sumYiSqrtW += Ys[i]*sqrtW;
354  sumYiCbrtW += Ys[i]*cbrtW;
355  }
356 
357  Cps = max(Cps, rootVSmall);
358 
359  rhos *= td.pc()/(RR*T);
360  rhos = max(rhos, rootVSmall);
361 
362  mus /= sumYiSqrtW;
363  mus = max(mus, rootVSmall);
364 
365  kappas /= sumYiCbrtW;
366  kappas = max(kappas, rootVSmall);
367 
368  Prs = Cps*mus/kappas;
369 }
370 
371 
372 template<class ParcelType>
373 template<class TrackCloudType>
375 (
376  TrackCloudType& cloud,
377  trackingData& td,
378  const scalar dt
379 )
380 {
381  typedef typename TrackCloudType::thermoCloudType thermoCloudType;
382  const CompositionModel<thermoCloudType>& composition =
383  cloud.composition();
384 
385 
386  // Define local properties at beginning of time step
387  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
388 
389  const scalar np0 = this->nParticle_;
390  const scalar d0 = this->d_;
391  const vector& U0 = this->U_;
392  const scalar T0 = this->T_;
393  const scalar mass0 = this->mass();
394 
395 
396  // Calc surface values
397  scalar Ts, rhos, mus, Prs, kappas;
398  this->calcSurfaceValues(cloud, td, T0, Ts, rhos, mus, Prs, kappas);
399  scalar Res = this->Re(rhos, U0, td.Uc(), d0, mus);
400 
401 
402  // Sources
403  // ~~~~~~~
404 
405  // Explicit momentum source for particle
406  vector Su = Zero;
407 
408  // Linearised momentum source coefficient
409  scalar Spu = 0.0;
410 
411  // Momentum transfer from the particle to the carrier phase
412  vector dUTrans = Zero;
413 
414  // Explicit enthalpy source for particle
415  scalar Sh = 0.0;
416 
417  // Linearised enthalpy source coefficient
418  scalar Sph = 0.0;
419 
420  // Sensible enthalpy transfer from the particle to the carrier phase
421  scalar dhsTrans = 0.0;
422 
423 
424  // 1. Compute models that contribute to mass transfer - U, T held constant
425  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
426 
427  // Phase change
428  // ~~~~~~~~~~~~
429 
430  // Mass transfer due to phase change
431  scalarField dMassPC(Y_.size(), 0.0);
432 
433  // Molar flux of species emitted from the particle (kmol/m^2/s)
434  scalar Ne = 0.0;
435 
436  // Sum Ni*Cpi*Wi of emission species
437  scalar NCpW = 0.0;
438 
439  // Surface concentrations of emitted species
440  scalarField Cs(composition.carrier().species().size(), 0.0);
441 
442  // Calc mass and enthalpy transfer due to phase change
443  calcPhaseChange
444  (
445  cloud,
446  td,
447  dt,
448  Res,
449  Prs,
450  Ts,
451  mus/rhos,
452  d0,
453  T0,
454  mass0,
455  0,
456  1.0,
457  Y_,
458  dMassPC,
459  Sh,
460  Ne,
461  NCpW,
462  Cs
463  );
464 
465 
466  // 2. Update the parcel properties due to change in mass
467  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
468 
469  scalarField dMass(dMassPC);
470  scalar mass1 = updateMassFraction(mass0, dMass, Y_);
471 
472  this->Cp_ = composition.Cp(0, Y_, td.pc(), T0);
473 
474  // Update particle density or diameter
475  if (cloud.constProps().constantVolume())
476  {
477  this->rho_ = mass1/this->volume();
478  }
479  else
480  {
481  this->d_ = cbrt(mass1/this->rho_*6.0/pi);
482  }
483 
484  // Remove the particle when mass falls below minimum threshold
485  if (np0*mass1 < cloud.constProps().minParcelMass())
486  {
487  td.keepParticle = false;
488 
489  if (cloud.solution().coupled())
490  {
491  scalar dm = np0*mass0;
492 
493  // Absorb parcel into carrier phase
494  forAll(Y_, i)
495  {
496  scalar dmi = dm*Y_[i];
497  label gid = composition.localToCarrierId(0, i);
498  scalar hs = composition.carrier().hsi(gid, td.pc(), T0);
499 
500  cloud.rhoTrans(gid)[this->cell()] += dmi;
501  cloud.hsTransRef()[this->cell()] += dmi*hs;
502  }
503  cloud.UTransRef()[this->cell()] += dm*U0;
504 
505  cloud.phaseChange().addToPhaseChangeMass(np0*mass1);
506  }
507 
508  return;
509  }
510 
511  // Correct surface values due to emitted species
512  correctSurfaceValues(cloud, td, Ts, Cs, rhos, mus, Prs, kappas);
513  Res = this->Re(rhos, U0, td.Uc(), this->d(), mus);
514 
515 
516  // 3. Compute heat- and momentum transfers
517  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
518 
519  // Heat transfer
520  // ~~~~~~~~~~~~~
521 
522  // Calculate new particle temperature
523  this->T_ =
524  this->calcHeatTransfer
525  (
526  cloud,
527  td,
528  dt,
529  Res,
530  Prs,
531  kappas,
532  NCpW,
533  Sh,
534  dhsTrans,
535  Sph
536  );
537 
538  this->Cp_ = composition.Cp(0, Y_, td.pc(), T0);
539 
540 
541  // Motion
542  // ~~~~~~
543 
544  // Calculate new particle velocity
545  this->U_ =
546  this->calcVelocity(cloud, td, dt, Res, mus, mass1, Su, dUTrans, Spu);
547 
548 
549  // 4. Accumulate carrier phase source terms
550  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
551 
552  if (cloud.solution().coupled())
553  {
554  // Transfer mass lost to carrier mass, momentum and enthalpy sources
555  forAll(dMass, i)
556  {
557  scalar dm = np0*dMass[i];
558  label gid = composition.localToCarrierId(0, i);
559  scalar hs = composition.carrier().hsi(gid, td.pc(), T0);
560 
561  cloud.rhoTrans(gid)[this->cell()] += dm;
562  cloud.UTransRef()[this->cell()] += dm*U0;
563  cloud.hsTransRef()[this->cell()] += dm*hs;
564  }
565 
566  // Update momentum transfer
567  cloud.UTransRef()[this->cell()] += np0*dUTrans;
568  cloud.UCoeffRef()[this->cell()] += np0*Spu;
569 
570  // Update sensible enthalpy transfer
571  cloud.hsTransRef()[this->cell()] += np0*dhsTrans;
572  cloud.hsCoeffRef()[this->cell()] += np0*Sph;
573 
574  // Update radiation fields
575  if (cloud.radiation())
576  {
577  const scalar ap = this->areaP();
578  const scalar T4 = pow4(T0);
579  cloud.radAreaP()[this->cell()] += dt*np0*ap;
580  cloud.radT4()[this->cell()] += dt*np0*T4;
581  cloud.radAreaPT4()[this->cell()] += dt*np0*ap*T4;
582  }
583  }
584 }
585 
586 
587 // * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * * * //
588 
589 #include "ReactingParcelIO.C"
590 
591 // ************************************************************************* //
scalar hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:11
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Templated reacting parcel composition model class Consists of carrier species (via thermo package),...
const liquidMixtureProperties & liquids() const
Return the global (additional) liquids.
const fluidMulticomponentThermo & carrier() const
Return the carrier components (wrapper function)
label localToCarrierId(const label phaseI, const label id, const bool allowNotFound=false) const
Return carrier id of component given local id.
virtual scalar Cp(const label phaseI, const scalarField &Y, const scalar p, const scalar T) const
Return specific heat capacity for the phase phaseI.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Templated phase change model class.
virtual scalar dh(const label idc, const label idl, const scalar p, const scalar T) const
Return the enthalpy per unit mass.
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.
void addToPhaseChangeMass(const scalar dMass)
Add to phase change mass.
virtual scalar Tvap(const scalarField &X) const
Return vapourisation temperature.
virtual scalar TMax(const scalar p, const scalarField &X) const
Return maximum/limiting temperature.
Reacting parcel class with one/two-way coupling with the continuous phase.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
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.
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.
ParcelType::trackingData trackingData
Use base tracking data.
scalar updateMassFraction(const scalar mass0, const scalarField &dMass, scalarField &Y) const
Update mass fraction.
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
ReactingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label facei)
Construct from mesh, coordinates and topology.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
Base-class for multi-component fluid thermodynamic properties.
virtual scalar mui(const label speciei, const scalar p, const scalar T) const =0
Dynamic viscosity [kg/m/s].
const PtrList< liquidProperties > & properties() const
Return the liquid properties.
scalarField X(const scalarField &Y) const
Returns the mole fractions corresponding to the given mass fractions.
virtual scalar Cpi(const label speciei, const scalar p, const scalar T) const =0
Heat capacity at constant pressure [J/kg/K].
virtual const speciesTable & species() const =0
The table of species.
virtual scalar kappai(const label speciei, const scalar p, const scalar T) const =0
Thermal conductivity [W/m/K].
virtual scalar WiValue(const label speciei) const =0
Molecular weight [kg/kmol].
virtual PtrList< volScalarField > & Y()=0
Access the mass-fraction fields.
virtual scalar hsi(const label speciei, const scalar p, const scalar T) const =0
Sensible enthalpy [J/kg].
#define WarningInFunction
Report a warning using Foam::Warning.
mathematical constants.
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
static const zero Zero
Definition: zero.H:97
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
scalarList W(const fluidMulticomponentThermo &thermo)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedScalar pow4(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar cbrt(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
volScalarField & p
PtrList< volScalarField > & Y
scalar T0
Definition: createFields.H:22