/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | Website: https://openfoam.org \\ / A nd | Copyright (C) 2018-2023 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 "PhaseTransferPhaseSystem.H" #include "phaseTransferModel.H" #include "fvmSup.H" // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * // template Foam::autoPtr Foam::PhaseTransferPhaseSystem::totalDmdtfs() const { autoPtr totalDmdtfsPtr ( new phaseSystem::dmdtfTable ); phaseSystem::dmdtfTable& totalDmdtfs = totalDmdtfsPtr(); forAllConstIter ( phaseTransferModelTable, phaseTransferModels_, phaseTransferModelIter ) { const phaseInterface& interface = phaseTransferModelIter()->interface(); totalDmdtfs.insert(interface, phaseSystem::dmdtf(interface).ptr()); if (phaseTransferModelIter()->mixture()) { *totalDmdtfs[interface] += *dmdtfs_[interface]; } forAllConstIter ( HashPtrTable, *dmidtfs_[interface], dmidtfIter ) { *totalDmdtfs[interface] += *dmidtfIter(); } } return totalDmdtfsPtr; } // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * // template void Foam::PhaseTransferPhaseSystem::addDmdtYfs ( const phaseSystem::dmdtfTable& dmdtfs, phaseSystem::specieTransferTable& eqns ) const { forAllConstIter(phaseSystem::dmdtfTable, dmdtfs, dmdtfIter) { const phaseInterface interface(*this, dmdtfIter.key()); const volScalarField& dmdtf = *dmdtfIter(); const volScalarField dmdtf12(negPart(dmdtf)); const volScalarField dmdtf21(posPart(dmdtf)); const phaseModel& phase1 = interface.phase1(); const phaseModel& phase2 = interface.phase2(); forAll(phase1.Y(), Yi1) { const volScalarField& Y1 = phase1.Y()[Yi1]; const volScalarField& Y2 = phase2.Y(Y1.member()); *eqns[Y1.name()] += dmdtf21*Y2 + fvm::Sp(dmdtf12, Y1); *eqns[Y2.name()] -= dmdtf12*Y1 + fvm::Sp(dmdtf21, Y2); } } } template void Foam::PhaseTransferPhaseSystem::addDmidtYf ( const phaseSystem::dmidtfTable& dmidtfs, phaseSystem::specieTransferTable& eqns ) const { forAllConstIter(phaseSystem::dmidtfTable, dmidtfs, dmidtfIter) { const phaseInterface interface(*this, dmidtfIter.key()); const phaseModel& phase1 = interface.phase1(); const phaseModel& phase2 = interface.phase2(); forAllConstIter(HashPtrTable, *dmidtfIter(), dmidtfJter) { const word& member = dmidtfJter.key(); const volScalarField& dmidtf = *dmidtfJter(); if (!phase1.pure()) { const volScalarField& Y1 = phase1.Y(member); *eqns[Y1.name()] += dmidtf; } if (!phase2.pure()) { const volScalarField& Y2 = phase2.Y(member); *eqns[Y2.name()] -= dmidtf; } } } } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // template Foam::PhaseTransferPhaseSystem::PhaseTransferPhaseSystem ( const fvMesh& mesh ) : BasePhaseSystem(mesh) { this->generateInterfacialModels(phaseTransferModels_); forAllConstIter ( phaseTransferModelTable, phaseTransferModels_, phaseTransferModelIter ) { const phaseInterface& interface = phaseTransferModelIter()->interface(); this->template validateMassTransfer(interface); if (phaseTransferModelIter()->mixture()) { dmdtfs_.insert ( interface, new volScalarField ( IOobject ( IOobject::groupName ( "phaseTransfer:dmdtf", interface.name() ), this->mesh().time().name(), this->mesh() ), this->mesh(), dimensionedScalar(dimDensity/dimTime, 0) ) ); d2mdtdpfs_.insert ( interface, new volScalarField ( IOobject ( IOobject::groupName ( "phaseTransfer:d2mdtdpf", interface.name() ), this->mesh().time().name(), this->mesh() ), this->mesh(), dimensionedScalar(dimDensity/dimTime/dimPressure, 0) ) ); } dmidtfs_.insert(interface, new HashPtrTable()); const hashedWordList species(phaseTransferModelIter()->species()); forAllConstIter(hashedWordList, species, specieIter) { const word& specie = *specieIter; dmidtfs_[interface]->insert ( specie, new volScalarField ( IOobject ( IOobject::groupName ( IOobject::groupName ( "phaseTransfer:dmidtf", specie ), interface.name() ), this->mesh().time().name(), this->mesh() ), this->mesh(), dimensionedScalar(dimDensity/dimTime, 0) ) ); } } } // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // template Foam::PhaseTransferPhaseSystem::~PhaseTransferPhaseSystem() {} // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * // template Foam::tmp Foam::PhaseTransferPhaseSystem::dmdtf ( const phaseInterfaceKey& key ) const { tmp tDmdtf = BasePhaseSystem::dmdtf(key); if (phaseTransferModels_.found(key)) { if (phaseTransferModels_[key]->mixture()) { tDmdtf.ref() += *dmdtfs_[key]; } forAllConstIter ( HashPtrTable, *dmidtfs_[key], dmidtfIter ) { tDmdtf.ref() += *dmidtfIter(); } } return tDmdtf; } template Foam::PtrList Foam::PhaseTransferPhaseSystem::dmdts() const { PtrList dmdts(BasePhaseSystem::dmdts()); autoPtr totalDmdtfsPtr = this->totalDmdtfs(); const phaseSystem::dmdtfTable& totalDmdtfs = totalDmdtfsPtr(); forAllConstIter(phaseSystem::dmdtfTable, totalDmdtfs, totalDmdtfIter) { const phaseInterface interface(*this, totalDmdtfIter.key()); addField(interface.phase1(), "dmdt", *totalDmdtfIter(), dmdts); addField(interface.phase2(), "dmdt", - *totalDmdtfIter(), dmdts); } return dmdts; } template Foam::PtrList Foam::PhaseTransferPhaseSystem::d2mdtdps() const { PtrList d2mdtdps(BasePhaseSystem::d2mdtdps()); forAllConstIter(phaseSystem::dmdtfTable, d2mdtdpfs_, d2mdtdpfIter) { const phaseInterface interface(*this, d2mdtdpfIter.key()); addField(interface.phase1(), "d2mdtdp", *d2mdtdpfIter(), d2mdtdps); addField(interface.phase2(), "d2mdtdp", - *d2mdtdpfIter(), d2mdtdps); } return d2mdtdps; } template Foam::autoPtr Foam::PhaseTransferPhaseSystem::momentumTransfer() { autoPtr eqnsPtr = BasePhaseSystem::momentumTransfer(); phaseSystem::momentumTransferTable& eqns = eqnsPtr(); this->addDmdtUfs(totalDmdtfs(), eqns); return eqnsPtr; } template Foam::autoPtr Foam::PhaseTransferPhaseSystem::momentumTransferf() { autoPtr eqnsPtr = BasePhaseSystem::momentumTransferf(); phaseSystem::momentumTransferTable& eqns = eqnsPtr(); this->addDmdtUfs(totalDmdtfs(), eqns); return eqnsPtr; } template Foam::autoPtr Foam::PhaseTransferPhaseSystem::heatTransfer() const { autoPtr eqnsPtr = BasePhaseSystem::heatTransfer(); phaseSystem::heatTransferTable& eqns = eqnsPtr(); this->addDmdtHefs(dmdtfs_, eqns); this->addDmidtHefs(dmidtfs_, eqns); return eqnsPtr; } template Foam::autoPtr Foam::PhaseTransferPhaseSystem::specieTransfer() const { autoPtr eqnsPtr ( new phaseSystem::specieTransferTable() ); phaseSystem::specieTransferTable& eqns = eqnsPtr(); // Create a mass transfer matrix for each species of each phase forAll(this->phaseModels_, phasei) { const phaseModel& phase = this->phaseModels_[phasei]; const PtrList& Yi = phase.Y(); forAll(Yi, i) { eqns.insert ( Yi[i].name(), new fvScalarMatrix(Yi[i], dimMass/dimTime) ); } } this->addDmdtYfs(dmdtfs_, eqns); this->addDmidtYf(dmidtfs_, eqns); return eqnsPtr; } template void Foam::PhaseTransferPhaseSystem::correct() { BasePhaseSystem::correct(); // Reset all the mass transfer rates to zero forAllConstIter ( phaseTransferModelTable, phaseTransferModels_, phaseTransferModelIter ) { const phaseInterface& interface = phaseTransferModelIter()->interface(); if (phaseTransferModelIter()->mixture()) { *dmdtfs_[interface] = Zero; *d2mdtdpfs_[interface] = Zero; } const hashedWordList species(phaseTransferModelIter()->species()); forAllConstIter(hashedWordList, species, specieIter) { const word& specie = *specieIter; *(*dmidtfs_[interface])[specie] = Zero; } } // Evaluate the models and sum the results into the mass transfer tables forAllIter ( phaseTransferModelTable, phaseTransferModels_, phaseTransferModelIter ) { const phaseInterface& interface = phaseTransferModelIter()->interface(); if (phaseTransferModelIter()->mixture()) { *dmdtfs_[interface] += phaseTransferModelIter()->dmdtf(); *d2mdtdpfs_[interface] += phaseTransferModelIter()->d2mdtdpf(); } const HashPtrTable dmidtf ( phaseTransferModelIter()->dmidtf() ); forAllConstIter(HashPtrTable, dmidtf, dmidtfIter) { *(*dmidtfs_[interface])[dmidtfIter.key()] += *dmidtfIter(); } } } template bool Foam::PhaseTransferPhaseSystem::read() { if (BasePhaseSystem::read()) { bool readOK = true; // Models ... return readOK; } else { return false; } } // ************************************************************************* //