Files
openfoam/applications/solvers/multiphase/twoPhaseEulerFoam/twoPhaseSystem/phaseModel/phaseModel.C
Mark Olesen 660f3e5492 ENH: cleanup autoPtr class (issue #639)
Improve alignment of its behaviour with std::unique_ptr

  - element_type typedef
  - release() method - identical to ptr() method
  - get() method to get the pointer without checking and without releasing it.
  - operator*() for dereferencing

Method name changes

  - renamed rawPtr() to get()
  - renamed rawRef() to ref(), removed unused const version.

Removed methods/operators

  - assignment from a raw pointer was deleted (was rarely used).
    Can be convenient, but uncontrolled and potentially unsafe.
    Do allow assignment from a literal nullptr though, since this
    can never leak (and also corresponds to the unique_ptr API).

Additional methods

  - clone() method: forwards to the clone() method of the underlying
    data object with argument forwarding.

  - reset(autoPtr&&) as an alternative to operator=(autoPtr&&)

STYLE: avoid implicit conversion from autoPtr to object type in many places

- existing implementation has the following:

     operator const T&() const { return operator*(); }

  which means that the following code works:

       autoPtr<mapPolyMesh> map = ...;
       updateMesh(*map);    // OK: explicit dereferencing
       updateMesh(map());   // OK: explicit dereferencing
       updateMesh(map);     // OK: implicit dereferencing

  for clarity it may preferable to avoid the implicit dereferencing

- prefer operator* to operator() when deferenced a return value
  so it is clearer that a pointer is involve and not a function call
  etc    Eg,   return *meshPtr_;  vs.  return meshPtr_();
2018-02-26 12:00:00 +01:00

273 lines
6.8 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "phaseModel.H"
#include "twoPhaseSystem.H"
#include "diameterModel.H"
#include "fvMatrix.H"
#include "PhaseCompressibleTurbulenceModel.H"
#include "dragModel.H"
#include "heatTransferModel.H"
#include "fixedValueFvsPatchFields.H"
#include "fixedValueFvPatchFields.H"
#include "slipFvPatchFields.H"
#include "partialSlipFvPatchFields.H"
#include "fvcFlux.H"
#include "surfaceInterpolate.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::phaseModel::phaseModel
(
const twoPhaseSystem& fluid,
const dictionary& phaseProperties,
const word& phaseName
)
:
volScalarField
(
IOobject
(
IOobject::groupName("alpha", phaseName),
fluid.mesh().time().timeName(),
fluid.mesh(),
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
fluid.mesh(),
dimensionedScalar("alpha", dimless, 0)
),
fluid_(fluid),
name_(phaseName),
phaseDict_
(
phaseProperties.subDict(name_)
),
residualAlpha_
(
"residualAlpha",
dimless,
fluid.subDict(phaseName).lookup("residualAlpha")
),
alphaMax_(phaseDict_.lookupOrDefault("alphaMax", 1.0)),
thermo_(rhoThermo::New(fluid.mesh(), name_)),
U_
(
IOobject
(
IOobject::groupName("U", name_),
fluid.mesh().time().timeName(),
fluid.mesh(),
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
fluid.mesh()
),
alphaPhi_
(
IOobject
(
IOobject::groupName("alphaPhi", name_),
fluid.mesh().time().timeName(),
fluid.mesh()
),
fluid.mesh(),
dimensionedScalar("0", dimensionSet(0, 3, -1, 0, 0), 0)
),
alphaRhoPhi_
(
IOobject
(
IOobject::groupName("alphaRhoPhi", name_),
fluid.mesh().time().timeName(),
fluid.mesh()
),
fluid.mesh(),
dimensionedScalar("0", dimensionSet(1, 0, -1, 0, 0), 0)
)
{
alphaPhi_.setOriented();
alphaRhoPhi_.setOriented();
thermo_->validate("phaseModel " + name_, "h", "e");
const word phiName = IOobject::groupName("phi", name_);
IOobject phiHeader
(
phiName,
fluid_.mesh().time().timeName(),
fluid_.mesh(),
IOobject::NO_READ
);
if (phiHeader.typeHeaderOk<surfaceScalarField>(true))
{
Info<< "Reading face flux field " << phiName << endl;
phiPtr_.reset
(
new surfaceScalarField
(
IOobject
(
phiName,
fluid_.mesh().time().timeName(),
fluid_.mesh(),
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
fluid_.mesh()
)
);
}
else
{
Info<< "Calculating face flux field " << phiName << endl;
wordList phiTypes
(
U_.boundaryField().size(),
calculatedFvPatchScalarField::typeName
);
forAll(U_.boundaryField(), i)
{
if
(
isA<fixedValueFvPatchVectorField>(U_.boundaryField()[i])
|| isA<slipFvPatchVectorField>(U_.boundaryField()[i])
|| isA<partialSlipFvPatchVectorField>(U_.boundaryField()[i])
)
{
phiTypes[i] = fixedValueFvsPatchScalarField::typeName;
}
}
phiPtr_.reset
(
new surfaceScalarField
(
IOobject
(
phiName,
fluid_.mesh().time().timeName(),
fluid_.mesh(),
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
fvc::flux(U_),
phiTypes
)
);
}
dPtr_ = diameterModel::New
(
phaseDict_,
*this
);
turbulence_ =
PhaseCompressibleTurbulenceModel<phaseModel>::New
(
*this,
thermo_->rho(),
U_,
alphaRhoPhi_,
phi(),
*this
);
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::phaseModel::~phaseModel()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
const Foam::phaseModel& Foam::phaseModel::otherPhase() const
{
return fluid_.otherPhase(*this);
}
Foam::tmp<Foam::volScalarField> Foam::phaseModel::d() const
{
return dPtr_().d();
}
Foam::PhaseCompressibleTurbulenceModel<Foam::phaseModel>&
Foam::phaseModel::turbulence()
{
return *turbulence_;
}
const Foam::PhaseCompressibleTurbulenceModel<Foam::phaseModel>&
Foam::phaseModel::turbulence() const
{
return *turbulence_;
}
void Foam::phaseModel::correct()
{
return dPtr_->correct();
}
bool Foam::phaseModel::read(const dictionary& phaseProperties)
{
phaseDict_ = phaseProperties.subDict(name_);
return dPtr_->read(phaseDict_);
}
void Foam::phaseModel::correctInflowOutflow(surfaceScalarField& alphaPhi) const
{
surfaceScalarField::Boundary& alphaPhiBf = alphaPhi.boundaryFieldRef();
const volScalarField::Boundary& alphaBf = boundaryField();
const surfaceScalarField::Boundary& phiBf = phi().boundaryField();
forAll(alphaPhiBf, patchi)
{
fvsPatchScalarField& alphaPhip = alphaPhiBf[patchi];
if (!alphaPhip.coupled())
{
alphaPhip = phiBf[patchi]*alphaBf[patchi];
}
}
}
// ************************************************************************* //