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  this->coeffDict_.lookup("mixtureDiffusionCoefficients")
171  >> mixtureDiffusionCoefficients_;
172 
173  if (mixtureDiffusionCoefficients_)
174  {
175  const dictionary& Ddict = this->coeffDict_.subDict("Dm");
176 
177  forAll(species, i)
178  {
179  DmFuncs_.set
180  (
181  i,
183  (
184  species[i],
185  dimPressure,
188  Ddict
189  ).ptr()
190  );
191  }
192  }
193  else
194  {
195  const dictionary& Ddict = this->coeffDict_.subDict("D");
196 
197  // Read the array of specie binary mass diffusion coefficient
198  // functions
199  forAll(species, i)
200  {
201  DFuncs_[i].setSize(species.size());
202 
203  forAll(species, j)
204  {
205  if (j >= i)
206  {
207  const word nameij(species[i] + '-' + species[j]);
208  const word nameji(species[j] + '-' + species[i]);
209 
210  word Dname;
211 
212  if (Ddict.found(nameij) && Ddict.found(nameji))
213  {
214  if (i != j)
215  {
217  << "Binary mass diffusion coefficients "
218  "for Both " << nameij
219  << " and " << nameji << " provided, using "
220  << nameij << endl;
221  }
222 
223  Dname = nameij;
224  }
225  else if (Ddict.found(nameij))
226  {
227  Dname = nameij;
228  }
229  else if (Ddict.found(nameji))
230  {
231  Dname = nameji;
232  }
233  else
234  {
236  << "Binary mass diffusion coefficient for pair "
237  << nameij << " or " << nameji << " not provided"
238  << exit(FatalIOError);
239  }
240 
241  DFuncs_[i].set
242  (
243  j,
245  (
246  Dname,
247  dimPressure,
250  Ddict
251  ).ptr()
252  );
253  }
254  }
255  }
256  }
257 
258  // Optionally read the List of specie Soret thermal diffusion
259  // coefficient functions
260  if (this->coeffDict_.found("DT"))
261  {
262  const dictionary& DTdict = this->coeffDict_.subDict("DT");
263 
264  forAll(species, i)
265  {
266  DTFuncs_.set
267  (
268  i,
270  (
271  species[i],
272  dimPressure,
275  DTdict
276  ).ptr()
277  );
278  }
279  }
280 
281  return true;
282  }
283  else
284  {
285  return false;
286  }
287 }
288 
289 
290 template<class BasicThermophysicalTransportModel>
292 (
293  const volScalarField& Yi
294 ) const
295 {
296  return volScalarField::New
297  (
298  "DEff",
299  this->momentumTransport().rho()
300  *Dm()[this->thermo().specieIndex(Yi)]
301  );
302 }
303 
304 
305 template<class BasicThermophysicalTransportModel>
307 (
308  const volScalarField& Yi,
309  const label patchi
310 ) const
311 {
312  return
313  this->momentumTransport().rho().boundaryField()[patchi]
314  *Dm()[this->thermo().specieIndex(Yi)].boundaryField()[patchi];
315 }
316 
317 
318 template<class BasicThermophysicalTransportModel>
320 {
322  (
324  (
326  (
327  "q",
328  this->momentumTransport().alphaRhoPhi().group()
329  ),
330  -fvc::interpolate(this->alpha()*this->kappaEff())
331  *fvc::snGrad(this->thermo().T())
332  )
333  );
334 
335  const PtrList<volScalarField>& Y = this->thermo().Y();
336  const volScalarField& p = this->thermo().p();
337  const volScalarField& T = this->thermo().T();
338 
339  if (Y.size())
340  {
341  surfaceScalarField sumJ
342  (
344  (
345  "sumJ",
346  Y[0].mesh(),
348  )
349  );
350 
351  surfaceScalarField sumJh
352  (
354  (
355  "sumJh",
356  Y[0].mesh(),
358  )
359  );
360 
361  forAll(Y, i)
362  {
363  if (i != this->thermo().defaultSpecie())
364  {
365  const volScalarField hi(this->thermo().hsi(i, p, T));
366 
367  const surfaceScalarField ji(this->j(Y[i]));
368  sumJ += ji;
369 
370  sumJh += ji*fvc::interpolate(hi);
371  }
372  }
373 
374  {
375  const label i = this->thermo().defaultSpecie();
376 
377  const volScalarField hi(this->thermo().hsi(i, p, T));
378 
379  sumJh -= sumJ*fvc::interpolate(hi);
380  }
381 
382  tmpq.ref() += sumJh;
383  }
384 
385  return tmpq;
386 }
387 
388 
389 template<class BasicThermophysicalTransportModel>
391 (
393 ) const
394 {
395  tmp<fvScalarMatrix> tmpDivq
396  (
397  fvm::Su
398  (
399  -fvc::laplacian(this->alpha()*this->kappaEff(), this->thermo().T()),
400  he
401  )
402  );
403 
404  const PtrList<volScalarField>& Y = this->thermo().Y();
405  const volScalarField& p = this->thermo().p();
406  const volScalarField& T = this->thermo().T();
407 
408  tmpDivq.ref() -=
409  fvm::laplacianCorrection(this->alpha()*this->alphaEff(), he);
410 
411  surfaceScalarField sumJ
412  (
414  (
415  "sumJ",
416  he.mesh(),
418  )
419  );
420 
421  surfaceScalarField sumJh
422  (
424  (
425  "sumJh",
426  he.mesh(),
427  dimensionedScalar(sumJ.dimensions()*he.dimensions(), 0)
428  )
429  );
430 
431  forAll(Y, i)
432  {
433  if (i != this->thermo().defaultSpecie())
434  {
435  const volScalarField hi(this->thermo().hsi(i, p, T));
436 
437  const surfaceScalarField ji(this->j(Y[i]));
438  sumJ += ji;
439 
440  sumJh += ji*fvc::interpolate(hi);
441  }
442  }
443 
444  {
445  const label i = this->thermo().defaultSpecie();
446 
447  const volScalarField hi(this->thermo().hsi(i, p, T));
448 
449  sumJh -= sumJ*fvc::interpolate(hi);
450  }
451 
452  tmpDivq.ref() += fvc::div(sumJh*he.mesh().magSf());
453 
454  return tmpDivq;
455 }
456 
457 
458 template<class BasicThermophysicalTransportModel>
460 (
461  const volScalarField& Yi
462 ) const
463 {
464  if (DTFuncs_.size())
465  {
466  const volScalarField& p = this->thermo().p();
467  const volScalarField& T = this->thermo().T();
468 
469  return
470  BasicThermophysicalTransportModel::j(Yi)
472  (
473  evaluate
474  (
475  DTFuncs_[this->thermo().specieIndex(Yi)],
477  p,
478  T
479  )
480  )
482  }
483  else
484  {
485  return BasicThermophysicalTransportModel::j(Yi);
486  }
487 }
488 
489 
490 template<class BasicThermophysicalTransportModel>
492 (
493  volScalarField& Yi
494 ) const
495 {
496  if (DTFuncs_.size())
497  {
498  const volScalarField& p = this->thermo().p();
499  const volScalarField& T = this->thermo().T();
500 
501  return
502  BasicThermophysicalTransportModel::divj(Yi)
503  - fvc::div
504  (
506  (
507  evaluate
508  (
509  DTFuncs_[this->thermo().specieIndex(Yi)],
511  p,
512  T
513  )
514  )
516  *T.mesh().magSf()
517  );
518  }
519  else
520  {
521  return BasicThermophysicalTransportModel::divj(Yi);
522  }
523 }
524 
525 
526 template<class BasicThermophysicalTransportModel>
528 {
529  BasicThermophysicalTransportModel::predict();
530  updateDm();
531 }
532 
533 
534 template<class BasicThermophysicalTransportModel>
536 {
537  return true;
538 }
539 
540 
541 template<class BasicThermophysicalTransportModel>
543 (
544  const polyTopoChangeMap& map
545 )
546 {
547  // Delete the cached Dm, will be re-created in predict
548  Dm_.clear();
549 }
550 
551 
552 template<class BasicThermophysicalTransportModel>
554 (
555  const polyMeshMap& map
556 )
557 {
558  // Delete the cached Dm, will be re-created in predict
559  Dm_.clear();
560 }
561 
562 
563 template<class BasicThermophysicalTransportModel>
565 (
566  const polyDistributionMap& map
567 )
568 {
569  // Delete the cached Dm, will be re-created in predict
570  Dm_.clear();
571 }
572 
573 
574 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
575 
576 } // End namespace Foam
577 
578 // ************************************************************************* //
#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:535
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
Definition: Fickian.C:292
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Definition: Fickian.C:391
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
Definition: Fickian.C:492
virtual void predict()
Update the diffusion coefficients.
Definition: Fickian.C:527
virtual void mapMesh(const polyMeshMap &map)
Update from another mesh using the given map.
Definition: Fickian.C:554
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:460
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
Definition: Fickian.C:319
virtual void distribute(const polyDistributionMap &map)
Redistribute or update using the given distribution map.
Definition: Fickian.C:565
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:543
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.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, 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
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
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 keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
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: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: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:257
const dimensionSet dimKinematicViscosity
const dimensionSet dimless
const dimensionSet dimTemperature
const dimensionSet dimTime
scalarList W(const fluidMulticomponentThermo &thermo)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
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)
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