Fickian.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 "Fickian.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>
43 {
44  const PtrList<volScalarField>& Y = this->thermo().Y();
45  const volScalarField& p = this->thermo().p();
46  const volScalarField& T = this->thermo().T();
47 
48  Dm_.setSize(Y.size());
49 
50  if (mixtureDiffusionCoefficients_)
51  {
52  forAll(Y, i)
53  {
54  Dm_.set(i, evaluate(DmFuncs_[i], dimKinematicViscosity, p, T));
55  }
56  }
57  else
58  {
59  const volScalarField Wm(this->thermo().W());
60  volScalarField sumXbyD
61  (
63  (
64  "sumXbyD",
65  T.mesh(),
66  dimless/dimKinematicViscosity/Wm.dimensions()
67  )
68  );
69 
70  forAll(Dm_, i)
71  {
72  sumXbyD = Zero;
73 
74  forAll(Y, j)
75  {
76  if (j != i)
77  {
78  sumXbyD +=
79  Y[j]
80  /(
81  this->thermo().Wi(j)
82  *(
83  i < j
84  ? evaluate
85  (
86  DFuncs_[i][j],
88  p,
89  T
90  )
91  : evaluate
92  (
93  DFuncs_[j][i],
95  p,
96  T
97  )
98  )
99  );
100  }
101  }
102 
103  Dm_.set
104  (
105  i,
106  (1/Wm - Y[i]/this->thermo().Wi(i))
107  /max(sumXbyD, dimensionedScalar(sumXbyD.dimensions(), small))
108  );
109  }
110  }
111 }
112 
113 
114 template<class BasicThermophysicalTransportModel>
117 {
118  if (!Dm_.size())
119  {
120  updateDm();
121  }
122 
123  return Dm_;
124 }
125 
126 
127 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
128 
129 template<class BasicThermophysicalTransportModel>
131 (
132  const word& type,
133  const momentumTransportModel& momentumTransport,
134  const thermoModel& thermo
135 )
136 :
137  BasicThermophysicalTransportModel
138  (
139  type,
140  momentumTransport,
141  thermo
142  ),
143 
145 
146  mixtureDiffusionCoefficients_(true),
147 
148  DFuncs_(this->thermo().species().size()),
149 
150  DmFuncs_(this->thermo().species().size()),
151 
152  DTFuncs_
153  (
154  this->coeffDict().found("DT")
155  ? this->thermo().species().size()
156  : 0
157  )
158 {}
159 
160 
161 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162 
163 template<class BasicThermophysicalTransportModel>
165 {
167  {
168  const speciesTable& species = this->thermo().species();
169 
170  const dictionary& coeffDict = this->coeffDict();
171 
172  coeffDict.lookup("mixtureDiffusionCoefficients")
173  >> mixtureDiffusionCoefficients_;
174 
175  if (mixtureDiffusionCoefficients_)
176  {
177  const dictionary& Ddict = coeffDict.subDict("Dm");
178 
179  forAll(species, i)
180  {
181  DmFuncs_.set
182  (
183  i,
185  (
186  species[i],
187  dimPressure,
190  Ddict
191  ).ptr()
192  );
193  }
194  }
195  else
196  {
197  const dictionary& Ddict = coeffDict.subDict("D");
198 
199  // Read the array of specie binary mass diffusion coefficient
200  // functions
201  forAll(species, i)
202  {
203  DFuncs_[i].setSize(species.size());
204 
205  forAll(species, j)
206  {
207  if (j >= i)
208  {
209  const word nameij(species[i] + '-' + species[j]);
210  const word nameji(species[j] + '-' + species[i]);
211 
212  word Dname;
213 
214  if (Ddict.found(nameij) && Ddict.found(nameji))
215  {
216  if (i != j)
217  {
219  << "Binary mass diffusion coefficients "
220  "for Both " << nameij
221  << " and " << nameji << " provided, using "
222  << nameij << endl;
223  }
224 
225  Dname = nameij;
226  }
227  else if (Ddict.found(nameij))
228  {
229  Dname = nameij;
230  }
231  else if (Ddict.found(nameji))
232  {
233  Dname = nameji;
234  }
235  else
236  {
238  << "Binary mass diffusion coefficient for pair "
239  << nameij << " or " << nameji << " not provided"
240  << exit(FatalIOError);
241  }
242 
243  DFuncs_[i].set
244  (
245  j,
247  (
248  Dname,
249  dimPressure,
252  Ddict
253  ).ptr()
254  );
255  }
256  }
257  }
258  }
259 
260  // Optionally read the List of specie Soret thermal diffusion
261  // coefficient functions
262  if (coeffDict.found("DT"))
263  {
264  const dictionary& DTdict = coeffDict.subDict("DT");
265 
266  forAll(species, i)
267  {
268  DTFuncs_.set
269  (
270  i,
272  (
273  species[i],
274  dimPressure,
277  DTdict
278  ).ptr()
279  );
280  }
281  }
282 
283  return true;
284  }
285  else
286  {
287  return false;
288  }
289 }
290 
291 
292 template<class BasicThermophysicalTransportModel>
294 (
295  const volScalarField& Yi
296 ) const
297 {
298  return volScalarField::New
299  (
300  "DEff",
301  this->momentumTransport().rho()
302  *Dm()[this->thermo().specieIndex(Yi)]
303  );
304 }
305 
306 
307 template<class BasicThermophysicalTransportModel>
309 (
310  const volScalarField& Yi,
311  const label patchi
312 ) const
313 {
314  return
315  this->momentumTransport().rho().boundaryField()[patchi]
316  *Dm()[this->thermo().specieIndex(Yi)].boundaryField()[patchi];
317 }
318 
319 
320 template<class BasicThermophysicalTransportModel>
322 {
324  (
326  (
328  (
329  "q",
330  this->momentumTransport().alphaRhoPhi().group()
331  ),
332  -fvc::interpolate(this->alpha()*this->kappaEff())
333  *fvc::snGrad(this->thermo().T())
334  )
335  );
336 
337  const PtrList<volScalarField>& Y = this->thermo().Y();
338  const volScalarField& p = this->thermo().p();
339  const volScalarField& T = this->thermo().T();
340 
341  if (Y.size())
342  {
343  surfaceScalarField sumJ
344  (
346  (
347  "sumJ",
348  Y[0].mesh(),
350  )
351  );
352 
353  surfaceScalarField sumJh
354  (
356  (
357  "sumJh",
358  Y[0].mesh(),
360  )
361  );
362 
363  forAll(Y, i)
364  {
365  if (i != this->thermo().defaultSpecie())
366  {
367  const volScalarField hi(this->thermo().hsi(i, p, T));
368 
369  const surfaceScalarField ji(this->j(Y[i]));
370 
371  sumJ += ji;
372 
373  sumJh += ji*fvc::interpolate(hi);
374  }
375  }
376 
377  {
378  const label i = this->thermo().defaultSpecie();
379 
380  const volScalarField hi(this->thermo().hsi(i, p, T));
381 
382  sumJh -= sumJ*fvc::interpolate(hi);
383  }
384 
385  tmpq.ref() += sumJh;
386  }
387 
388  return tmpq;
389 }
390 
391 
392 template<class BasicThermophysicalTransportModel>
394 (
395  const label patchi
396 ) const
397 {
398  tmp<scalarField> tmpq
399  (
400  - (
401  this->alpha().boundaryField()[patchi]
402  *this->kappaEff(patchi)
403  *this->thermo().T().boundaryField()[patchi].snGrad()
404  )
405  );
406 
407  const PtrList<volScalarField>& Y = this->thermo().Y();
408  const scalarField& p = this->thermo().p().boundaryField()[patchi];
409  const scalarField& T = this->thermo().T().boundaryField()[patchi];
410 
411  if (Y.size())
412  {
413  scalarField sumJ(tmpq->size(), scalar(0));
414  scalarField sumJh(tmpq->size(), scalar(0));
415 
416  forAll(Y, i)
417  {
418  if (i != this->thermo().defaultSpecie())
419  {
420  const scalarField hi(this->thermo().hsi(i, p, T));
421 
422  const scalarField ji(this->j(Y[i], patchi));
423 
424  sumJ += ji;
425 
426  sumJh += ji*hi;
427  }
428  }
429 
430  {
431  const label i = this->thermo().defaultSpecie();
432 
433  const scalarField hi(this->thermo().hsi(i, p, T));
434 
435  sumJh -= sumJ*hi;
436  }
437 
438  tmpq.ref() += sumJh;
439  }
440 
441  return tmpq;
442 }
443 
444 
445 template<class BasicThermophysicalTransportModel>
447 (
449 ) const
450 {
451  tmp<fvScalarMatrix> tmpDivq
452  (
453  fvm::Su
454  (
455  -fvc::laplacian(this->alpha()*this->kappaEff(), this->thermo().T()),
456  he
457  )
458  );
459 
460  const PtrList<volScalarField>& Y = this->thermo().Y();
461  const volScalarField& p = this->thermo().p();
462  const volScalarField& T = this->thermo().T();
463 
464  tmpDivq.ref() -=
465  fvm::laplacianCorrection(this->alpha()*this->alphaEff(), he);
466 
467  surfaceScalarField sumJ
468  (
470  (
471  "sumJ",
472  he.mesh(),
474  )
475  );
476 
477  surfaceScalarField sumJh
478  (
480  (
481  "sumJh",
482  he.mesh(),
483  dimensionedScalar(sumJ.dimensions()*he.dimensions(), 0)
484  )
485  );
486 
487  forAll(Y, i)
488  {
489  if (i != this->thermo().defaultSpecie())
490  {
491  const volScalarField hi(this->thermo().hsi(i, p, T));
492 
493  const surfaceScalarField ji(this->j(Y[i]));
494 
495  sumJ += ji;
496 
497  sumJh += ji*fvc::interpolate(hi);
498  }
499  }
500 
501  {
502  const label i = this->thermo().defaultSpecie();
503 
504  const volScalarField hi(this->thermo().hsi(i, p, T));
505 
506  sumJh -= sumJ*fvc::interpolate(hi);
507  }
508 
509  tmpDivq.ref() += fvc::div(sumJh*he.mesh().magSf());
510 
511  return tmpDivq;
512 }
513 
514 
515 template<class BasicThermophysicalTransportModel>
517 (
518  const volScalarField& Yi
519 ) const
520 {
521  if (DTFuncs_.size())
522  {
523  const volScalarField& p = this->thermo().p();
524  const volScalarField& T = this->thermo().T();
525 
526  return
527  BasicThermophysicalTransportModel::j(Yi)
529  (
530  evaluate
531  (
532  DTFuncs_[this->thermo().specieIndex(Yi)],
534  p,
535  T
536  )
537  )
539  }
540  else
541  {
542  return BasicThermophysicalTransportModel::j(Yi);
543  }
544 }
545 
546 
547 template<class BasicThermophysicalTransportModel>
549 (
550  const volScalarField& Yi,
551  const label patchi
552 ) const
553 {
554  if (DTFuncs_.size())
555  {
556  const scalarField& p = this->thermo().p().boundaryField()[patchi];
557  const scalarField& T = this->thermo().T().boundaryField()[patchi];
558 
559  return
560  BasicThermophysicalTransportModel::j(Yi, patchi)
561  - DTFuncs_[this->thermo().specieIndex(Yi)].value(p, T)
562  *this->thermo().T().boundaryField()[patchi].snGrad()/T;
563  }
564  else
565  {
566  return BasicThermophysicalTransportModel::j(Yi, patchi);
567  }
568 }
569 
570 
571 template<class BasicThermophysicalTransportModel>
573 (
574  volScalarField& Yi
575 ) const
576 {
577  if (DTFuncs_.size())
578  {
579  const volScalarField& p = this->thermo().p();
580  const volScalarField& T = this->thermo().T();
581 
582  return
583  BasicThermophysicalTransportModel::divj(Yi)
584  - fvc::div
585  (
587  (
588  evaluate
589  (
590  DTFuncs_[this->thermo().specieIndex(Yi)],
592  p,
593  T
594  )
595  )
597  *T.mesh().magSf()
598  );
599  }
600  else
601  {
602  return BasicThermophysicalTransportModel::divj(Yi);
603  }
604 }
605 
606 
607 template<class BasicThermophysicalTransportModel>
609 {
610  BasicThermophysicalTransportModel::predict();
611  updateDm();
612 }
613 
614 
615 template<class BasicThermophysicalTransportModel>
617 {
618  return true;
619 }
620 
621 
622 template<class BasicThermophysicalTransportModel>
624 (
625  const polyTopoChangeMap& map
626 )
627 {
628  // Delete the cached Dm, will be re-created in predict
629  Dm_.clear();
630 }
631 
632 
633 template<class BasicThermophysicalTransportModel>
635 (
636  const polyMeshMap& map
637 )
638 {
639  // Delete the cached Dm, will be re-created in predict
640  Dm_.clear();
641 }
642 
643 
644 template<class BasicThermophysicalTransportModel>
646 (
647  const polyDistributionMap& map
648 )
649 {
650  // Delete the cached Dm, will be re-created in predict
651  Dm_.clear();
652 }
653 
654 
655 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
656 
657 } // End namespace Foam
658 
659 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
const dimensionSet & dimensions() const
Return dimensions.
virtual bool movePoints()
Update for mesh motion.
Definition: Fickian.C:616
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
Definition: Fickian.C:294
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Definition: Fickian.C:447
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
Definition: Fickian.C:573
virtual void predict()
Update the diffusion coefficients.
Definition: Fickian.C:608
virtual void mapMesh(const polyMeshMap &map)
Update from another mesh using the given map.
Definition: Fickian.C:635
virtual tmp< surfaceScalarField > j(const volScalarField &Yi) const
Return the specie flux for the given specie mass-fraction [kg/m^2/s].
Definition: Fickian.C:517
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
Definition: Fickian.C:321
virtual void distribute(const polyDistributionMap &map)
Redistribute or update using the given distribution map.
Definition: Fickian.C:646
BasicThermophysicalTransportModel::thermoModel thermoModel
Definition: Fickian.H:105
BasicThermophysicalTransportModel::momentumTransportModel momentumTransportModel
Definition: Fickian.H:102
Fickian(const word &type, const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from a momentum transport model and a thermo model.
Definition: Fickian.C:131
virtual void topoChange(const polyTopoChangeMap &map)
Update topology using the given map.
Definition: Fickian.C:624
virtual bool read()
Read thermophysicalTransport dictionary.
Definition: Fickian.C:164
Run-time selectable function of two variables.
Definition: Function2.H:98
static autoPtr< Function2< Type > > New(const word &name, const Function2s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function2New.C:32
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)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
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
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
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
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 dimensionedScalar Wi(const label speciei) const =0
Molecular weight [kg/kmol].
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
Base class for multi-component Fickian based temperature gradient heat flux models with optional Sore...
A class for handling words, derived from string.
Definition: word.H:62
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
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
const dimensionSet dimless
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTemperature
const dimensionSet dimTime
scalarList W(const fluidMulticomponentThermo &thermo)
void evaluate(GeometricField< Type, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, GeoMesh > &x)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
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()
volScalarField & p
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31