/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2015 OpenFOAM Foundation \\/ 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" // * * * * * * * * * * * * * 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); } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // template Foam::InterfaceCompositionModel::InterfaceCompositionModel ( const dictionary& dict, const phasePair& pair ) : interfaceCompositionModel(dict, pair), thermo_ ( pair.phase1().mesh().lookupObject ( IOobject::groupName(basicThermo::dictName, pair.phase1().name()) ) ), otherThermo_ ( pair.phase2().mesh().lookupObject ( IOobject::groupName(basicThermo::dictName, pair.phase2().name()) ) ), Le_("Le", dimless, dict) {} // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // template Foam::InterfaceCompositionModel:: ~InterfaceCompositionModel() {} // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * // template Foam::tmp Foam::InterfaceCompositionModel::dY ( const word& speciesName, const volScalarField& Tf ) const { return Yf(speciesName, Tf) - thermo_.composition().Y() [ thermo_.composition().species()[speciesName] ]; } template Foam::tmp Foam::InterfaceCompositionModel::D ( const word& speciesName ) const { const typename Thermo::thermoType& localThermo = getLocalThermo ( speciesName, thermo_ ); const volScalarField& p(thermo_.p()); const volScalarField& T(thermo_.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()); forAll(p, cellI) { D[cellI] = localThermo.alphah(p[cellI], T[cellI]) /localThermo.rho(p[cellI], T[cellI]); } D /= Le_; return tmpD; } template Foam::tmp Foam::InterfaceCompositionModel::L ( const word& speciesName, const volScalarField& Tf ) const { const typename Thermo::thermoType& localThermo = getLocalThermo ( speciesName, thermo_ ); const typename OtherThermo::thermoType& otherLocalThermo = getLocalThermo ( speciesName, otherThermo_ ); const volScalarField& p(thermo_.p()); const volScalarField& otherP(otherThermo_.p()); tmp tmpL ( new volScalarField ( IOobject ( IOobject::groupName("L", pair_.name()), p.time().timeName(), p.mesh() ), p.mesh(), dimensionedScalar("zero", dimEnergy/dimMass, 0) ) ); volScalarField& L(tmpL()); forAll(p, cellI) { L[cellI] = localThermo.Ha(p[cellI], Tf[cellI]) - otherLocalThermo.Ha(otherP[cellI], Tf[cellI]); } return tmpL; } template void Foam::InterfaceCompositionModel::addMDotL ( const volScalarField& K, const volScalarField& Tf, volScalarField& mDotL, volScalarField& mDotLPrime ) const { forAllConstIter(hashedWordList, this->speciesNames_, iter) { volScalarField rhoKDL ( thermo_.rhoThermo::rho() *K *D(*iter) *L(*iter, Tf) ); mDotL += rhoKDL*dY(*iter, Tf); mDotLPrime += rhoKDL*YfPrime(*iter, Tf); } } // ************************************************************************* //