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-2023 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 basicSpecieMixture& composition = this->thermo().composition();
45  const PtrList<volScalarField>& Y = composition.Y();
46  const volScalarField& p = this->thermo().p();
47  const volScalarField& T = this->thermo().T();
48 
49  Dm_.setSize(Y.size());
50 
51  if (mixtureDiffusionCoefficients_)
52  {
53  forAll(Y, i)
54  {
55  Dm_.set(i, evaluate(DmFuncs_[i], dimViscosity, p, T));
56  }
57  }
58  else
59  {
60  const volScalarField Wm(this->thermo().W());
61  volScalarField sumXbyD
62  (
64  (
65  "sumXbyD",
66  T.mesh(),
67  dimless/dimViscosity/Wm.dimensions()
68  )
69  );
70 
71  forAll(Dm_, i)
72  {
73  sumXbyD = Zero;
74 
75  forAll(Y, j)
76  {
77  if (j != i)
78  {
79  sumXbyD +=
80  Y[j]
81  /(
83  (
84  "Wj",
85  Wm.dimensions(),
86  composition.Wi(j)
87  )
88  *(
89  i < j
90  ? evaluate(DFuncs_[i][j], dimViscosity, p, T)
91  : evaluate(DFuncs_[j][i], dimViscosity, p, T)
92  )
93  );
94  }
95  }
96 
97  Dm_.set
98  (
99  i,
100  (
101  1/Wm
102  - Y[i]
103  /dimensionedScalar("Wi", Wm.dimensions(), composition.Wi(i))
104  )/max(sumXbyD, dimensionedScalar(sumXbyD.dimensions(), small))
105  );
106  }
107  }
108 }
109 
110 
111 template<class BasicThermophysicalTransportModel>
114 {
115  if (!Dm_.size())
116  {
117  updateDm();
118  }
119 
120  return Dm_;
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
125 
126 template<class BasicThermophysicalTransportModel>
128 (
129  const word& type,
130  const momentumTransportModel& momentumTransport,
132 )
133 :
134  BasicThermophysicalTransportModel
135  (
136  type,
137  momentumTransport,
138  thermo
139  ),
140 
141  UpdateableMeshObject(*this, thermo.mesh()),
142 
143  mixtureDiffusionCoefficients_(true),
144 
145  DFuncs_(this->thermo().composition().species().size()),
146 
147  DmFuncs_(this->thermo().composition().species().size()),
148 
149  DTFuncs_
150  (
151  this->coeffDict_.found("DT")
152  ? this->thermo().composition().species().size()
153  : 0
154  )
155 {}
156 
157 
158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
159 
160 template<class BasicThermophysicalTransportModel>
162 {
163  if
164  (
166  )
167  {
169  const speciesTable& species = composition.species();
170 
171  this->coeffDict_.lookup("mixtureDiffusionCoefficients")
172  >> mixtureDiffusionCoefficients_;
173 
174  if (mixtureDiffusionCoefficients_)
175  {
176  const dictionary& Ddict = this->coeffDict_.subDict("Dm");
177 
178  forAll(species, i)
179  {
180  DmFuncs_.set
181  (
182  i,
183  Function2<scalar>::New(species[i], Ddict).ptr()
184  );
185  }
186  }
187  else
188  {
189  const dictionary& Ddict = this->coeffDict_.subDict("D");
190 
191  // Read the array of specie binary mass diffusion coefficient
192  // functions
193  forAll(species, i)
194  {
195  DFuncs_[i].setSize(species.size());
196 
197  forAll(species, j)
198  {
199  if (j >= i)
200  {
201  const word nameij(species[i] + '-' + species[j]);
202  const word nameji(species[j] + '-' + species[i]);
203 
204  word Dname;
205 
206  if (Ddict.found(nameij) && Ddict.found(nameji))
207  {
208  if (i != j)
209  {
211  << "Binary mass diffusion coefficients "
212  "for Both " << nameij
213  << " and " << nameji << " provided, using "
214  << nameij << endl;
215  }
216 
217  Dname = nameij;
218  }
219  else if (Ddict.found(nameij))
220  {
221  Dname = nameij;
222  }
223  else if (Ddict.found(nameji))
224  {
225  Dname = nameji;
226  }
227  else
228  {
230  << "Binary mass diffusion coefficient for pair "
231  << nameij << " or " << nameji << " not provided"
232  << exit(FatalIOError);
233  }
234 
235  DFuncs_[i].set
236  (
237  j,
238  Function2<scalar>::New(Dname, Ddict).ptr()
239  );
240  }
241  }
242  }
243  }
244 
245  // Optionally read the List of specie Soret thermal diffusion
246  // coefficient functions
247  if (this->coeffDict_.found("DT"))
248  {
249  const dictionary& DTdict = this->coeffDict_.subDict("DT");
250 
251  forAll(species, i)
252  {
253  DTFuncs_.set
254  (
255  i,
256  Function2<scalar>::New(species[i], DTdict).ptr()
257  );
258  }
259  }
260 
261  return true;
262  }
263  else
264  {
265  return false;
266  }
267 }
268 
269 
270 template<class BasicThermophysicalTransportModel>
272 (
273  const volScalarField& Yi
274 ) const
275 {
277 
278  return volScalarField::New
279  (
280  "DEff",
281  this->momentumTransport().rho()
282  *Dm()[composition.index(Yi)]
283  );
284 }
285 
286 
287 template<class BasicThermophysicalTransportModel>
289 (
290  const volScalarField& Yi,
291  const label patchi
292 ) const
293 {
295 
296  return
297  this->momentumTransport().rho().boundaryField()[patchi]
298  *Dm()[composition.index(Yi)].boundaryField()[patchi];
299 }
300 
301 
302 template<class BasicThermophysicalTransportModel>
304 {
306  (
308  (
310  (
311  "q",
312  this->momentumTransport().alphaRhoPhi().group()
313  ),
314  -fvc::interpolate(this->alpha()*this->kappaEff())
315  *fvc::snGrad(this->thermo().T())
316  )
317  );
318 
321 
322  if (Y.size())
323  {
324  surfaceScalarField sumJ
325  (
327  (
328  "sumJ",
329  Y[0].mesh(),
331  )
332  );
333 
334  surfaceScalarField sumJh
335  (
337  (
338  "sumJh",
339  Y[0].mesh(),
341  )
342  );
343 
344  forAll(Y, i)
345  {
346  if (i != composition.defaultSpecie())
347  {
348  const volScalarField hi
349  (
350  composition.Hs(i, this->thermo().p(), this->thermo().T())
351  );
352 
353  const surfaceScalarField ji(this->j(Y[i]));
354  sumJ += ji;
355 
356  sumJh += ji*fvc::interpolate(hi);
357  }
358  }
359 
360  {
361  const label i = composition.defaultSpecie();
362 
363  const volScalarField hi
364  (
365  composition.Hs(i, this->thermo().p(), this->thermo().T())
366  );
367 
368  sumJh -= sumJ*fvc::interpolate(hi);
369  }
370 
371  tmpq.ref() += sumJh;
372  }
373 
374  return tmpq;
375 }
376 
377 
378 template<class BasicThermophysicalTransportModel>
380 (
382 ) const
383 {
384  tmp<fvScalarMatrix> tmpDivq
385  (
386  fvm::Su
387  (
388  -fvc::laplacian(this->alpha()*this->kappaEff(), this->thermo().T()),
389  he
390  )
391  );
392 
395 
396  tmpDivq.ref() -=
397  fvm::laplacianCorrection(this->alpha()*this->alphaEff(), he);
398 
399  surfaceScalarField sumJ
400  (
402  (
403  "sumJ",
404  he.mesh(),
406  )
407  );
408 
409  surfaceScalarField sumJh
410  (
412  (
413  "sumJh",
414  he.mesh(),
415  dimensionedScalar(sumJ.dimensions()*he.dimensions(), 0)
416  )
417  );
418 
419  forAll(Y, i)
420  {
421  if (i != composition.defaultSpecie())
422  {
423  const volScalarField hi
424  (
425  composition.Hs(i, this->thermo().p(), this->thermo().T())
426  );
427 
428  const surfaceScalarField ji(this->j(Y[i]));
429  sumJ += ji;
430 
431  sumJh += ji*fvc::interpolate(hi);
432  }
433  }
434 
435  {
436  const label i = composition.defaultSpecie();
437 
438  const volScalarField hi
439  (
440  composition.Hs(i, this->thermo().p(), this->thermo().T())
441  );
442 
443  sumJh -= sumJ*fvc::interpolate(hi);
444  }
445 
446  tmpDivq.ref() += fvc::div(sumJh*he.mesh().magSf());
447 
448  return tmpDivq;
449 }
450 
451 
452 template<class BasicThermophysicalTransportModel>
454 (
455  const volScalarField& Yi
456 ) const
457 {
458  if (DTFuncs_.size())
459  {
461  const volScalarField& p = this->thermo().p();
462  const volScalarField& T = this->thermo().T();
463 
464  return
465  BasicThermophysicalTransportModel::j(Yi)
467  (
468  evaluate
469  (
470  DTFuncs_[composition.index(Yi)],
472  p,
473  T
474  )
475  )
477  }
478  else
479  {
480  return BasicThermophysicalTransportModel::j(Yi);
481  }
482 }
483 
484 
485 template<class BasicThermophysicalTransportModel>
487 (
488  volScalarField& Yi
489 ) const
490 {
491  if (DTFuncs_.size())
492  {
494  const volScalarField& p = this->thermo().p();
495  const volScalarField& T = this->thermo().T();
496 
497  return
498  BasicThermophysicalTransportModel::divj(Yi)
499  - fvc::div
500  (
502  (
503  evaluate
504  (
505  DTFuncs_[composition.index(Yi)],
507  p,
508  T
509  )
510  )
512  *T.mesh().magSf()
513  );
514  }
515  else
516  {
517  return BasicThermophysicalTransportModel::divj(Yi);
518  }
519 }
520 
521 
522 template<class BasicThermophysicalTransportModel>
524 {
525  BasicThermophysicalTransportModel::predict();
526  updateDm();
527 }
528 
529 
530 template<class BasicThermophysicalTransportModel>
532 {
533  return true;
534 }
535 
536 
537 template<class BasicThermophysicalTransportModel>
539 (
540  const polyTopoChangeMap& map
541 )
542 {
543  // Delete the cached Dm, will be re-created in predict
544  Dm_.clear();
545 }
546 
547 
548 template<class BasicThermophysicalTransportModel>
550 (
551  const polyMeshMap& map
552 )
553 {
554  // Delete the cached Dm, will be re-created in predict
555  Dm_.clear();
556 }
557 
558 
559 template<class BasicThermophysicalTransportModel>
561 (
562  const polyDistributionMap& map
563 )
564 {
565  // Delete the cached Dm, will be re-created in predict
566  Dm_.clear();
567 }
568 
569 
570 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
571 
572 } // End namespace Foam
573 
574 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const dimensionSet & dimensions() const
Return dimensions.
virtual bool movePoints()
Update for mesh motion.
Definition: Fickian.C:531
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
Definition: Fickian.C:272
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Definition: Fickian.C:380
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
Definition: Fickian.C:487
virtual void predict()
Update the diffusion coefficients.
Definition: Fickian.C:523
virtual void mapMesh(const polyMeshMap &map)
Update from another mesh using the given map.
Definition: Fickian.C:550
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:454
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
Definition: Fickian.C:303
virtual void distribute(const polyDistributionMap &map)
Redistribute or update using the given distribution map.
Definition: Fickian.C:561
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:128
virtual void topoChange(const polyTopoChangeMap &map)
Update topology using the given map.
Definition: Fickian.C:539
virtual bool read()
Read thermophysicalTransport dictionary.
Definition: Fickian.C:161
Run-time selectable function of two variables.
Definition: Function2.H:64
static autoPtr< Function2< Type > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function2New.C:32
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
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
void setSize(const label)
Reset size of List.
Definition: List.C:281
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
Specialisation of basicMixture for a mixture consisting of a number for molecular species.
virtual scalar rho(const label speciei, const scalar p, const scalar T) const =0
Density [kg/m^3].
virtual scalar Hs(const label speciei, const scalar p, const scalar T) const =0
Sensible enthalpy [J/kg].
virtual scalar Wi(const label speciei) const =0
Molecular weight of the given specie [kg/kmol].
const speciesTable & species() const
Return the table of species.
label defaultSpecie() const
Return the index of the default specie.
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
label index(const volScalarField &Yi) const
Return the specie index of the given mass-fraction field.
virtual const volScalarField & T() const =0
Temperature [K].
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:998
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
virtual basicSpecieMixture & composition()=0
Return the composition of the multi-component mixture.
virtual volScalarField & p()=0
Pressure [Pa].
A wordList with hashed indices for faster lookup by name.
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:181
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
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
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.
const dimensionSet dimViscosity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimless
const dimensionSet dimTime
scalarList W(const fluidMulticomponentThermo &thermo)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
const dimensionSet dimMass
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
const dimensionSet dimArea
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
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()
basicSpecieMixture & composition
volScalarField & p
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31