31 #include "surfaceInterpolate.H"
41 template<
class BasicThermophysicalTransportModel>
42 void MaxwellStefan<BasicThermophysicalTransportModel>::
43 transformDiffusionCoefficient()
const
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]);
78 A(is, is) -= X[j]*Wm/(DD(i, j)*
W[i]);
83 X[i]*(Wm/(DD(i, j)*
W[j]) - Wm/(DD(i, d)*
W[d]));
85 B(is, js) = X[i]*(Wm/
W[j] - Wm/
W[d]);
108 template<
class BasicThermophysicalTransportModel>
109 void MaxwellStefan<BasicThermophysicalTransportModel>::
110 transformDiffusionCoefficientFields()
const
120 Y[i] = (*YPtrs[i])[
pi];
125 DD(i, j) = (*DijPtrs[i][j])[
pi];
130 transformDiffusionCoefficient();
148 (*DijPtrs[i][j])[
pi] =
D(is, js);
161 template<
class BasicThermophysicalTransportModel>
162 void MaxwellStefan<BasicThermophysicalTransportModel>::transform
164 List<PtrList<volScalarField>>& Dij
167 const PtrList<volScalarField>&
Y = this->
thermo().
Y();
173 YPtrs[i] = &
Y[i].primitiveField();
176 DijPtrs[i][i] = &Dii_[i].primitiveFieldRef();
183 DijPtrs[i][j] = &Dij[i][j].primitiveFieldRef();
190 transformDiffusionCoefficientFields();
197 YPtrs[i] = &
Y[i].boundaryField()[
patchi];
200 DijPtrs[i][i] = &Dii_[i].boundaryFieldRef()[
patchi];
207 DijPtrs[i][j] = &Dij[i][j].boundaryFieldRef()[
patchi];
214 transformDiffusionCoefficientFields();
219 template<
class BasicThermophysicalTransportModel>
229 Dii_.setSize(
Y.size());
230 jexp_.setSize(
Y.size());
251 "jexp" +
Y[i].
name(),
261 Dij[i].setSize(
Y.size());
275 Dij[i].set(j, Dij[j][i].
clone());
293 if (i != d && i != j)
321 template<
class BasicThermophysicalTransportModel>
322 const PtrList<volScalarField>&
323 MaxwellStefan<BasicThermophysicalTransportModel>::Dii()
const
334 template<
class BasicThermophysicalTransportModel>
335 const PtrList<surfaceScalarField>&
336 MaxwellStefan<BasicThermophysicalTransportModel>::jexp()
const
349 template<
class BasicThermophysicalTransportModel>
357 BasicThermophysicalTransportModel
366 DFuncs_(this->
thermo().species().size()),
370 this->coeffDict().
found(
"DT")
371 ? this->
thermo().species().size()
375 W(this->
thermo().species().size()),
391 W[i] = this->
thermo().Wi(i).value();
398 template<
class BasicThermophysicalTransportModel>
408 const dictionary& coeffDict = this->coeffDict();
414 DFuncs_[i].setSize(species.
size());
420 const word nameij(species[i] +
'-' + species[j]);
421 const word nameji(species[j] +
'-' + species[i]);
425 if (Ddict.
found(nameij) && Ddict.
found(nameji))
430 <<
"Binary mass diffusion coefficients "
431 "for both " << nameij <<
" and " << nameji
432 <<
" provided, using " << nameij <<
endl;
437 else if (Ddict.
found(nameij))
441 else if (Ddict.
found(nameji))
448 <<
"Binary mass diffusion coefficients for pair "
449 << nameij <<
" or " << nameji <<
" not provided"
471 if (coeffDict.
found(
"DT"))
501 template<
class BasicThermophysicalTransportModel>
510 this->momentumTransport().
rho()*Dii()[this->
thermo().specieIndex(Yi)]
515 template<
class BasicThermophysicalTransportModel>
523 this->momentumTransport().rho().boundaryField()[
patchi]
528 template<
class BasicThermophysicalTransportModel>
539 this->momentumTransport().alphaRhoPhi().
group()
603 template<
class BasicThermophysicalTransportModel>
658 template<
class BasicThermophysicalTransportModel>
730 template<
class BasicThermophysicalTransportModel>
738 if (this->
thermo().specieIndex(Yi) == d)
749 this->momentumTransport().alphaRhoPhi().
group()
771 BasicThermophysicalTransportModel::j(Yi)
777 template<
class BasicThermophysicalTransportModel>
786 if (this->
thermo().specieIndex(Yi) == d)
809 BasicThermophysicalTransportModel::j(Yi,
patchi)
815 template<
class BasicThermophysicalTransportModel>
822 BasicThermophysicalTransportModel::divj(Yi)
827 template<
class BasicThermophysicalTransportModel>
830 BasicThermophysicalTransportModel::predict();
835 template<
class BasicThermophysicalTransportModel>
842 template<
class BasicThermophysicalTransportModel>
854 template<
class BasicThermophysicalTransportModel>
866 template<
class BasicThermophysicalTransportModel>
#define forAll(list, i)
Loop across all elements in list.
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Run-time selectable function of two variables.
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...
void size(const label)
Override size to be inconsistent with allocated storage.
Base class for multi-component Maxwell Stefan generalised Fick's law diffusion coefficients based tem...
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...
label size() const
Return the number of elements in the UPtrList.
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.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Dimension set for the base types.
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.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
A class for handling words, derived from string.
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.
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.
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.
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)
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const VolField< Type > &)
tmp< fvMatrix< Type > > laplacianCorrection(const VolField< scalar > &gamma, const VolField< Type > &vf)
static const coefficient D("D", dimTemperature, 257.14)
static const coefficient B("B", dimless, 18.678)
static const coefficient A("A", dimPressure, 611.21)
errorManipArg< error, int > exit(error &err, const int errNo=1)
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.
const dimensionSet dimEnergy
const dimensionSet dimPressure
Ostream & endl(Ostream &os)
Add newline and flush stream.
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
dimensionSet transform(const dimensionSet &)
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.
scalarList Y0(nSpecie, 0.0)
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo