/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | Website: https://openfoam.org \\ / A nd | Copyright (C) 2015-2024 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 "phaseModel.H" #include "phaseSystem.H" #include "diameterModel.H" // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // namespace Foam { defineTypeNameAndDebug(phaseModel, 0); defineRunTimeSelectionTable(phaseModel, phaseSystem); } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // Foam::phaseModel::phaseModel ( const phaseSystem& fluid, const word& phaseName, const bool referencePhase, const label index ) : // volScalarField // ( // referencePhase // ? volScalarField // ( // IOobject // ( // IOobject::groupName("alpha", phaseName), // fluid.mesh().time().name(), // fluid.mesh(), // IOobject::NO_READ, // IOobject::AUTO_WRITE // ), // fluid.mesh(), // dimensionedScalar(dimless, 0) // ) // : volScalarField // ( // IOobject // ( // IOobject::groupName("alpha", phaseName), // fluid.mesh().time().name(), // fluid.mesh(), // IOobject::MUST_READ, // IOobject::AUTO_WRITE // ), // fluid.mesh() // ) // ), // Replaced by the following because the Clang compiler does not use // std::move to transfer the result of the ternary operator volScalarField ( IOobject ( IOobject::groupName("alpha", phaseName), fluid.mesh().time().name(), fluid.mesh(), IOobject::NO_READ, IOobject::AUTO_WRITE ), referencePhase ? volScalarField::New ( IOobject::groupName("alpha", phaseName), fluid.mesh(), dimensionedScalar(dimless, 0) ) : tmp ( new volScalarField ( IOobject ( IOobject::groupName("alpha", phaseName), fluid.mesh().time().name(), fluid.mesh(), IOobject::MUST_READ, IOobject::AUTO_WRITE, false ), fluid.mesh() ) ) ), fluid_(fluid), name_(phaseName), index_(index), residualAlpha_ ( "residualAlpha", dimless, fluid.subDict(phaseName).lookup("residualAlpha") ), alphaMax_(fluid.subDict(phaseName).lookupOrDefault("alphaMax", 1.0)) { diameterModel_ = diameterModel::New(fluid.subDict(phaseName), *this); } Foam::autoPtr Foam::phaseModel::clone() const { NotImplemented; return autoPtr(nullptr); } // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // Foam::phaseModel::~phaseModel() {} // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // const Foam::word& Foam::phaseModel::name() const { return name_; } const Foam::word& Foam::phaseModel::keyword() const { return name_; } Foam::label Foam::phaseModel::index() const { return index_; } const Foam::phaseSystem& Foam::phaseModel::fluid() const { return fluid_; } const Foam::dimensionedScalar& Foam::phaseModel::residualAlpha() const { return residualAlpha_; } Foam::scalar Foam::phaseModel::alphaMax() const { return alphaMax_; } Foam::tmp Foam::phaseModel::d() const { return diameterModel_().d(); } const Foam::diameterModel& Foam::phaseModel::diameter() const { return diameterModel_(); } void Foam::phaseModel::correct() { diameterModel_->correct(); } void Foam::phaseModel::correctContinuityError(const volScalarField& source) {} void Foam::phaseModel::correctKinematics() {} void Foam::phaseModel::correctThermo() {} void Foam::phaseModel::correctReactions() {} void Foam::phaseModel::correctSpecies() {} void Foam::phaseModel::predictMomentumTransport() {} void Foam::phaseModel::predictThermophysicalTransport() {} void Foam::phaseModel::correctMomentumTransport() {} void Foam::phaseModel::correctThermophysicalTransport() {} void Foam::phaseModel::correctUf() {} bool Foam::phaseModel::read() { return diameterModel_->read(fluid_.subDict(name_)); } void Foam::phaseModel::correctInflowOutflow(surfaceScalarField& alphaPhi) const { surfaceScalarField::Boundary& alphaPhiBf = alphaPhi.boundaryFieldRef(); const volScalarField::Boundary& alphaBf = boundaryField(); const surfaceScalarField::Boundary& phiBf = phiRef().boundaryField(); forAll(alphaPhiBf, patchi) { fvsPatchScalarField& alphaPhip = alphaPhiBf[patchi]; if (!alphaPhip.coupled()) { alphaPhip = phiBf[patchi]*alphaBf[patchi]; } } } // ************************************************************************* //