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-2024 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40 
41 template<class BasicThermophysicalTransportModel>
42 void MaxwellStefan<BasicThermophysicalTransportModel>::
43 transformDiffusionCoefficient() const
44 {
45  const label d = this->thermo().defaultSpecie();
46 
47  // Calculate the molecular weight of the mixture and the mole fractions
48  scalar Wm = 0;
49 
50  forAll(W, i)
51  {
52  X[i] = Y[i]/W[i];
53  Wm += X[i];
54  }
55 
56  Wm = 1/Wm;
57  X *= Wm;
58 
59  // i counter for the specie sub-system without the default specie
60  label is = 0;
61 
62  // Calculate the A and B matrices from the binary mass diffusion
63  // coefficients and specie mole fractions
64  forAll(X, i)
65  {
66  if (i != d)
67  {
68  A(is, is) = -X[i]*Wm/(DD(i, d)*W[d]);
69  B(is, is) = -(X[i]*Wm/W[d] + (1 - X[i])*Wm/W[i]);
70 
71  // j counter for the specie sub-system without the default specie
72  label js = 0;
73 
74  forAll(X, j)
75  {
76  if (j != i)
77  {
78  A(is, is) -= X[j]*Wm/(DD(i, j)*W[i]);
79 
80  if (j != d)
81  {
82  A(is, js) =
83  X[i]*(Wm/(DD(i, j)*W[j]) - Wm/(DD(i, d)*W[d]));
84 
85  B(is, js) = X[i]*(Wm/W[j] - Wm/W[d]);
86  }
87  }
88 
89  if (j != d)
90  {
91  js++;
92  }
93  }
94 
95  is++;
96  }
97  }
98 
99  // LU decompose A and invert
100  A.decompose();
101  A.inv(invA);
102 
103  // Calculate the generalised Fick's law diffusion coefficients
104  multiply(D, invA, B);
105 }
106 
107 
108 template<class BasicThermophysicalTransportModel>
109 void MaxwellStefan<BasicThermophysicalTransportModel>::
110 transformDiffusionCoefficientFields() const
111 {
112  const label d = this->thermo().defaultSpecie();
113 
114  // For each cell or patch face
115  forAll(*(YPtrs[0]), pi)
116  {
117  forAll(W, i)
118  {
119  // Map YPtrs -> Y
120  Y[i] = (*YPtrs[i])[pi];
121 
122  // Map DijPtrs -> DD
123  forAll(W, j)
124  {
125  DD(i, j) = (*DijPtrs[i][j])[pi];
126  }
127  }
128 
129  // Transform DD -> D
130  transformDiffusionCoefficient();
131 
132  // i counter for the specie sub-system without the default specie
133  label is = 0;
134 
135  forAll(W, i)
136  {
137  if (i != d)
138  {
139  // j counter for the specie sub-system
140  // without the default specie
141  label js = 0;
142 
143  // Map D -> DijPtrs
144  forAll(W, j)
145  {
146  if (j != d)
147  {
148  (*DijPtrs[i][j])[pi] = D(is, js);
149 
150  js++;
151  }
152  }
153 
154  is++;
155  }
156  }
157  }
158 }
159 
160 
161 template<class BasicThermophysicalTransportModel>
162 void MaxwellStefan<BasicThermophysicalTransportModel>::transform
163 (
164  List<PtrList<volScalarField>>& Dij
165 ) const
166 {
167  const PtrList<volScalarField>& Y = this->thermo().Y();
168  const volScalarField& Y0 = Y[0];
169 
170  forAll(W, i)
171  {
172  // Map this->thermo().Y() internal fields -> YPtrs
173  YPtrs[i] = &Y[i].primitiveField();
174 
175  // Map Dii_ internal fields -> DijPtrs
176  DijPtrs[i][i] = &Dii_[i].primitiveFieldRef();
177 
178  // Map Dij internal fields -> DijPtrs
179  forAll(W, j)
180  {
181  if (j != i)
182  {
183  DijPtrs[i][j] = &Dij[i][j].primitiveFieldRef();
184  }
185  }
186  }
187 
188  // Transform binary mass diffusion coefficients internal field DijPtrs ->
189  // generalised Fick's law diffusion coefficients DijPtrs
190  transformDiffusionCoefficientFields();
191 
192  forAll(Y0.boundaryField(), patchi)
193  {
194  forAll(W, i)
195  {
196  // Map this->thermo().Y() patch fields -> YPtrs
197  YPtrs[i] = &Y[i].boundaryField()[patchi];
198 
199  // Map Dii_ patch fields -> DijPtrs
200  DijPtrs[i][i] = &Dii_[i].boundaryFieldRef()[patchi];
201 
202  // Map Dij patch fields -> DijPtrs
203  forAll(W, j)
204  {
205  if (j != i)
206  {
207  DijPtrs[i][j] = &Dij[i][j].boundaryFieldRef()[patchi];
208  }
209  }
210  }
211 
212  // Transform binary mass diffusion coefficients patch field DijPtrs ->
213  // generalised Fick's law diffusion coefficients DijPtrs
214  transformDiffusionCoefficientFields();
215  }
216 }
217 
218 
219 template<class BasicThermophysicalTransportModel>
221 {
222  const label d = this->thermo().defaultSpecie();
223 
224  const PtrList<volScalarField>& Y = this->thermo().Y();
225  const volScalarField& p = this->thermo().p();
226  const volScalarField& T = this->thermo().T();
227  const volScalarField& rho = this->momentumTransport().rho();
228 
229  Dii_.setSize(Y.size());
230  jexp_.setSize(Y.size());
231 
233 
234  // Evaluate the specie binary mass diffusion coefficient functions
235  // and initialise the explicit part of the specie mass flux fields
236  forAll(Y, i)
237  {
238  if (i != d)
239  {
240  if (jexp_.set(i))
241  {
242  jexp_[i] = Zero;
243  }
244  else
245  {
246  jexp_.set
247  (
248  i,
250  (
251  "jexp" + Y[i].name(),
252  Y[i].mesh(),
253  dimensionedScalar(dimensionSet(1, -2, -1, 0, 0), 0)
254  )
255  );
256  }
257  }
258 
259  Dii_.set(i, evaluate(DFuncs_[i][i], dimKinematicViscosity, p, T));
260 
261  Dij[i].setSize(Y.size());
262 
263  forAll(Y, j)
264  {
265  if (j > i)
266  {
267  Dij[i].set
268  (
269  j,
270  evaluate(DFuncs_[i][j], dimKinematicViscosity, p, T)
271  );
272  }
273  else if (j < i)
274  {
275  Dij[i].set(j, Dij[j][i].clone());
276  }
277  }
278  }
279 
280  //- Transform the binary mass diffusion coefficients into the
281  // the generalised Fick's law diffusion coefficients
282  transform(Dij);
283 
284  // Accumulate the explicit part of the specie mass flux fields
285  forAll(Y, j)
286  {
287  if (j != d)
288  {
289  const surfaceScalarField snGradYj(fvc::snGrad(Y[j]));
290 
291  forAll(Y, i)
292  {
293  if (i != d && i != j)
294  {
295  jexp_[i] -= fvc::interpolate(rho*Dij[i][j])*snGradYj;
296  }
297  }
298  }
299  }
300 
301  // Optionally add the Soret thermal diffusion contribution to the
302  // explicit part of the specie mass flux fields
303  if (DTFuncs_.size())
304  {
306 
307  forAll(Y, i)
308  {
309  if (i != d)
310  {
311  jexp_[i] -= fvc::interpolate
312  (
313  evaluate(DTFuncs_[i], dimDynamicViscosity, p, T)
314  )*gradTbyT;
315  }
316  }
317  }
318 }
319 
320 
321 template<class BasicThermophysicalTransportModel>
322 const PtrList<volScalarField>&
323 MaxwellStefan<BasicThermophysicalTransportModel>::Dii() const
324 {
325  if (!Dii_.size())
326  {
327  updateDii();
328  }
329 
330  return Dii_;
331 }
332 
333 
334 template<class BasicThermophysicalTransportModel>
335 const PtrList<surfaceScalarField>&
336 MaxwellStefan<BasicThermophysicalTransportModel>::jexp() const
337 {
338  if (!jexp_.size())
339  {
340  updateDii();
341  }
342 
343  return jexp_;
344 }
345 
346 
347 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
348 
349 template<class BasicThermophysicalTransportModel>
351 (
352  const word& type,
353  const momentumTransportModel& momentumTransport,
354  const thermoModel& thermo
355 )
356 :
357  BasicThermophysicalTransportModel
358  (
359  type,
360  momentumTransport,
361  thermo
362  ),
363 
365 
366  DFuncs_(this->thermo().species().size()),
367 
368  DTFuncs_
369  (
370  this->coeffDict().found("DT")
371  ? this->thermo().species().size()
372  : 0
373  ),
374 
375  W(this->thermo().species().size()),
376 
377  YPtrs(W.size()),
378  DijPtrs(W.size()),
379 
380  Y(W.size()),
381  X(W.size()),
382  DD(W.size()),
383  A(W.size() - 1),
384  B(A.m()),
385  invA(A.m()),
386  D(W.size())
387 {
388  // Set the molecular weights of the species
389  forAll(W, i)
390  {
391  W[i] = this->thermo().Wi(i).value();
392  }
393 }
394 
395 
396 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
397 
398 template<class BasicThermophysicalTransportModel>
400 {
401  if
402  (
404  )
405  {
406  const speciesTable& species = this->thermo().species();
407 
408  const dictionary& coeffDict = this->coeffDict();
409  const dictionary& Ddict = coeffDict.subDict("D");
410 
411  // Read the array of specie binary mass diffusion coefficient functions
412  forAll(species, i)
413  {
414  DFuncs_[i].setSize(species.size());
415 
416  forAll(species, j)
417  {
418  if (j >= i)
419  {
420  const word nameij(species[i] + '-' + species[j]);
421  const word nameji(species[j] + '-' + species[i]);
422 
423  word Dname;
424 
425  if (Ddict.found(nameij) && Ddict.found(nameji))
426  {
427  if (i != j)
428  {
430  << "Binary mass diffusion coefficients "
431  "for both " << nameij << " and " << nameji
432  << " provided, using " << nameij << endl;
433  }
434 
435  Dname = nameij;
436  }
437  else if (Ddict.found(nameij))
438  {
439  Dname = nameij;
440  }
441  else if (Ddict.found(nameji))
442  {
443  Dname = nameji;
444  }
445  else
446  {
448  << "Binary mass diffusion coefficients for pair "
449  << nameij << " or " << nameji << " not provided"
450  << exit(FatalIOError);
451  }
452 
453  DFuncs_[i].set
454  (
455  j,
457  (
458  Dname,
459  dimPressure,
462  Ddict
463  ).ptr()
464  );
465  }
466  }
467  }
468 
469  // Optionally read the List of specie Soret thermal diffusion
470  // coefficient functions
471  if (coeffDict.found("DT"))
472  {
473  const dictionary& DTdict = coeffDict.subDict("DT");
474 
475  forAll(species, i)
476  {
477  DTFuncs_.set
478  (
479  i,
481  (
482  species[i],
483  dimPressure,
486  DTdict
487  ).ptr()
488  );
489  }
490  }
491 
492  return true;
493  }
494  else
495  {
496  return false;
497  }
498 }
499 
500 
501 template<class BasicThermophysicalTransportModel>
503 (
504  const volScalarField& Yi
505 ) const
506 {
507  return volScalarField::New
508  (
509  "DEff",
510  this->momentumTransport().rho()*Dii()[this->thermo().specieIndex(Yi)]
511  );
512 }
513 
514 
515 template<class BasicThermophysicalTransportModel>
517 (
518  const volScalarField& Yi,
519  const label patchi
520 ) const
521 {
522  return
523  this->momentumTransport().rho().boundaryField()[patchi]
524  *Dii()[this->thermo().specieIndex(Yi)].boundaryField()[patchi];
525 }
526 
527 
528 template<class BasicThermophysicalTransportModel>
531 {
533  (
535  (
537  (
538  "q",
539  this->momentumTransport().alphaRhoPhi().group()
540  ),
541  -fvc::interpolate(this->alpha()*this->kappaEff())
542  *fvc::snGrad(this->thermo().T())
543  )
544  );
545 
546  const label d = this->thermo().defaultSpecie();
547 
548  const PtrList<volScalarField>& Y = this->thermo().Y();
549  const volScalarField& p = this->thermo().p();
550  const volScalarField& T = this->thermo().T();
551 
552  if (Y.size())
553  {
554  surfaceScalarField sumJ
555  (
557  (
558  "sumJ",
559  Y[0].mesh(),
561  )
562  );
563 
564  surfaceScalarField sumJh
565  (
567  (
568  "sumJh",
569  Y[0].mesh(),
571  )
572  );
573 
574  forAll(Y, i)
575  {
576  if (i != d)
577  {
578  const volScalarField hi(this->thermo().hsi(i, p, T));
579 
580  const surfaceScalarField ji(this->j(Y[i]));
581 
582  sumJ += ji;
583 
584  sumJh += ji*fvc::interpolate(hi);
585  }
586  }
587 
588  {
589  const label i = d;
590 
591  const volScalarField hi(this->thermo().hsi(i, p, T));
592 
593  sumJh -= sumJ*fvc::interpolate(hi);
594  }
595 
596  tmpq.ref() += sumJh;
597  }
598 
599  return tmpq;
600 }
601 
602 
603 template<class BasicThermophysicalTransportModel>
605 (
606  const label patchi
607 ) const
608 {
609  tmp<scalarField> tmpq
610  (
611  - (
612  this->alpha().boundaryField()[patchi]
613  *this->kappaEff(patchi)
614  *this->thermo().T().boundaryField()[patchi]
615  )
616  );
617 
618  const label d = this->thermo().defaultSpecie();
619 
620  const PtrList<volScalarField>& Y = this->thermo().Y();
621  const volScalarField& p = this->thermo().p();
622  const volScalarField& T = this->thermo().T();
623 
624  if (Y.size())
625  {
626  scalarField sumJ(tmpq->size(), scalar(0));
627  scalarField sumJh(tmpq->size(), scalar(0));
628 
629  forAll(Y, i)
630  {
631  if (i != d)
632  {
633  const scalarField hi(this->thermo().hsi(i, p, T));
634 
635  const scalarField ji(this->j(Y[i], patchi));
636 
637  sumJ += ji;
638 
639  sumJh += ji*hi;
640  }
641  }
642 
643  {
644  const label i = d;
645 
646  const scalarField hi(this->thermo().hsi(i, p, T));
647 
648  sumJh -= sumJ*hi;
649  }
650 
651  tmpq.ref() += sumJh;
652  }
653 
654  return tmpq;
655 }
656 
657 
658 template<class BasicThermophysicalTransportModel>
660 (
662 ) const
663 {
664  tmp<fvScalarMatrix> tmpDivq
665  (
666  fvm::Su
667  (
668  -fvc::laplacian(this->alpha()*this->kappaEff(), this->thermo().T()),
669  he
670  )
671  );
672 
673  const label d = this->thermo().defaultSpecie();
674 
675  const PtrList<volScalarField>& Y = this->thermo().Y();
676  const volScalarField& p = this->thermo().p();
677  const volScalarField& T = this->thermo().T();
678 
679  tmpDivq.ref() -=
680  fvm::laplacianCorrection(this->alpha()*this->alphaEff(), he);
681 
682  surfaceScalarField sumJ
683  (
685  (
686  "sumJ",
687  he.mesh(),
689  )
690  );
691 
692  surfaceScalarField sumJh
693  (
695  (
696  "sumJh",
697  he.mesh(),
698  dimensionedScalar(sumJ.dimensions()*he.dimensions(), 0)
699  )
700  );
701 
702  forAll(Y, i)
703  {
704  if (i != d)
705  {
706  const volScalarField hi(this->thermo().hsi(i, p, T));
707 
708  const surfaceScalarField ji(this->j(Y[i]));
709 
710  sumJ += ji;
711 
712  sumJh += ji*fvc::interpolate(hi);
713  }
714  }
715 
716  {
717  const label i = d;
718 
719  const volScalarField hi(this->thermo().hsi(i, p, T));
720 
721  sumJh -= sumJ*fvc::interpolate(hi);
722  }
723 
724  tmpDivq.ref() += fvc::div(sumJh*he.mesh().magSf());
725 
726  return tmpDivq;
727 }
728 
729 
730 template<class BasicThermophysicalTransportModel>
732 (
733  const volScalarField& Yi
734 ) const
735 {
736  const label d = this->thermo().defaultSpecie();
737 
738  if (this->thermo().specieIndex(Yi) == d)
739  {
740  const PtrList<volScalarField>& Y = this->thermo().Y();
741 
743  (
745  (
747  (
748  "j" + name(d),
749  this->momentumTransport().alphaRhoPhi().group()
750  ),
751  Yi.mesh(),
753  )
754  );
755 
756  surfaceScalarField& jd = tjd.ref();
757 
758  forAll(Y, i)
759  {
760  if (i != d)
761  {
762  jd -= this->j(Y[i]);
763  }
764  }
765 
766  return tjd;
767  }
768  else
769  {
770  return
771  BasicThermophysicalTransportModel::j(Yi)
772  + jexp()[this->thermo().specieIndex(Yi)];
773  }
774 }
775 
776 
777 template<class BasicThermophysicalTransportModel>
779 (
780  const volScalarField& Yi,
781  const label patchi
782 ) const
783 {
784  const label d = this->thermo().defaultSpecie();
785 
786  if (this->thermo().specieIndex(Yi) == d)
787  {
788  const PtrList<volScalarField>& Y = this->thermo().Y();
789 
790  tmp<scalarField> tjd
791  (
792  new scalarField(Yi.boundaryField()[patchi].size(), scalar(0))
793  );
794  scalarField& jd = tjd.ref();
795 
796  forAll(Y, i)
797  {
798  if (i != d)
799  {
800  jd -= this->j(Y[i], patchi);
801  }
802  }
803 
804  return tjd;
805  }
806  else
807  {
808  return
809  BasicThermophysicalTransportModel::j(Yi, patchi)
810  + jexp()[this->thermo().specieIndex(Yi)].boundaryField()[patchi];
811  }
812 }
813 
814 
815 template<class BasicThermophysicalTransportModel>
817 (
818  volScalarField& Yi
819 ) const
820 {
821  return
822  BasicThermophysicalTransportModel::divj(Yi)
823  + fvc::div(jexp()[this->thermo().specieIndex(Yi)]*Yi.mesh().magSf());
824 }
825 
826 
827 template<class BasicThermophysicalTransportModel>
829 {
830  BasicThermophysicalTransportModel::predict();
831  updateDii();
832 }
833 
834 
835 template<class BasicThermophysicalTransportModel>
837 {
838  return true;
839 }
840 
841 
842 template<class BasicThermophysicalTransportModel>
844 (
845  const polyTopoChangeMap& map
846 )
847 {
848  // Delete the cached Dii and jexp, will be re-created in predict
849  Dii_.clear();
850  jexp_.clear();
851 }
852 
853 
854 template<class BasicThermophysicalTransportModel>
856 (
857  const polyMeshMap& map
858 )
859 {
860  // Delete the cached Dii and jexp, will be re-created in predict
861  Dii_.clear();
862  jexp_.clear();
863 }
864 
865 
866 template<class BasicThermophysicalTransportModel>
868 (
869  const polyDistributionMap& map
870 )
871 {
872  // Delete the cached Dii and jexp, will be re-created in predict
873  Dii_.clear();
874  jexp_.clear();
875 }
876 
877 
878 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
879 
880 } // End namespace Foam
881 
882 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Run-time selectable function of two variables.
Definition: Function2.H:98
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
static word groupName(Name name, const word &group)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Base class for multi-component Maxwell Stefan generalised Fick's law diffusion coefficients based tem...
Definition: MaxwellStefan.H:76
virtual bool movePoints()
Update for mesh motion.
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
virtual void predict()
Update the diffusion coefficients and flux corrections.
virtual void mapMesh(const polyMeshMap &map)
Update from another mesh using the given map.
virtual tmp< surfaceScalarField > j(const volScalarField &Yi) const
Return the specie flux for the given specie mass-fraction [kg/m^2/s].
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
virtual void distribute(const polyDistributionMap &map)
Redistribute or update using the given distribution map.
MaxwellStefan(const word &type, const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from a momentum transport model and a thermo model.
BasicThermophysicalTransportModel::thermoModel thermoModel
BasicThermophysicalTransportModel::momentumTransportModel momentumTransportModel
virtual void topoChange(const polyTopoChangeMap &map)
Update topology using the given map.
virtual bool read()
Read thermophysicalTransport dictionary.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
virtual const volScalarField & T() const =0
Temperature [K].
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:849
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
Dimension set for the base types.
Definition: dimensionSet.H:125
virtual const volScalarField & p() const =0
Pressure [Pa].
A wordList with hashed indices for faster lookup by name.
virtual const speciesTable & species() const =0
The table of species.
virtual PtrList< volScalarField > & Y()=0
Access the mass-fraction fields.
virtual label defaultSpecie() const =0
The index of the default specie.
label specieIndex(const volScalarField &Yi) const
Access the specie index of the given mass-fraction field.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
volSymmTensorField invA(inv(I *UEqn.A()+drag->Dcu()))
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
Calculate the divergence of the given field.
Calculate the laplacian of the given field.
Calculate the snGrad of the given volField.
Calculate the matrix for implicit and explicit sources.
label patchi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
#define WarningInFunction
Report a warning using Foam::Warning.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
const char *const group
Group name for atomic constants.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const VolField< Type > &)
tmp< fvMatrix< Type > > laplacianCorrection(const VolField< scalar > &gamma, const VolField< Type > &vf)
Definition: fvmLaplacian.C:340
static const coefficient D("D", dimTemperature, 257.14)
static const coefficient B("B", dimless, 18.678)
static const coefficient A("A", dimPressure, 611.21)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
const dimensionSet dimDynamicViscosity
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
const dimensionSet dimEnergy
const dimensionSet dimPressure
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimKinematicViscosity
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
SurfaceField< scalar > surfaceScalarField
const dimensionSet dimTemperature
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTime
scalarList W(const fluidMulticomponentThermo &thermo)
void evaluate(GeometricField< Type, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, GeoMesh > &x)
void multiply(LagrangianPatchField< Type > &f, const LagrangianPatchField< scalar > &f1, const LagrangianPatchField< Type > &f2)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:501
T clone(const T &t)
Definition: List.H:55
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
const dimensionSet dimArea
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
thermo he()
scalarList Y0(nSpecie, 0.0)
volScalarField & p
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31