MaxwellStefan.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) 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 
26 #include "MaxwellStefan.H"
27 #include "fvcDiv.H"
28 #include "fvcLaplacian.H"
29 #include "fvcSnGrad.H"
30 #include "fvmSup.H"
31 #include "surfaceInterpolate.H"
32 #include "Function2Evaluate.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 template<class BasicThermophysicalTransportModel>
43 (
44  const word& type,
45  const momentumTransportModel& momentumTransport,
46  const thermoModel& thermo
47 )
48 :
49  BasicThermophysicalTransportModel
50  (
51  type,
52  momentumTransport,
53  thermo
54  ),
55 
56  DFuncs_(this->thermo().composition().species().size()),
57 
58  DTFuncs_
59  (
60  this->coeffDict_.found("DT")
61  ? this->thermo().composition().species().size()
62  : 0
63  ),
64 
65  Dii_(this->thermo().composition().species().size()),
66  jexp_(this->thermo().composition().species().size()),
67 
68  W(this->thermo().composition().species().size()),
69 
70  YPtrs(W.size()),
71  DijPtrs(W.size()),
72 
73  Y(W.size()),
74  X(W.size()),
75  DD(W.size()),
76  A(W.size() - 1),
77  B(A.m()),
78  invA(A.m()),
79  D(W.size())
80 {
82 
83  // Set the molecular weights of the species
84  forAll(W, i)
85  {
86  W[i] = composition.Wi(i);
87  }
88 }
89 
90 
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 
93 template<class BasicThermophysicalTransportModel>
95 {
96  if
97  (
99  )
100  {
102  const speciesTable& species = composition.species();
103 
104  const dictionary& Ddict = this->coeffDict_.subDict("D");
105 
106  // Read the array of specie binary mass diffusion coefficient functions
107  forAll(species, i)
108  {
109  DFuncs_[i].setSize(species.size());
110 
111  forAll(species, j)
112  {
113  if (j >= i)
114  {
115  const word nameij(species[i] + '-' + species[j]);
116  const word nameji(species[j] + '-' + species[i]);
117 
118  word Dname;
119 
120  if (Ddict.found(nameij) && Ddict.found(nameji))
121  {
122  if (i != j)
123  {
125  << "Binary mass diffusion coefficients "
126  "for both " << nameij << " and " << nameji
127  << " provided, using " << nameij << endl;
128  }
129 
130  Dname = nameij;
131  }
132  else if (Ddict.found(nameij))
133  {
134  Dname = nameij;
135  }
136  else if (Ddict.found(nameji))
137  {
138  Dname = nameji;
139  }
140  else
141  {
143  << "Binary mass diffusion coefficients for pair "
144  << nameij << " or " << nameji << " not provided"
145  << exit(FatalIOError);
146  }
147 
148  DFuncs_[i].set
149  (
150  j,
151  Function2<scalar>::New(Dname, Ddict).ptr()
152  );
153  }
154  }
155  }
156 
157  // Optionally read the List of specie Soret thermal diffusion
158  // coefficient functions
159  if (this->coeffDict_.found("DT"))
160  {
161  const dictionary& DTdict = this->coeffDict_.subDict("DT");
162 
163  forAll(species, i)
164  {
165  DTFuncs_.set
166  (
167  i,
168  Function2<scalar>::New(species[i], DTdict).ptr()
169  );
170  }
171  }
172 
173  return true;
174  }
175  else
176  {
177  return false;
178  }
179 }
180 
181 
182 template<class BasicThermophysicalTransportModel>
184 (
185  const volScalarField& Yi
186 ) const
187 {
189 
190  return volScalarField::New
191  (
192  "DEff",
193  this->momentumTransport().rho()*Dii_[composition.index(Yi)]
194  );
195 }
196 
197 
198 template<class BasicThermophysicalTransportModel>
200 (
201  const volScalarField& Yi,
202  const label patchi
203 ) const
204 {
206 
207  return
208  this->momentumTransport().rho().boundaryField()[patchi]
209  *Dii_[composition.index(Yi)].boundaryField()[patchi];
210 }
211 
212 
213 template<class BasicThermophysicalTransportModel>
216 {
218  (
220  (
222  (
223  "q",
224  this->momentumTransport().alphaRhoPhi().group()
225  ),
226  -fvc::interpolate(this->alpha()*this->kappaEff())
227  *fvc::snGrad(this->thermo().T())
228  )
229  );
230 
232  const label d = composition.defaultSpecie();
233 
234  const PtrList<volScalarField>& Y = composition.Y();
235 
236  if (Y.size())
237  {
238  surfaceScalarField sumJ
239  (
241  (
242  "sumJ",
243  Y[0].mesh(),
245  )
246  );
247 
248  surfaceScalarField sumJh
249  (
251  (
252  "sumJh",
253  Y[0].mesh(),
255  )
256  );
257 
258  forAll(Y, i)
259  {
260  if (i != d)
261  {
262  const volScalarField hi
263  (
264  composition.Hs(i, this->thermo().p(), this->thermo().T())
265  );
266 
267  const surfaceScalarField ji(this->j(Y[i]));
268  sumJ += ji;
269 
270  sumJh += ji*fvc::interpolate(hi);
271  }
272  }
273 
274  {
275  const label i = d;
276 
277  const volScalarField hi
278  (
279  composition.Hs(i, this->thermo().p(), this->thermo().T())
280  );
281 
282  sumJh -= sumJ*fvc::interpolate(hi);
283  }
284 
285  tmpq.ref() += sumJh;
286  }
287 
288  return tmpq;
289 }
290 
291 
292 template<class BasicThermophysicalTransportModel>
294 (
295  volScalarField& he
296 ) const
297 {
298  tmp<fvScalarMatrix> tmpDivq
299  (
300  fvm::Su
301  (
302  -fvc::laplacian(this->alpha()*this->kappaEff(), this->thermo().T()),
303  he
304  )
305  );
306 
308  const label d = composition.defaultSpecie();
309 
310  const PtrList<volScalarField>& Y = composition.Y();
311 
312  tmpDivq.ref() -=
313  correction(fvm::laplacian(this->alpha()*this->alphaEff(), he));
314 
315  surfaceScalarField sumJ
316  (
318  (
319  "sumJ",
320  he.mesh(),
322  )
323  );
324 
325  surfaceScalarField sumJh
326  (
328  (
329  "sumJh",
330  he.mesh(),
331  dimensionedScalar(sumJ.dimensions()*he.dimensions(), 0)
332  )
333  );
334 
335  forAll(Y, i)
336  {
337  if (i != d)
338  {
339  const volScalarField hi
340  (
341  composition.Hs(i, this->thermo().p(), this->thermo().T())
342  );
343 
344  const surfaceScalarField ji(this->j(Y[i]));
345  sumJ += ji;
346 
347  sumJh += ji*fvc::interpolate(hi);
348  }
349  }
350 
351  {
352  const label i = d;
353 
354  const volScalarField hi
355  (
356  composition.Hs(i, this->thermo().p(), this->thermo().T())
357  );
358 
359  sumJh -= sumJ*fvc::interpolate(hi);
360  }
361 
362  tmpDivq.ref() += fvc::div(sumJh*he.mesh().magSf());
363 
364  return tmpDivq;
365 }
366 
367 
368 template<class BasicThermophysicalTransportModel>
370 (
371  const volScalarField& Yi
372 ) const
373 {
375  const label d = composition.defaultSpecie();
376 
377  if (composition.index(Yi) == d)
378  {
379  const PtrList<volScalarField>& Y = composition.Y();
380 
382  (
384  (
386  (
387  "j" + name(d),
388  this->momentumTransport().alphaRhoPhi().group()
389  ),
390  Yi.mesh(),
392  )
393  );
394 
395  surfaceScalarField& jd = tjd.ref();
396 
397  forAll(Y, i)
398  {
399  if (i != d)
400  {
401  jd -= this->j(Y[i]);
402  }
403  }
404 
405  return tjd;
406  }
407  else
408  {
409  return
410  BasicThermophysicalTransportModel::j(Yi)
411  + jexp_[composition.index(Yi)];
412  }
413 }
414 
415 
416 template<class BasicThermophysicalTransportModel>
418 (
419  volScalarField& Yi
420 ) const
421 {
423  return
424  BasicThermophysicalTransportModel::divj(Yi)
425  + fvc::div(jexp_[composition.index(Yi)]*Yi.mesh().magSf());
426 }
427 
428 
429 template<class BasicThermophysicalTransportModel>
432 {
434  const label d = composition.defaultSpecie();
435 
436  // Calculate the molecular weight of the mixture and the mole fractions
437  scalar Wm = 0;
438 
439  forAll(W, i)
440  {
441  X[i] = Y[i]/W[i];
442  Wm += X[i];
443  }
444 
445  Wm = 1/Wm;
446  X *= Wm;
447 
448  // i counter for the specie sub-system without the default specie
449  label is = 0;
450 
451  // Calculate the A and B matrices from the binary mass diffusion
452  // coefficients and specie mole fractions
453  forAll(X, i)
454  {
455  if (i != d)
456  {
457  A(is, is) = -X[i]*Wm/(DD(i, d)*W[d]);
458  B(is, is) = -(X[i]*Wm/W[d] + (1 - X[i])*Wm/W[i]);
459 
460  // j counter for the specie sub-system without the default specie
461  label js = 0;
462 
463  forAll(X, j)
464  {
465  if (j != i)
466  {
467  A(is, is) -= X[j]*Wm/(DD(i, j)*W[i]);
468 
469  if (j != d)
470  {
471  A(is, js) =
472  X[i]*(Wm/(DD(i, j)*W[j]) - Wm/(DD(i, d)*W[d]));
473 
474  B(is, js) = X[i]*(Wm/W[j] - Wm/W[d]);
475  }
476  }
477 
478  if (j != d)
479  {
480  js++;
481  }
482  }
483 
484  is++;
485  }
486  }
487 
488  // LU decompose A and invert
489  A.decompose();
490  A.inv(invA);
491 
492  // Calculate the generalised Fick's law diffusion coefficients
493  multiply(D, invA, B);
494 }
495 
496 
497 template<class BasicThermophysicalTransportModel>
500 {
502  const label d = composition.defaultSpecie();
503 
504  // For each cell or patch face
505  forAll(*(YPtrs[0]), pi)
506  {
507  forAll(W, i)
508  {
509  // Map YPtrs -> Y
510  Y[i] = (*YPtrs[i])[pi];
511 
512  // Map DijPtrs -> DD
513  forAll(W, j)
514  {
515  DD(i, j) = (*DijPtrs[i][j])[pi];
516  }
517  }
518 
519  // Transform DD -> D
520  transformDiffusionCoefficient();
521 
522  // i counter for the specie sub-system without the default specie
523  label is = 0;
524 
525  forAll(W, i)
526  {
527  if (i != d)
528  {
529  // j counter for the specie sub-system
530  // without the default specie
531  label js = 0;
532 
533  // Map D -> DijPtrs
534  forAll(W, j)
535  {
536  if (j != d)
537  {
538  (*DijPtrs[i][j])[pi] = D(is, js);
539 
540  js++;
541  }
542  }
543 
544  is++;
545  }
546  }
547  }
548 }
549 
550 
551 template<class BasicThermophysicalTransportModel>
553 (
555 )
556 {
558  const PtrList<volScalarField>& Y = composition.Y();
559  const volScalarField& Y0 = Y[0];
560 
561  forAll(W, i)
562  {
563  // Map composition.Y() internal fields -> YPtrs
564  YPtrs[i] = &Y[i].primitiveField();
565 
566  // Map Dii_ internal fields -> DijPtrs
567  DijPtrs[i][i] = &Dii_[i].primitiveFieldRef();
568 
569  // Map Dij internal fields -> DijPtrs
570  forAll(W, j)
571  {
572  if (j != i)
573  {
574  DijPtrs[i][j] = &Dij[i][j].primitiveFieldRef();
575  }
576  }
577  }
578 
579  // Transform binary mass diffusion coefficients internal field DijPtrs ->
580  // generalised Fick's law diffusion coefficients DijPtrs
581  transformDiffusionCoefficientFields();
582 
584  {
585  forAll(W, i)
586  {
587  // Map composition.Y() patch fields -> YPtrs
588  YPtrs[i] = &Y[i].boundaryField()[patchi];
589 
590  // Map Dii_ patch fields -> DijPtrs
591  DijPtrs[i][i] = &Dii_[i].boundaryFieldRef()[patchi];
592 
593  // Map Dij patch fields -> DijPtrs
594  forAll(W, j)
595  {
596  if (j != i)
597  {
598  DijPtrs[i][j] = &Dij[i][j].boundaryFieldRef()[patchi];
599  }
600  }
601  }
602 
603  // Transform binary mass diffusion coefficients patch field DijPtrs ->
604  // generalised Fick's law diffusion coefficients DijPtrs
605  transformDiffusionCoefficientFields();
606  }
607 }
608 
609 
610 template<class BasicThermophysicalTransportModel>
612 {
614 
616  const label d = composition.defaultSpecie();
617 
618  const PtrList<volScalarField>& Y = composition.Y();
619  const volScalarField& p = this->thermo().p();
620  const volScalarField& T = this->thermo().T();
621  const volScalarField& rho = this->momentumTransport().rho();
622 
624 
625  // Evaluate the specie binary mass diffusion coefficient functions
626  // and initialise the explicit part of the specie mass flux fields
627  forAll(Y, i)
628  {
629  if (i != d)
630  {
631  if (jexp_.set(i))
632  {
633  jexp_[i] = Zero;
634  }
635  else
636  {
637  jexp_.set
638  (
639  i,
641  (
642  "jexp" + Y[i].name(),
643  Y[i].mesh(),
644  dimensionedScalar(dimensionSet(1, -2, -1, 0, 0), 0)
645  )
646  );
647  }
648  }
649 
650  Dii_.set(i, evaluate(DFuncs_[i][i], dimViscosity, p, T));
651 
652  Dij[i].setSize(Y.size());
653 
654  forAll(Y, j)
655  {
656  if (j > i)
657  {
658  Dij[i].set(j, evaluate(DFuncs_[i][j], dimViscosity, p, T));
659  }
660  else if (j < i)
661  {
662  Dij[i].set(j, Dij[j][i].clone());
663  }
664  }
665  }
666 
667  //- Transform the binary mass diffusion coefficients into the
668  // the generalised Fick's law diffusion coefficients
669  transform(Dij);
670 
671  // Accumulate the explicit part of the specie mass flux fields
672  forAll(Y, j)
673  {
674  if (j != d)
675  {
676  const surfaceScalarField snGradYj(fvc::snGrad(Y[j]));
677 
678  forAll(Y, i)
679  {
680  if (i != d && i != j)
681  {
682  jexp_[i] -= fvc::interpolate(rho*Dij[i][j])*snGradYj;
683  }
684  }
685  }
686  }
687 
688  // Optionally add the Soret thermal diffusion contribution to the
689  // explicit part of the specie mass flux fields
690  if (DTFuncs_.size())
691  {
692  const surfaceScalarField gradTbyT(fvc::snGrad(T)/fvc::interpolate(T));
693 
694  forAll(Y, i)
695  {
696  if (i != d)
697  {
698  jexp_[i] -= fvc::interpolate
699  (
700  evaluate(DTFuncs_[i], dimDynamicViscosity, p, T)
701  )*gradTbyT;
702  }
703  }
704  }
705 }
706 
707 
708 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
709 
710 } // End namespace Foam
711 
712 // ************************************************************************* //
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
const char *const group
Group name for atomic constants.
const dimensionSet dimViscosity
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
const dimensionSet dimArea
fluidReactionThermo & thermo
Definition: createFields.H:28
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
basicSpecieMixture & composition
label index(const volScalarField &Yi) const
Return the specie index of the given mass-fraction field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual void correct()
Update the diffusion coefficients and flux corrections.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual scalar Wi(const label speciei) const =0
Molecular weight of the given specie [kg/kmol].
Calculate the snGrad of the given volField.
virtual scalar rho(const label speciei, const scalar p, const scalar T) const =0
Density [kg/m^3].
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
Specialisation of basicMixture for a mixture consisting of a number for molecular species...
fvMesh & mesh
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
virtual volScalarField & p()=0
Pressure [Pa].
virtual basicSpecieMixture & composition()=0
Return the composition of the multi-component mixture.
scalarList Y0(nSpecie, 0.0)
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
const dimensionSet dimTime
const dimensionSet & dimensions() const
Return dimensions.
volSymmTensorField invA(inv(I *UEqn.A()+drag->Dcu()))
Dimension set for the base types.
Definition: dimensionSet.H:121
T clone(const T &t)
Definition: List.H:54
MaxwellStefan(const word &type, const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from a momentum transport model and a thermo model.
Definition: MaxwellStefan.C:43
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
Calculate the laplacian of the given field.
virtual bool read()
Read thermophysicalTransport dictionary.
Definition: MaxwellStefan.C:94
static const zero Zero
Definition: zero.H:97
const dimensionSet dimDynamicViscosity
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
Calculate the divergence of the given field.
const Mesh & mesh() const
Return mesh.
const dimensionSet dimEnergy
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
virtual const volScalarField & T() const =0
Temperature [K].
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
A wordList with hashed indices for faster lookup by name.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
const scalarList W(::W(thermo))
PtrList< volScalarField > & Y
virtual scalar Hs(const label speciei, const scalar p, const scalar T) const =0
Sensible enthalpy [J/kg].
Run-time selectable function of two variables.
Definition: Function2.H:52
virtual tmp< surfaceScalarField > j(const volScalarField &Yi) const
Return the specie flux for the given specie mass-fraction [kg/m^2/s].
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:483
const speciesTable & species() const
Return the table of species.
label defaultSpecie() const
Return the index of the default specie.
IOerror FatalIOError