ReactingMultiphaseParcel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
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] = massSolid/massNew;
131 
132  scalar Ytotal = sum(YMix);
133 
134  YMix[GAS] /= Ytotal;
135  YMix[LIQ] /= Ytotal;
136  YMix[SLD] /= Ytotal;
137 
138  return massNew;
139 }
140 
141 
142 template<class ParcelType>
143 template<class TrackCloudType>
145 (
146  TrackCloudType& cloud,
147  const scalarField& dMassGas,
148  const scalarField& dMassLiquid,
149  const scalarField& dMassSolid,
150  const label idG,
151  const label idL,
152  const label idS,
153  const scalar p,
154  const scalar T
155 )
156 {
157  const auto& props = cloud.composition().phaseProps()[idG];
158  const auto& thermo = cloud.composition().thermo();
159 
160  scalarField dVolGas(dMassGas.size(), Zero);
161  forAll(dMassGas, i)
162  {
163  label cid = props.carrierIds()[i];
164  dVolGas[i] = -dMassGas[i]/thermo.carrier().rho(cid, p, T);
165  }
166 
167  scalarField dVolLiquid(dMassLiquid.size(), Zero);
168  forAll(dMassLiquid, i)
169  {
170  dVolLiquid[i] =
171  -dMassLiquid[i]/thermo.liquids().properties()[i].rho(p, T);
172  }
173 
174  scalarField dVolSolid(dMassSolid.size(), Zero);
175  forAll(dMassSolid, i)
176  {
177  dVolSolid[i] = -dMassSolid[i]/thermo.solids().properties()[i].rho();
178  }
179 
180  return (sum(dVolGas) + sum(dVolLiquid) + sum(dMassSolid));
181 }
182 
183 
184 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
185 
186 template<class ParcelType>
187 template<class TrackCloudType>
189 (
190  TrackCloudType& cloud,
191  trackingData& td
192 )
193 {
194  ParcelType::setCellValues(cloud, td);
195 }
196 
197 
198 template<class ParcelType>
199 template<class TrackCloudType>
201 (
202  TrackCloudType& cloud,
203  trackingData& td,
204  const scalar dt
205 )
206 {
207  // Re-use correction from reacting parcel
208  ParcelType::cellValueSourceCorrection(cloud, td, dt);
209 }
210 
211 
212 template<class ParcelType>
213 template<class TrackCloudType>
215 (
216  TrackCloudType& cloud,
217  trackingData& td,
218  const scalar dt
219 )
220 {
221  typedef typename TrackCloudType::reactingCloudType reactingCloudType;
223  cloud.composition();
224 
225  // Define local properties at beginning of timestep
226  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
227 
228  const scalar np0 = this->nParticle_;
229  const scalar d0 = this->d_;
230  const vector& U0 = this->U_;
231  const scalar T0 = this->T_;
232  const scalar mass0 = this->mass();
233  const scalar rho0 = this->rho_;
234 
235 
236  const scalar pc = td.pc();
237 
238  const scalarField& YMix = this->Y_;
239  const label idG = composition.idGas();
240  const label idL = composition.idLiquid();
241  const label idS = composition.idSolid();
242 
243 
244  // Calc surface values
245  scalar Ts, rhos, mus, Prs, kappas;
246  this->calcSurfaceValues(cloud, td, T0, Ts, rhos, mus, Prs, kappas);
247  scalar Res = this->Re(rhos, U0, td.Uc(), d0, mus);
248 
249 
250  // Sources
251  //~~~~~~~~
252 
253  // Explicit momentum source for particle
254  vector Su = Zero;
255 
256  // Linearised momentum source coefficient
257  scalar Spu = 0.0;
258 
259  // Momentum transfer from the particle to the carrier phase
260  vector dUTrans = Zero;
261 
262  // Explicit enthalpy source for particle
263  scalar Sh = 0.0;
264 
265  // Linearised enthalpy source coefficient
266  scalar Sph = 0.0;
267 
268  // Sensible enthalpy transfer from the particle to the carrier phase
269  scalar dhsTrans = 0.0;
270 
271 
272  // 1. Compute models that contribute to mass transfer - U, T held constant
273  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
274 
275  // Phase change in liquid phase
276  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
277 
278  // Mass transfer due to phase change
279  scalarField dMassPC(YLiquid_.size(), Zero);
280 
281  // Molar flux of species emitted from the particle (kmol/m^2/s)
282  scalar Ne = 0.0;
283 
284  // Sum Ni*Cpi*Wi of emission species
285  scalar NCpW = 0.0;
286 
287  // Surface concentrations of emitted species
288  scalarField Cs(composition.carrier().species().size(), Zero);
289 
290  // Calc mass and enthalpy transfer due to phase change
291  this->calcPhaseChange
292  (
293  cloud,
294  td,
295  dt,
296  Res,
297  Prs,
298  Ts,
299  mus/rhos,
300  d0,
301  T0,
302  mass0,
303  rho0,
304  idL,
305  YMix[LIQ],
306  YLiquid_,
307  YMix[SLD]*YSolid_,
308  dMassPC,
309  Sh,
310  Ne,
311  NCpW,
312  Cs
313  );
314 
315 
316  // Devolatilisation
317  // ~~~~~~~~~~~~~~~~
318 
319  // Mass transfer due to devolatilisation
320  scalarField dMassDV(YGas_.size(), Zero);
321 
322  // Calc mass and enthalpy transfer due to devolatilisation
323  calcDevolatilisation
324  (
325  cloud,
326  td,
327  dt,
328  this->age_,
329  Ts,
330  d0,
331  T0,
332  mass0,
333  this->mass0_,
334  YMix[GAS]*YGas_,
335  YMix[LIQ]*YLiquid_,
336  YMix[SLD]*YSolid_,
337  canCombust_,
338  dMassDV,
339  Sh,
340  Ne,
341  NCpW,
342  Cs
343  );
344 
345  // Surface reactions
346  // ~~~~~~~~~~~~~~~~~
347 
348  // Change in carrier phase composition due to surface reactions
349  scalarField dMassSRGas(YGas_.size(), Zero);
350  scalarField dMassSRLiquid(YLiquid_.size(), Zero);
351  scalarField dMassSRSolid(YSolid_.size(), Zero);
352  scalarField dMassSRCarrier(composition.carrier().species().size(), Zero);
353 
354  // Calc mass and enthalpy transfer due to surface reactions
355  calcSurfaceReactions
356  (
357  cloud,
358  td,
359  dt,
360  d0,
361  Res,
362  mus/rhos,
363  T0,
364  mass0,
365  canCombust_,
366  Ne,
367  YMix,
368  YGas_,
369  YLiquid_,
370  YSolid_,
371  dMassSRGas,
372  dMassSRLiquid,
373  dMassSRSolid,
374  dMassSRCarrier,
375  Sh,
376  dhsTrans
377  );
378 
379  // 2. Update the parcel properties due to change in mass
380  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
381 
382  scalarField dMassGas(dMassDV + dMassSRGas);
383  scalarField dMassLiquid(dMassPC + dMassSRLiquid);
384  scalarField dMassSolid(dMassSRSolid);
385 
386  scalar mass1 = mass0 - sum(dMassGas) - sum(dMassLiquid) - sum(dMassSolid);
387 
388  // Remove the particle when mass falls below minimum threshold
389  if (np0*mass1 < cloud.constProps().minParcelMass())
390  {
391  td.keepParticle = false;
392 
393  if (cloud.solution().coupled())
394  {
395  scalar dm = np0*mass0;
396 
397  // Absorb parcel into carrier phase
398  forAll(YGas_, i)
399  {
400  label gid = composition.localToCarrierId(GAS, i);
401  cloud.rhoTrans(gid)[this->cell()] += dm*YMix[GAS]*YGas_[i];
402  }
403  forAll(YLiquid_, i)
404  {
405  label gid = composition.localToCarrierId(LIQ, i);
406  cloud.rhoTrans(gid)[this->cell()] += dm*YMix[LIQ]*YLiquid_[i];
407  }
408 
409  // No mapping between solid components and carrier phase
410  /*
411  forAll(YSolid_, i)
412  {
413  label gid = composition.localToCarrierId(SLD, i);
414  cloud.rhoTrans(gid)[this->cell()] += dm*YMix[SLD]*YSolid_[i];
415  }
416  */
417 
418  cloud.UTrans()[this->cell()] += dm*U0;
419 
420  cloud.hsTrans()[this->cell()] +=
421  dm*HsEff(cloud, td, pc, T0, idG, idL, idS);
422 
423  cloud.phaseChange().addToPhaseChangeMass(np0*mass1);
424  }
425 
426  return;
427  }
428 
429  (void)updateMassFractions(mass0, dMassGas, dMassLiquid, dMassSolid);
430 
431  if
432  (
433  cloud.constProps().volUpdateType()
434  == constantProperties::volumeUpdateType::mUndefined
435  )
436  {
437  if (cloud.constProps().constantVolume())
438  {
439  this->rho_ = mass1/this->volume();
440  }
441  else
442  {
443  this->d_ = cbrt(mass1/this->rho_*6/pi);
444  }
445  }
446  else
447  {
448  switch (cloud.constProps().volUpdateType())
449  {
450  case constantProperties::volumeUpdateType::mConstRho :
451  {
452  this->d_ = cbrt(mass1/this->rho_*6/pi);
453  break;
454  }
455  case constantProperties::volumeUpdateType::mConstVol :
456  {
457  this->rho_ = mass1/this->volume();
458  break;
459  }
460  case constantProperties::volumeUpdateType::mUpdateRhoAndVol :
461  {
462  scalar deltaVol =
463  updatedDeltaVolume
464  (
465  cloud,
466  dMassGas,
467  dMassLiquid,
468  dMassSolid,
469  idG,
470  idL,
471  idS,
472  pc,
473  T0
474  );
475 
476  this->rho_ = mass1/(this->volume() + deltaVol);
477  this->d_ = cbrt(mass1/this->rho_*6/pi);
478  break;
479  }
480  }
481  }
482  // Correct surface values due to emitted species
483  this->correctSurfaceValues(cloud, td, Ts, Cs, rhos, mus, Prs, kappas);
484  Res = this->Re(rhos, U0, td.Uc(), this->d_, mus);
485 
486 
487  // 3. Compute heat- and momentum transfers
488  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
489 
490  // Heat transfer
491  // ~~~~~~~~~~~~~
492 
493  // Calculate new particle temperature
494  this->T_ =
495  this->calcHeatTransfer
496  (
497  cloud,
498  td,
499  dt,
500  Res,
501  Prs,
502  kappas,
503  NCpW,
504  Sh,
505  dhsTrans,
506  Sph
507  );
508 
509 
510  this->Cp_ = CpEff(cloud, td, pc, this->T_, idG, idL, idS);
511 
512 
513  // Motion
514  // ~~~~~~
515 
516  // Calculate new particle velocity
517  this->U_ =
518  this->calcVelocity(cloud, td, dt, Res, mus, mass1, Su, dUTrans, Spu);
519 
520 
521  // 4. Accumulate carrier phase source terms
522  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
523 
524  if (cloud.solution().coupled())
525  {
526  // Transfer mass lost to carrier mass, momentum and enthalpy sources
527  forAll(YGas_, i)
528  {
529  scalar dm = np0*dMassGas[i];
530  label gid = composition.localToCarrierId(GAS, i);
531  scalar hs = composition.carrier().Hs(gid, pc, T0);
532  cloud.rhoTrans(gid)[this->cell()] += dm;
533  cloud.UTrans()[this->cell()] += dm*U0;
534  cloud.hsTrans()[this->cell()] += dm*hs;
535  }
536  forAll(YLiquid_, i)
537  {
538  scalar dm = np0*dMassLiquid[i];
539  label gid = composition.localToCarrierId(LIQ, i);
540  scalar hs = composition.carrier().Hs(gid, pc, T0);
541  cloud.rhoTrans(gid)[this->cell()] += dm;
542  cloud.UTrans()[this->cell()] += dm*U0;
543  cloud.hsTrans()[this->cell()] += dm*hs;
544  }
545 
546  // No mapping between solid components and carrier phase
547  /*
548  forAll(YSolid_, i)
549  {
550  scalar dm = np0*dMassSolid[i];
551  label gid = composition.localToCarrierId(SLD, i);
552  scalar hs = composition.carrier().Hs(gid, pc, T0);
553  cloud.rhoTrans(gid)[this->cell()] += dm;
554  cloud.UTrans()[this->cell()] += dm*U0;
555  cloud.hsTrans()[this->cell()] += dm*hs;
556  }
557  */
558 
559  forAll(dMassSRCarrier, i)
560  {
561  scalar dm = np0*dMassSRCarrier[i];
562  scalar hs = composition.carrier().Hs(i, pc, T0);
563  cloud.rhoTrans(i)[this->cell()] += dm;
564  cloud.UTrans()[this->cell()] += dm*U0;
565  cloud.hsTrans()[this->cell()] += dm*hs;
566  }
567 
568  // Update momentum transfer
569  cloud.UTrans()[this->cell()] += np0*dUTrans;
570  cloud.UCoeff()[this->cell()] += np0*Spu;
571 
572  // Update sensible enthalpy transfer
573  cloud.hsTrans()[this->cell()] += np0*dhsTrans;
574  cloud.hsCoeff()[this->cell()] += np0*Sph;
575 
576  // Update radiation fields
577  if (cloud.radiation())
578  {
579  const scalar ap = this->areaP();
580  const scalar T4 = pow4(T0);
581  cloud.radAreaP()[this->cell()] += dt*np0*ap;
582  cloud.radT4()[this->cell()] += dt*np0*T4;
583  cloud.radAreaPT4()[this->cell()] += dt*np0*ap*T4;
584  }
585  }
586 }
587 
588 
589 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
590 
591 template<class ParcelType>
592 template<class TrackCloudType>
594 (
595  TrackCloudType& cloud,
596  trackingData& td,
597  const scalar dt,
598  const scalar age,
599  const scalar Ts,
600  const scalar d,
601  const scalar T,
602  const scalar mass,
603  const scalar mass0,
604  const scalarField& YGasEff,
605  const scalarField& YLiquidEff,
606  const scalarField& YSolidEff,
607  label& canCombust,
608  scalarField& dMassDV,
609  scalar& Sh,
610  scalar& N,
611  scalar& NCpW,
612  scalarField& Cs
613 ) const
614 {
615  // Check that model is active
616  if (!cloud.devolatilisation().active())
617  {
618  return;
619  }
620 
621  // Initialise demand-driven constants
622  (void)cloud.constProps().TDevol();
623  (void)cloud.constProps().LDevol();
624 
625  // Check that the parcel temperature is within necessary limits for
626  // devolatilisation to occur
627  if (T < cloud.constProps().TDevol() || canCombust == -1)
628  {
629  return;
630  }
631 
632  typedef typename TrackCloudType::reactingCloudType reactingCloudType;
634  cloud.composition();
635 
636 
637  // Total mass of volatiles evolved
638  cloud.devolatilisation().calculate
639  (
640  dt,
641  age,
642  mass0,
643  mass,
644  T,
645  YGasEff,
646  YLiquidEff,
647  YSolidEff,
648  canCombust,
649  dMassDV
650  );
651 
652  scalar dMassTot = sum(dMassDV);
653 
654  cloud.devolatilisation().addToDevolatilisationMass
655  (
656  this->nParticle_*dMassTot
657  );
658 
659  Sh -= dMassTot*cloud.constProps().LDevol()/dt;
660 
661  // Update molar emissions
662  if (cloud.heatTransfer().BirdCorrection())
663  {
664  // Molar average molecular weight of carrier mix
665  const scalar Wc = max(SMALL, td.rhoc()*RR*td.Tc()/td.pc());
666 
667  // Note: hardcoded gaseous diffusivities for now
668  // TODO: add to carrier thermo
669  const scalar beta = sqr(cbrt(15.0) + cbrt(15.0));
670 
671  forAll(dMassDV, i)
672  {
673  const label id = composition.localToCarrierId(GAS, i);
674  const scalar Cp = composition.carrier().Cp(id, td.pc(), Ts);
675  const scalar W = composition.carrier().W(id);
676  const scalar Ni = dMassDV[i]/(this->areaS(d)*dt*W);
677 
678  // Dab calc'd using API vapour mass diffusivity function
679  const scalar Dab =
680  3.6059e-3*(pow(1.8*Ts, 1.75))
681  *sqrt(1.0/W + 1.0/Wc)
682  /(td.pc()*beta);
683 
684  N += Ni;
685  NCpW += Ni*Cp*W;
686  Cs[id] += Ni*d/(2.0*Dab);
687  }
688  }
689 }
690 
691 
692 template<class ParcelType>
693 template<class TrackCloudType>
695 (
696  TrackCloudType& cloud,
697  trackingData& td,
698  const scalar dt,
699  const scalar d,
700  const scalar Re,
701  const scalar nu,
702  const scalar T,
703  const scalar mass,
704  const label canCombust,
705  const scalar N,
706  const scalarField& YMix,
707  const scalarField& YGas,
708  const scalarField& YLiquid,
709  const scalarField& YSolid,
710  scalarField& dMassSRGas,
711  scalarField& dMassSRLiquid,
712  scalarField& dMassSRSolid,
713  scalarField& dMassSRCarrier,
714  scalar& Sh,
715  scalar& dhsTrans
716 ) const
717 {
718  // Check that model is active
719  if (!cloud.surfaceReaction().active())
720  {
721  return;
722  }
723 
724  // Initialise demand-driven constants
725  (void)cloud.constProps().hRetentionCoeff();
726  (void)cloud.constProps().TMax();
727 
728  // Check that model is active
729  if (canCombust != 1)
730  {
731  return;
732  }
733 
734 
735  // Update surface reactions
736  const scalar hReaction = cloud.surfaceReaction().calculate
737  (
738  dt,
739  Re,
740  nu,
741  this->cell(),
742  d,
743  T,
744  td.Tc(),
745  td.pc(),
746  td.rhoc(),
747  mass,
748  YGas,
749  YLiquid,
750  YSolid,
751  YMix,
752  N,
753  dMassSRGas,
754  dMassSRLiquid,
755  dMassSRSolid,
756  dMassSRCarrier
757  );
758 
759  cloud.surfaceReaction().addToSurfaceReactionMass
760  (
761  this->nParticle_
762  *(sum(dMassSRGas) + sum(dMassSRLiquid) + sum(dMassSRSolid))
763  );
764 
765  const scalar xsi = min(T/cloud.constProps().TMax(), 1.0);
766  const scalar coeff =
767  (1.0 - xsi*xsi)*cloud.constProps().hRetentionCoeff();
768 
769  Sh += coeff*hReaction/dt;
770 
771  dhsTrans += (1.0 - coeff)*hReaction;
772 }
773 
774 
775 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
776 
777 template<class ParcelType>
779 (
781 )
782 :
783  ParcelType(p),
784  YGas_(p.YGas_),
785  YLiquid_(p.YLiquid_),
786  YSolid_(p.YSolid_),
787  canCombust_(p.canCombust_)
788 {}
789 
790 
791 template<class ParcelType>
793 (
794  const ReactingMultiphaseParcel<ParcelType>& p,
795  const polyMesh& mesh
796 )
797 :
798  ParcelType(p, mesh),
799  YGas_(p.YGas_),
800  YLiquid_(p.YLiquid_),
801  YSolid_(p.YSolid_),
802  canCombust_(p.canCombust_)
803 {}
804 
805 
806 // * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * * * //
807 
809 
810 // ************************************************************************* //
Foam::constant::thermodynamic::RR
const scalar RR
Universal gas constant: default in [J/(kmol K)].
Definition: thermodynamicConstants.C:46
mathematicalConstants.H
Foam::ReactingMultiphaseParcel::updatedDeltaVolume
scalar updatedDeltaVolume(TrackCloudType &cloud, const scalarField &dMassGas, const scalarField &dMassLiquid, const scalarField &dMassSolid, const label idG, const label idL, const label idS, const scalar p, const scalar T)
Return change of volume due to mass exchange.
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::ReactingMultiphaseParcel::cellValueSourceCorrection
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
Definition: ReactingMultiphaseParcel.C:201
Foam::ReactingMultiphaseParcel::calc
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
Definition: ReactingMultiphaseParcel.C:215
Foam::ReactingMultiphaseParcel::trackingData
ParcelType::trackingData trackingData
Use base tracking data.
Definition: ReactingMultiphaseParcel.H:129
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::ReactingMultiphaseParcel::SLD
static const label SLD
Definition: ReactingMultiphaseParcel.H:79
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::fieldTypes::volume
const wordList volume
Standard volume field types (scalar, vector, tensor, etc)
composition
basicSpecieMixture & composition
Definition: createFieldRefs.H:6
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Su
zeroField Su
Definition: alphaSuSp.H:1
Foam::CompositionModel
Templated reacting parcel composition model class Consists of carrier species (via thermo package),...
Definition: ReactingCloud.H:58
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::Field< scalar >
ReactingMultiphaseParcelIO.C
Foam::ReactingMultiphaseParcel
Multiphase variant of the reacting parcel class with one/two-way coupling with the continuous phase.
Definition: ReactingMultiphaseParcel.H:55
ReactingMultiphaseParcel.H
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::ReactingMultiphaseParcel::calcSurfaceReactions
void calcSurfaceReactions(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar d, const scalar Re, const scalar nu, 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.
Definition: ReactingMultiphaseParcel.C:695
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::ReactingMultiphaseParcel::LIQ
static const label LIQ
Definition: ReactingMultiphaseParcel.H:78
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:57
rho0
scalar rho0
Definition: readInitialConditions.H:89
Foam::ReactingMultiphaseParcel::GAS
static const label GAS
Definition: ReactingMultiphaseParcel.H:77
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical
Mathematical constants.
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::Vector< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::ReactingMultiphaseParcel::setCellValues
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
Definition: ReactingMultiphaseParcel.C:189
Foam::Re
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Cs
const scalarField & Cs
Definition: solveBulkSurfactant.H:10
Foam::ReactingMultiphaseParcel::calcDevolatilisation
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.
Definition: ReactingMultiphaseParcel.C:594
Foam::ReactingMultiphaseParcel::ReactingMultiphaseParcel
ReactingMultiphaseParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, position and topology.
Definition: ReactingMultiphaseParcelI.H:71
N
const Vector< label > N(dict.get< Vector< label >>("N"))
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
T0
scalar T0
Definition: createFields.H:22
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54