/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2017 OpenCFD Ltd. \\/ M anipulation | ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OpenFOAM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OpenFOAM. If not, see . \*---------------------------------------------------------------------------*/ #include "InterfaceCompositionModel.H" #include "phaseModel.H" #include "phasePair.H" #include "pureMixture.H" #include "multiComponentMixture.H" #include "rhoThermo.H" #include "zeroGradientFvPatchFields.H" // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * // template template const typename Foam::multiComponentMixture::thermoType& Foam::InterfaceCompositionModel::getLocalThermo ( const word& speciesName, const multiComponentMixture& globalThermo ) const { return globalThermo.getLocalThermo ( globalThermo.species() [ speciesName ] ); } template template const typename Foam::pureMixture::thermoType& Foam::InterfaceCompositionModel::getLocalThermo ( const word& speciesName, const pureMixture& globalThermo ) const { return globalThermo.cellMixture(0); } template template Foam::tmp Foam::InterfaceCompositionModel::getSpecieMassFraction ( const word& speciesName, const multiComponentMixture& mixture ) const { const fvMesh& mesh = fromThermo_.p().mesh(); tmp tY ( new volScalarField ( IOobject ( "tY", mesh.time().timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE ), mesh, dimensionedScalar("zero", dimless, 0), zeroGradientFvPatchScalarField::typeName ) ); volScalarField& Ys = tY.ref(); Ys = mixture.Y(speciesName); return tY; } template template Foam::tmp Foam::InterfaceCompositionModel::getSpecieMassFraction ( const word& speciesName, const pureMixture& mixture ) const { const fvMesh& mesh = fromThermo_.p().mesh(); tmp tY ( new volScalarField ( IOobject ( "tY", mesh.time().timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE ), mesh, dimensionedScalar("zero", dimless, 1), zeroGradientFvPatchScalarField::typeName ) ); return tY; } template template Foam::tmp Foam::InterfaceCompositionModel::MwMixture ( const pureMixture& mixture ) const { const fvMesh& mesh = fromThermo_.p().mesh(); tmp tM ( new volScalarField ( IOobject ( "tM", mesh.time().timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE ), mesh, dimensionedScalar ( "Mw", dimMass/dimMoles, 1e-3*mixture.cellMixture(0).W() ), zeroGradientFvPatchScalarField::typeName ) ); return tM; } template template Foam::tmp Foam::InterfaceCompositionModel::MwMixture ( const multiComponentMixture& mixture ) const { return refCast(mixture).W(); } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // template Foam::InterfaceCompositionModel::InterfaceCompositionModel ( const dictionary& dict, const phasePair& pair ) : interfaceCompositionModel(dict, pair), fromThermo_ ( pair.from().mesh().lookupObject ( IOobject::groupName ( basicThermo::dictName, pair.from().name() ) ) ), toThermo_ ( pair.to().mesh().lookupObject ( IOobject::groupName ( basicThermo::dictName, pair.to().name() ) ) ), Le_("Le", dimless, dict.lookupOrDefault("Le", 1.0)) {} // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * // template Foam::tmp Foam::InterfaceCompositionModel::D ( const word& speciesName ) const { const typename Thermo::thermoType& fromThermo = getLocalThermo ( speciesName, fromThermo_ ); const volScalarField& p(fromThermo_.p()); const volScalarField& T(fromThermo_.T()); tmp tmpD ( new volScalarField ( IOobject ( IOobject::groupName("D", pair_.name()), p.time().timeName(), p.mesh() ), p.mesh(), dimensionedScalar("zero", dimArea/dimTime, 0) ) ); volScalarField& D = tmpD.ref(); forAll(p, cellI) { D[cellI] = fromThermo.alphah(p[cellI], T[cellI]) /fromThermo.rho(p[cellI], T[cellI]); } D /= Le_; D.correctBoundaryConditions(); return tmpD; } template Foam::tmp Foam::InterfaceCompositionModel::L ( const word& speciesName, const volScalarField& Tf ) const { const typename Thermo::thermoType& fromThermo = getLocalThermo(speciesName, fromThermo_); const typename OtherThermo::thermoType& toThermo = getLocalThermo(speciesName, toThermo_); const volScalarField& p(fromThermo_.p()); tmp tmpL ( new volScalarField ( IOobject ( IOobject::groupName("L", pair_.name()), p.time().timeName(), p.mesh() ), p.mesh(), dimensionedScalar("zero", dimEnergy/dimMass, 0), zeroGradientFvPatchScalarField::typeName ) ); volScalarField& L = tmpL.ref(); // from Thermo (from) to Thermo (to) forAll(p, cellI) { L[cellI] = fromThermo.Hc() - toThermo.Hc(); } L.correctBoundaryConditions(); return tmpL; } template Foam::tmp Foam::InterfaceCompositionModel::dY ( const word& speciesName, const volScalarField& Tf ) const { NotImplemented; return tmp(); } template Foam::tmp Foam::InterfaceCompositionModel::Yf ( const word& speciesName, const volScalarField& Tf ) const { NotImplemented; return tmp(); } // ************************************************************************* //