/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2015-2017 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 "phaseSystem.H" #include "surfaceTensionModel.H" #include "aspectRatioModel.H" #include "surfaceInterpolate.H" #include "fvcDdt.H" // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // namespace Foam { defineTypeNameAndDebug(phaseSystem, 0); } const Foam::word Foam::phaseSystem::propertiesName("phaseProperties"); // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * // Foam::tmp Foam::phaseSystem::calcPhi ( const phaseModelList& phaseModels ) const { tmp tmpPhi ( new surfaceScalarField ( "phi", fvc::interpolate(phaseModels[0])*phaseModels[0].phi() ) ); for (label phasei=1; phasei ( new orderedPhasePair ( phaseModels_[key.first()], phaseModels_[key.second()] ) ) ); } // new unordered pair else { phasePairs_.insert ( key, autoPtr ( new phasePair ( phaseModels_[key.first()], phaseModels_[key.second()] ) ) ); } } } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // Foam::phaseSystem::phaseSystem ( const fvMesh& mesh ) : IOdictionary ( IOobject ( "phaseProperties", mesh.time().constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE ) ), mesh_(mesh), phaseModels_(lookup("phases"), phaseModel::iNew(*this)), phi_(calcPhi(phaseModels_)), dpdt_ ( IOobject ( "dpdt", mesh.time().timeName(), mesh ), mesh, dimensionedScalar("dpdt", dimPressure/dimTime, 0) ), MRF_(mesh_) { phi_.writeOpt() = IOobject::AUTO_WRITE; // Blending methods forAllConstIter(dictionary, subDict("blending"), iter) { blendingMethods_.insert ( iter().dict().dictName(), blendingMethod::New ( iter().dict(), phaseModels_.toc() ) ); } // Sub-models generatePairsAndSubModels("surfaceTension", surfaceTensionModels_); generatePairsAndSubModels("aspectRatio", aspectRatioModels_); correctKinematics(); } // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // Foam::phaseSystem::~phaseSystem() {} // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * // Foam::tmp Foam::phaseSystem::rho() const { tmp tmpRho ( phaseModels_[0]*phaseModels_[0].rho() ); for (label phasei=1; phasei Foam::phaseSystem::U() const { tmp tmpU ( phaseModels_[0]*phaseModels_[0].U() ); for (label phasei=1; phasei Foam::phaseSystem::E(const phasePairKey& key) const { if (aspectRatioModels_.found(key)) { return aspectRatioModels_[key]->E(); } else { return tmp ( new volScalarField ( IOobject ( aspectRatioModel::typeName + ":E", this->mesh_.time().timeName(), this->mesh_, IOobject::NO_READ, IOobject::NO_WRITE, false ), this->mesh_, dimensionedScalar("zero", dimless, 1) ) ); } } Foam::tmp Foam::phaseSystem::sigma(const phasePairKey& key) const { if (surfaceTensionModels_.found(key)) { return surfaceTensionModels_[key]->sigma(); } else { return tmp ( new volScalarField ( IOobject ( surfaceTensionModel::typeName + ":sigma", this->mesh_.time().timeName(), this->mesh_, IOobject::NO_READ, IOobject::NO_WRITE, false ), this->mesh_, dimensionedScalar("zero", surfaceTensionModel::dimSigma, 0) ) ); } } void Foam::phaseSystem::solve() {} void Foam::phaseSystem::correct() { forAll(phaseModels_, phasei) { phaseModels_[phasei].correct(); } } void Foam::phaseSystem::correctKinematics() { bool updateDpdt = false; forAll(phaseModels_, phasei) { phaseModels_[phasei].correctKinematics(); updateDpdt = updateDpdt || phaseModels_[phasei].thermo().dpdt(); } // Update the pressure time-derivative if required if (updateDpdt) { dpdt_ = fvc::ddt(phaseModels_.begin()().thermo().p()); } } void Foam::phaseSystem::correctThermo() { forAll(phaseModels_, phasei) { phaseModels_[phasei].correctThermo(); } } void Foam::phaseSystem::correctTurbulence() { forAll(phaseModels_, phasei) { phaseModels_[phasei].correctTurbulence(); } } void Foam::phaseSystem::correctEnergyTransport() { forAll(phaseModels_, phasei) { phaseModels_[phasei].correctEnergyTransport(); } } bool Foam::phaseSystem::read() { if (regIOobject::read()) { bool readOK = true; forAll(phaseModels_, phasei) { readOK &= phaseModels_[phasei].read(); } // models ... return readOK; } else { return false; } } // ************************************************************************* //