Skip to content
Snippets Groups Projects
BlendedInterfacialModel.C 10.1 KiB
Newer Older
  • Learn to ignore specific revisions
  • /*---------------------------------------------------------------------------*\
      =========                 |
      \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
       \\    /   O peration     |
    
    OpenFOAM bot's avatar
    OpenFOAM bot committed
        \\  /    A nd           | www.openfoam.com
    
    -------------------------------------------------------------------------------
    
    OpenFOAM bot's avatar
    OpenFOAM bot committed
        Copyright (C) 2014-2016 OpenFOAM Foundation
    
        Copyright (C) 2020 OpenCFD Ltd.
    
    -------------------------------------------------------------------------------
    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 "BlendedInterfacialModel.H"
    #include "fixedValueFvsPatchFields.H"
    
    
    // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
    
    template<class modelType>
    
    void Foam::BlendedInterfacialModel<modelType>::correctFixedFluxBCs
    (
    
        typename GeometricField::Boundary& fieldBf =
    
            field.boundaryFieldRef();
    
    
        forAll(pair_.phase1().phi().boundaryField(), patchi)
    
        {
            if
            (
                isA<fixedValueFvsPatchScalarField>
                (
    
                    pair_.phase1().phi().boundaryField()[patchi]
    
                fieldBf[patchi] = Zero;
    
            }
        }
    }
    
    
    // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
    
    template<class modelType>
    Foam::BlendedInterfacialModel<modelType>::BlendedInterfacialModel
    (
        const phasePair::dictTable& modelTable,
    
        const phasePair& pair,
        const orderedPhasePair& pair1In2,
    
        const orderedPhasePair& pair2In1,
        const bool correctFixedFluxBCs
    
    )
    :
        pair_(pair),
        pair1In2_(pair1In2),
        pair2In1_(pair2In1),
    
        blending_(blending),
        correctFixedFluxBCs_(correctFixedFluxBCs)
    
            (
                modelType::New
                (
                    modelTable[pair_],
                    pair_
                ).ptr()
            );
        }
    
        if (modelTable.found(pair1In2_))
        {
    
            (
                modelType::New
                (
                    modelTable[pair1In2_],
                    pair1In2_
                ).ptr()
            );
        }
    
        if (modelTable.found(pair2In1_))
        {
    
    
    
    // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
    
    template<class modelType>
    Foam::BlendedInterfacialModel<modelType>::~BlendedInterfacialModel()
    {}
    
    
    // * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * * //
    
    template<class modelType>
    Foam::tmp<Foam::volScalarField>
    Foam::BlendedInterfacialModel<modelType>::K() const
    {
    
        if (model_ || model1In2_)
    
            f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
    
        if (model_ || model2In1_)
    
            f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
    
                    pair_.phase1().mesh().time().timeName(),
    
                    pair_.phase1().mesh(),
                    IOobject::NO_READ,
                    IOobject::NO_WRITE,
                    false
    
                dimensionedScalar(modelType::dimK, Zero)
    
         && (model_ || model1In2_ || model2In1_)
    
    template<class modelType>
    Foam::tmp<Foam::surfaceScalarField>
    Foam::BlendedInterfacialModel<modelType>::Kf() const
    {
        tmp<surfaceScalarField> f1, f2;
    
    
        if (model_ || model1In2_)
    
        {
            f1 = fvc::interpolate
            (
                blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed())
            );
        }
    
    
        if (model_ || model2In1_)
    
        {
            f2 = fvc::interpolate
            (
                blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed())
            );
        }
    
        tmp<surfaceScalarField> x
        (
            new surfaceScalarField
            (
                IOobject
                (
                    modelType::typeName + ":Kf",
                    pair_.phase1().mesh().time().timeName(),
                    pair_.phase1().mesh(),
                    IOobject::NO_READ,
                    IOobject::NO_WRITE,
                    false
                ),
                pair_.phase1().mesh(),
    
                dimensionedScalar(modelType::dimK, Zero)
    
         && (model_ || model1In2_ || model2In1_)
    
    template<class modelType>
    template<class Type>
    
    Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh>>
    
    Foam::BlendedInterfacialModel<modelType>::F() const
    {
    
        if (model_ || model1In2_)
    
            f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
    
        if (model_ || model2In1_)
    
            f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
    
        auto x = tmp<GeometricField<Type, fvPatchField, volMesh>>::New
    
                modelType::typeName + ":F",
                pair_.phase1().mesh().time().timeName(),
    
                IOobject::NO_READ,
                IOobject::NO_WRITE,
                false
            ),
            pair_.phase1().mesh(),
            dimensioned<Type>(modelType::dimF, Zero)
    
            x.ref() -= model2In1_->F()*f2; // note : subtraction
    
         && (model_ || model1In2_ || model2In1_)
    
    template<class modelType>
    Foam::tmp<Foam::surfaceScalarField>
    Foam::BlendedInterfacialModel<modelType>::Ff() const
    {
        tmp<surfaceScalarField> f1, f2;
    
    
        if (model_ || model1In2_)
    
        {
            f1 = fvc::interpolate
            (
                blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed())
            );
        }
    
    
        if (model_ || model2In1_)
    
        {
            f2 = fvc::interpolate
            (
                blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed())
            );
        }
    
    
        auto x = tmp<surfaceScalarField>::New
    
                modelType::typeName + ":Ff",
                pair_.phase1().mesh().time().timeName(),
    
                IOobject::NO_READ,
                IOobject::NO_WRITE,
                false
            ),
            pair_.phase1().mesh(),
            dimensionedScalar(modelType::dimF*dimArea, Zero)
    
        x.ref().setOriented();
    
    
            x.ref() -= model2In1_->Ff()*f2; // note : subtraction
    
         && (model_ || model1In2_ || model2In1_)
    
    template<class modelType>
    Foam::tmp<Foam::volScalarField>
    Foam::BlendedInterfacialModel<modelType>::D() const
    {
        tmp<volScalarField> f1, f2;
    
    
        if (model_ || model1In2_)
    
        {
            f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
        }
    
    
        if (model_ || model2In1_)
    
        {
            f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
        }
    
        tmp<volScalarField> x
        (
            new volScalarField
            (
                IOobject
                (
    
                    pair_.phase1().mesh().time().timeName(),
    
                    pair_.phase1().mesh(),
                    IOobject::NO_READ,
                    IOobject::NO_WRITE,
                    false
    
                dimensionedScalar(modelType::dimD, Zero)
    
         && (model_ || model1In2_ || model2In1_)
    
    template<class modelType>
    bool Foam::BlendedInterfacialModel<modelType>::hasModel
    (
        const class phaseModel& phase
    ) const
    {
        return
    
        (
            &phase == &(pair_.phase1())
    
          ? bool(model1In2_)
          : bool(model2In1_)
    
    template<class modelType>
    const modelType& Foam::BlendedInterfacialModel<modelType>::phaseModel
    (
        const class phaseModel& phase
    ) const
    {
    
        return &phase == &(pair_.phase1()) ? *model1In2_ : *model2In1_;
    
    }
    
    
    // ************************************************************************* //