ReactingMultiphaseParcel.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-2021 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 
27 #include "CompositionModel.H"
28 #include "NoDevolatilisation.H"
29 #include "NoSurfaceReaction.H"
30 #include "mathematicalConstants.H"
31 
32 using namespace Foam::constant::mathematical;
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 template<class ParcelType>
38 
39 template<class ParcelType>
41 
42 template<class ParcelType>
44 
45 
46 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
47 
48 template<class ParcelType>
49 template<class TrackCloudType>
51 (
52  TrackCloudType& cloud,
53  trackingData& td,
54  const scalar p,
55  const scalar T,
56  const label idG,
57  const label idL,
58  const label idS
59 ) const
60 {
61  return
62  this->Y_[GAS]*cloud.composition().Cp(idG, YGas_, p, T)
63  + this->Y_[LIQ]*cloud.composition().Cp(idL, YLiquid_, p, T)
64  + this->Y_[SLD]*cloud.composition().Cp(idS, YSolid_, p, T);
65 }
66 
67 
68 template<class ParcelType>
69 template<class TrackCloudType>
71 (
72  TrackCloudType& cloud,
73  trackingData& td,
74  const scalar p,
75  const scalar T,
76  const label idG,
77  const label idL,
78  const label idS
79 ) const
80 {
81  return
82  this->Y_[GAS]*cloud.composition().Hs(idG, YGas_, p, T)
83  + this->Y_[LIQ]*cloud.composition().Hs(idL, YLiquid_, p, T)
84  + this->Y_[SLD]*cloud.composition().Hs(idS, YSolid_, p, T);
85 }
86 
87 
88 template<class ParcelType>
89 template<class TrackCloudType>
91 (
92  TrackCloudType& cloud,
93  trackingData& td,
94  const scalar p,
95  const scalar T,
96  const label idG,
97  const label idL,
98  const label idS
99 ) const
100 {
101  return
102  this->Y_[GAS]*cloud.composition().L(idG, YGas_, p, T)
103  + this->Y_[LIQ]*cloud.composition().L(idL, YLiquid_, p, T)
104  + this->Y_[SLD]*cloud.composition().L(idS, YSolid_, p, T);
105 }
106 
107 
108 template<class ParcelType>
110 (
111  const scalar mass0,
112  const scalarField& dMassGas,
113  const scalarField& dMassLiquid,
114  const scalarField& dMassSolid
115 )
116 {
117  scalarField& YMix = this->Y_;
118 
119  scalar massGas =
120  this->updateMassFraction(mass0*YMix[GAS], dMassGas, YGas_);
121  scalar massLiquid =
122  this->updateMassFraction(mass0*YMix[LIQ], dMassLiquid, YLiquid_);
123  scalar massSolid =
124  this->updateMassFraction(mass0*YMix[SLD], dMassSolid, YSolid_);
125 
126  scalar massNew = max(massGas + massLiquid + massSolid, rootVSmall);
127 
128  YMix[GAS] = massGas/massNew;
129  YMix[LIQ] = massLiquid/massNew;
130  YMix[SLD] = 1.0 - YMix[GAS] - YMix[LIQ];
131 
132  return massNew;
133 }
134 
135 
136 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
137 
138 template<class ParcelType>
139 template<class TrackCloudType>
141 (
142  TrackCloudType& cloud,
143  trackingData& td
144 )
145 {
146  ParcelType::setCellValues(cloud, td);
147 }
148 
149 
150 template<class ParcelType>
151 template<class TrackCloudType>
153 (
154  TrackCloudType& cloud,
155  trackingData& td,
156  const scalar dt
157 )
158 {
159  // Re-use correction from reacting parcel
160  ParcelType::cellValueSourceCorrection(cloud, td, dt);
161 }
162 
163 
164 template<class ParcelType>
165 template<class TrackCloudType>
167 (
168  TrackCloudType& cloud,
169  trackingData& td,
170  const scalar dt
171 )
172 {
173  typedef typename TrackCloudType::thermoCloudType thermoCloudType;
175  cloud.composition();
176 
177 
178  // Define local properties at beginning of timestep
179  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
180 
181  const scalar np0 = this->nParticle_;
182  const scalar d0 = this->d_;
183  const vector& U0 = this->U_;
184  const scalar T0 = this->T_;
185  const scalar mass0 = this->mass();
186 
187  const scalar pc = td.pc();
188 
189  const scalarField& YMix = this->Y_;
190  const label idG = composition.idGas();
191  const label idL = composition.idLiquid();
192  const label idS = composition.idSolid();
193 
194 
195  // Calc surface values
196  scalar Ts, rhos, mus, Prs, kappas;
197  this->calcSurfaceValues(cloud, td, T0, Ts, rhos, mus, Prs, kappas);
198  scalar Res = this->Re(rhos, U0, td.Uc(), d0, mus);
199 
200 
201  // Sources
202  //~~~~~~~~
203 
204  // Explicit momentum source for particle
205  vector Su = Zero;
206 
207  // Linearised momentum source coefficient
208  scalar Spu = 0.0;
209 
210  // Momentum transfer from the particle to the carrier phase
211  vector dUTrans = Zero;
212 
213  // Explicit enthalpy source for particle
214  scalar Sh = 0.0;
215 
216  // Linearised enthalpy source coefficient
217  scalar Sph = 0.0;
218 
219  // Sensible enthalpy transfer from the particle to the carrier phase
220  scalar dhsTrans = 0.0;
221 
222 
223  // 1. Compute models that contribute to mass transfer - U, T held constant
224  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
225 
226  // Phase change in liquid phase
227  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
228 
229  // Mass transfer due to phase change
230  scalarField dMassPC(YLiquid_.size(), 0.0);
231 
232  // Molar flux of species emitted from the particle (kmol/m^2/s)
233  scalar Ne = 0.0;
234 
235  // Sum Ni*Cpi*Wi of emission species
236  scalar NCpW = 0.0;
237 
238  // Surface concentrations of emitted species
239  scalarField Cs(composition.carrier().species().size(), 0.0);
240 
241  // Calc mass and enthalpy transfer due to phase change
242  this->calcPhaseChange
243  (
244  cloud,
245  td,
246  dt,
247  Res,
248  Prs,
249  Ts,
250  mus/rhos,
251  d0,
252  T0,
253  mass0,
254  idL,
255  YMix[LIQ],
256  YLiquid_,
257  dMassPC,
258  Sh,
259  Ne,
260  NCpW,
261  Cs
262  );
263 
264 
265  // Devolatilisation
266  // ~~~~~~~~~~~~~~~~
267 
268  // Mass transfer due to devolatilisation
269  scalarField dMassDV(YGas_.size(), 0.0);
270 
271  // Calc mass and enthalpy transfer due to devolatilisation
272  calcDevolatilisation
273  (
274  cloud,
275  td,
276  dt,
277  this->age_,
278  Ts,
279  d0,
280  T0,
281  mass0,
282  this->mass0_,
283  YMix[GAS]*YGas_,
284  YMix[LIQ]*YLiquid_,
285  YMix[SLD]*YSolid_,
286  canCombust_,
287  dMassDV,
288  Sh,
289  Ne,
290  NCpW,
291  Cs
292  );
293 
294 
295  // Surface reactions
296  // ~~~~~~~~~~~~~~~~~
297 
298  // Change in carrier phase composition due to surface reactions
299  scalarField dMassSRGas(YGas_.size(), 0.0);
300  scalarField dMassSRLiquid(YLiquid_.size(), 0.0);
301  scalarField dMassSRSolid(YSolid_.size(), 0.0);
302  scalarField dMassSRCarrier(composition.carrier().species().size(), 0.0);
303 
304  // Calc mass and enthalpy transfer due to surface reactions
305  calcSurfaceReactions
306  (
307  cloud,
308  td,
309  dt,
310  d0,
311  T0,
312  mass0,
313  canCombust_,
314  Ne,
315  YMix,
316  YGas_,
317  YLiquid_,
318  YSolid_,
319  dMassSRGas,
320  dMassSRLiquid,
321  dMassSRSolid,
322  dMassSRCarrier,
323  Sh,
324  dhsTrans
325  );
326 
327 
328  // 2. Update the parcel properties due to change in mass
329  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
330 
331  scalarField dMassGas(dMassDV + dMassSRGas);
332  scalarField dMassLiquid(dMassPC + dMassSRLiquid);
333  scalarField dMassSolid(dMassSRSolid);
334  scalar mass1 =
335  updateMassFractions(mass0, dMassGas, dMassLiquid, dMassSolid);
336 
337  this->Cp_ = CpEff(cloud, td, pc, T0, idG, idL, idS);
338 
339  // Update particle density or diameter
340  if (cloud.constProps().constantVolume())
341  {
342  this->rho_ = mass1/this->volume();
343  }
344  else
345  {
346  this->d_ = cbrt(mass1/this->rho_*6.0/pi);
347  }
348 
349  // Remove the particle when mass falls below minimum threshold
350  if (np0*mass1 < cloud.constProps().minParcelMass())
351  {
352  td.keepParticle = false;
353 
354  if (cloud.solution().coupled())
355  {
356  scalar dm = np0*mass0;
357 
358  // Absorb parcel into carrier phase
359  forAll(YGas_, i)
360  {
361  label gid = composition.localToCarrierId(GAS, i);
362  cloud.rhoTrans(gid)[this->cell()] += dm*YMix[GAS]*YGas_[i];
363  }
364  forAll(YLiquid_, i)
365  {
366  label gid = composition.localToCarrierId(LIQ, i);
367  cloud.rhoTrans(gid)[this->cell()] += dm*YMix[LIQ]*YLiquid_[i];
368  }
369 
370  // No mapping between solid components and carrier phase
371  /*
372  forAll(YSolid_, i)
373  {
374  label gid = composition.localToCarrierId(SLD, i);
375  cloud.rhoTrans(gid)[this->cell()] += dm*YMix[SLD]*YSolid_[i];
376  }
377  */
378 
379  cloud.UTransRef()[this->cell()] += dm*U0;
380 
381  cloud.hsTransRef()[this->cell()] +=
382  dm*HsEff(cloud, td, pc, T0, idG, idL, idS);
383 
384  cloud.phaseChange().addToPhaseChangeMass(np0*mass1);
385  }
386 
387  return;
388  }
389 
390  // Correct surface values due to emitted species
391  this->correctSurfaceValues(cloud, td, Ts, Cs, rhos, mus, Prs, kappas);
392  Res = this->Re(rhos, U0, td.Uc(), this->d_, mus);
393 
394 
395  // 3. Compute heat- and momentum transfers
396  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
397 
398  // Heat transfer
399  // ~~~~~~~~~~~~~
400 
401  // Calculate new particle temperature
402  this->T_ =
403  this->calcHeatTransfer
404  (
405  cloud,
406  td,
407  dt,
408  Res,
409  Prs,
410  kappas,
411  NCpW,
412  Sh,
413  dhsTrans,
414  Sph
415  );
416 
417 
418  this->Cp_ = CpEff(cloud, td, pc, this->T_, idG, idL, idS);
419 
420 
421  // Motion
422  // ~~~~~~
423 
424  // Calculate new particle velocity
425  this->U_ =
426  this->calcVelocity(cloud, td, dt, Res, mus, mass1, Su, dUTrans, Spu);
427 
428 
429  // 4. Accumulate carrier phase source terms
430  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
431 
432  if (cloud.solution().coupled())
433  {
434  // Transfer mass lost to carrier mass, momentum and enthalpy sources
435  forAll(YGas_, i)
436  {
437  scalar dm = np0*dMassGas[i];
438  label gid = composition.localToCarrierId(GAS, i);
439  scalar hs = composition.carrier().Hs(gid, pc, T0);
440  cloud.rhoTrans(gid)[this->cell()] += dm;
441  cloud.UTransRef()[this->cell()] += dm*U0;
442  cloud.hsTransRef()[this->cell()] += dm*hs;
443  }
444  forAll(YLiquid_, i)
445  {
446  scalar dm = np0*dMassLiquid[i];
447  label gid = composition.localToCarrierId(LIQ, i);
448  scalar hs = composition.carrier().Hs(gid, pc, T0);
449  cloud.rhoTrans(gid)[this->cell()] += dm;
450  cloud.UTransRef()[this->cell()] += dm*U0;
451  cloud.hsTransRef()[this->cell()] += dm*hs;
452  }
453 
454  // No mapping between solid components and carrier phase
455  /*
456  forAll(YSolid_, i)
457  {
458  scalar dm = np0*dMassSolid[i];
459  label gid = composition.localToCarrierId(SLD, i);
460  scalar hs = composition.carrier().Hs(gid, pc, T0);
461  cloud.rhoTrans(gid)[this->cell()] += dm;
462  cloud.UTransRef()[this->cell()] += dm*U0;
463  cloud.hsTransRef()[this->cell()] += dm*hs;
464  }
465  */
466 
467  forAll(dMassSRCarrier, i)
468  {
469  scalar dm = np0*dMassSRCarrier[i];
470  scalar hs = composition.carrier().Hs(i, pc, T0);
471  cloud.rhoTrans(i)[this->cell()] += dm;
472  cloud.UTransRef()[this->cell()] += dm*U0;
473  cloud.hsTransRef()[this->cell()] += dm*hs;
474  }
475 
476  // Update momentum transfer
477  cloud.UTransRef()[this->cell()] += np0*dUTrans;
478  cloud.UCoeffRef()[this->cell()] += np0*Spu;
479 
480  // Update sensible enthalpy transfer
481  cloud.hsTransRef()[this->cell()] += np0*dhsTrans;
482  cloud.hsCoeffRef()[this->cell()] += np0*Sph;
483 
484  // Update radiation fields
485  if (cloud.radiation())
486  {
487  const scalar ap = this->areaP();
488  const scalar T4 = pow4(T0);
489  cloud.radAreaP()[this->cell()] += dt*np0*ap;
490  cloud.radT4()[this->cell()] += dt*np0*T4;
491  cloud.radAreaPT4()[this->cell()] += dt*np0*ap*T4;
492  }
493  }
494 }
495 
496 
497 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
498 
499 template<class ParcelType>
500 template<class TrackCloudType>
502 (
503  TrackCloudType& cloud,
504  trackingData& td,
505  const scalar dt,
506  const scalar age,
507  const scalar Ts,
508  const scalar d,
509  const scalar T,
510  const scalar mass,
511  const scalar mass0,
512  const scalarField& YGasEff,
513  const scalarField& YLiquidEff,
514  const scalarField& YSolidEff,
515  label& canCombust,
516  scalarField& dMassDV,
517  scalar& Sh,
518  scalar& N,
519  scalar& NCpW,
520  scalarField& Cs
521 ) const
522 {
523  // Check that model is active
524  if
525  (
526  isType
527  <
529  <
530  typename TrackCloudType::reactingMultiphaseCloudType
531  >
532  >
533  (
534  cloud.devolatilisation()
535  )
536  )
537  {
538  if (canCombust != -1)
539  {
540  canCombust = 1;
541  }
542  return;
543  }
544 
545  // Initialise demand-driven constants
546  (void)cloud.constProps().TDevol();
547  (void)cloud.constProps().LDevol();
548 
549  // Check that the parcel temperature is within necessary limits for
550  // devolatilisation to occur
551  if (T < cloud.constProps().TDevol() || canCombust == -1)
552  {
553  return;
554  }
555 
556  typedef typename TrackCloudType::thermoCloudType thermoCloudType;
558  cloud.composition();
559 
560 
561  // Total mass of volatiles evolved
562  cloud.devolatilisation().calculate
563  (
564  dt,
565  age,
566  mass0,
567  mass,
568  T,
569  YGasEff,
570  YLiquidEff,
571  YSolidEff,
572  canCombust,
573  dMassDV
574  );
575 
576  scalar dMassTot = sum(dMassDV);
577 
578  cloud.devolatilisation().addToDevolatilisationMass
579  (
580  this->nParticle_*dMassTot
581  );
582 
583  Sh -= dMassTot*cloud.constProps().LDevol()/dt;
584 
585  // Update molar emissions
586  if (cloud.heatTransfer().BirdCorrection())
587  {
588  // Molar average molecular weight of carrier mix
589  const scalar Wc = max(small, td.rhoc()*RR*td.Tc()/td.pc());
590 
591  // Note: hardcoded gaseous diffusivities for now
592  // TODO: add to carrier thermo
593  const scalar beta = sqr(cbrt(15.0) + cbrt(15.0));
594 
595  forAll(dMassDV, i)
596  {
597  const label id = composition.localToCarrierId(GAS, i);
598  const scalar Cp = composition.carrier().Cp(id, td.pc(), Ts);
599  const scalar W = composition.carrier().Wi(id);
600  const scalar Ni = dMassDV[i]/(this->areaS(d)*dt*W);
601 
602  // Dab calc'd using API vapour mass diffusivity function
603  const scalar Dab =
604  3.6059e-3*(pow(1.8*Ts, 1.75))
605  *sqrt(1.0/W + 1.0/Wc)
606  /(td.pc()*beta);
607 
608  N += Ni;
609  NCpW += Ni*Cp*W;
610  Cs[id] += Ni*d/(2.0*Dab);
611  }
612  }
613 }
614 
615 
616 template<class ParcelType>
617 template<class TrackCloudType>
619 (
620  TrackCloudType& cloud,
621  trackingData& td,
622  const scalar dt,
623  const scalar d,
624  const scalar T,
625  const scalar mass,
626  const label canCombust,
627  const scalar N,
628  const scalarField& YMix,
629  const scalarField& YGas,
630  const scalarField& YLiquid,
631  const scalarField& YSolid,
632  scalarField& dMassSRGas,
633  scalarField& dMassSRLiquid,
634  scalarField& dMassSRSolid,
635  scalarField& dMassSRCarrier,
636  scalar& Sh,
637  scalar& dhsTrans
638 ) const
639 {
640  // Check that model is active
641  if
642  (
643  isType
644  <
646  <
647  typename TrackCloudType::reactingMultiphaseCloudType
648  >
649  >
650  (
651  cloud.surfaceReaction()
652  )
653  )
654  {
655  return;
656  }
657 
658  // Initialise demand-driven constants
659  (void)cloud.constProps().hRetentionCoeff();
660  (void)cloud.constProps().TMax();
661 
662  // Check that model is active
663  if (canCombust != 1)
664  {
665  return;
666  }
667 
668 
669  // Update surface reactions
670  const scalar hReaction = cloud.surfaceReaction().calculate
671  (
672  dt,
673  this->cell(),
674  d,
675  T,
676  td.Tc(),
677  td.pc(),
678  td.rhoc(),
679  mass,
680  YGas,
681  YLiquid,
682  YSolid,
683  YMix,
684  N,
685  dMassSRGas,
686  dMassSRLiquid,
687  dMassSRSolid,
688  dMassSRCarrier
689  );
690 
691  cloud.surfaceReaction().addToSurfaceReactionMass
692  (
693  this->nParticle_
694  *(sum(dMassSRGas) + sum(dMassSRLiquid) + sum(dMassSRSolid))
695  );
696 
697  const scalar xsi = min(T/cloud.constProps().TMax(), 1.0);
698  const scalar coeff =
699  (1.0 - xsi*xsi)*cloud.constProps().hRetentionCoeff();
700 
701  Sh += coeff*hReaction/dt;
702 
703  dhsTrans += (1.0 - coeff)*hReaction;
704 }
705 
706 
707 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
708 
709 template<class ParcelType>
711 (
713 )
714 :
715  ParcelType(p),
716  YGas_(p.YGas_),
717  YLiquid_(p.YLiquid_),
718  YSolid_(p.YSolid_),
719  canCombust_(p.canCombust_)
720 {}
721 
722 
723 template<class ParcelType>
725 (
727  const polyMesh& mesh
728 )
729 :
730  ParcelType(p, mesh),
731  YGas_(p.YGas_),
732  YLiquid_(p.YLiquid_),
733  YSolid_(p.YSolid_),
734  canCombust_(p.canCombust_)
735 {}
736 
737 
738 // * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * * * //
739 
741 
742 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
basicSpecieMixture & composition
label canCombust_
Flag to identify if the particle can devolatilise and combust.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual label idSolid() const =0
Solid id.
void calcDevolatilisation(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar age, const scalar Ts, const scalar d, const scalar T, const scalar mass, const scalar mass0, const scalarField &YGasEff, const scalarField &YLiquidEff, const scalarField &YSolidEff, label &canCombust, scalarField &dMassDV, scalar &Sh, scalar &N, scalar &NCpW, scalarField &Cs) const
Calculate Devolatilisation.
dimensionedScalar sqrt(const dimensionedScalar &ds)
label localToCarrierId(const label phaseI, const label id, const bool allowNotFound=false) const
Return carrier id of component given local id.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Dummy surface reaction model for &#39;none&#39;.
ReactingMultiphaseParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, position and topology.
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
ParcelType::trackingData trackingData
Use base tracking data.
scalarField YGas_
Mass fractions of gases [].
mathematical constants.
bool isType(const Type &t)
Check the typeid.
Definition: typeInfo.H:126
dimensionedScalar cbrt(const dimensionedScalar &ds)
static const zero Zero
Definition: zero.H:97
virtual label idLiquid() const =0
Liquid id.
void calcSurfaceReactions(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar d, const scalar T, const scalar mass, const label canCombust, const scalar N, const scalarField &YMix, const scalarField &YGas, const scalarField &YLiquid, const scalarField &YSolid, scalarField &dMassSRGas, scalarField &dMassSRLiquid, scalarField &dMassSRSolid, scalarField &dMassSRCarrier, scalar &Sh, scalar &dhsTrans) const
Calculate surface reactions.
scalarField YSolid_
Mass fractions of solids [].
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
scalarField YLiquid_
Mass fractions of liquids [].
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
const scalarList W(::W(thermo))
dimensionedScalar pow4(const dimensionedScalar &ds)
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
Multiphase variant of the reacting parcel class with one/two-way coupling with the continuous phase...
label N
Definition: createFields.H:9
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
virtual label idGas() const =0
Gas id.
Templated reacting parcel composition model class Consists of carrier species (via thermo package)...
Definition: ThermoCloud.H:61
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
const basicSpecieMixture & carrier() const
Return the carrier components (wrapper function)
scalar T0
Definition: createFields.H:22
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
Dummy devolatilisation model for &#39;none&#39;.