/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | Website: https://openfoam.org \\ / A nd | Copyright (C) 2015-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 "InterfaceCompositionPhaseChangePhaseSystem.H" #include "interfaceCompositionModel.H" #include "heatTransferModel.H" #include "diffusiveMassTransferModel.H" #include "fvmSup.H" // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * // template void Foam::InterfaceCompositionPhaseChangePhaseSystem:: correctDmdtfs() { forAllConstIter ( interfaceCompositionModelTable, interfaceCompositionModels_, interfaceCompositionModelIter ) { const sidedInterfaceCompositionModel& compositionModel = interfaceCompositionModelIter(); const phaseInterface& interface = compositionModel.interface(); *dmdtfs_[interface] = Zero; forAllConstIter(phaseInterface, interface, interfaceIter) { const phaseModel& phase = interfaceIter(); if (!compositionModel.haveModelInThe(phase)) continue; forAllConstIter ( hashedWordList, compositionModel.modelInThe(phase).species(), specieIter ) { const word& specie = *specieIter; *dmdtfs_[interface] += (interfaceIter.index() == 0 ? +1 : -1) *( *(*dmidtfSus_[interface])[specie] + *(*dmidtfSps_[interface])[specie]*phase.Y(specie) ); } } } } template Foam::autoPtr Foam::InterfaceCompositionPhaseChangePhaseSystem:: totalDmidtfs() const { autoPtr totalDmidtfsPtr ( new phaseSystem::dmidtfTable ); phaseSystem::dmidtfTable& totalDmidtfs = totalDmidtfsPtr(); forAllConstIter ( interfaceCompositionModelTable, interfaceCompositionModels_, interfaceCompositionModelIter ) { const sidedInterfaceCompositionModel& compositionModel = interfaceCompositionModelIter(); const phaseInterface& interface = compositionModel.interface(); if (!totalDmidtfs.found(interface)) { totalDmidtfs.insert(interface, new HashPtrTable()); } forAllConstIter(phaseInterface, interface, interfaceIter) { const phaseModel& phase = interfaceIter(); if (!compositionModel.haveModelInThe(phase)) continue; forAllConstIter ( hashedWordList, compositionModel.modelInThe(phase).species(), specieIter ) { const word& specie = *specieIter; tmp dmidtf ( (interfaceIter.index() == 0 ? +1 : -1) *( *(*dmidtfSus_[interface])[specie] + *(*dmidtfSps_[interface])[specie]*phase.Y(specie) ) ); if (totalDmidtfs[interface]->found(specie)) { *(*totalDmidtfs[interface])[specie] += dmidtf; } else { totalDmidtfs[interface]->insert(specie, dmidtf.ptr()); } } } } return totalDmidtfsPtr; } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // template Foam::InterfaceCompositionPhaseChangePhaseSystem:: InterfaceCompositionPhaseChangePhaseSystem ( const fvMesh& mesh ) : BasePhaseSystem(mesh), nInterfaceCorrectors_ ( this->template lookupOrDefault