1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171
|
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2014-2016 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/>.
Class
Foam::BlendedInterfacialModel
Description
SourceFiles
BlendedInterfacialModel.C
\*---------------------------------------------------------------------------*/
#ifndef BlendedInterfacialModel_H
#define BlendedInterfacialModel_H
#include "blendingMethod.H"
#include "phasePair.H"
#include "orderedPhasePair.H"
#include "geometricZeroField.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class BlendedInterfacialModel Declaration
\*---------------------------------------------------------------------------*/
template<class ModelType>
class BlendedInterfacialModel
{
// Private data
//- Reference to phase 1
const phaseModel& phase1_;
//- Reference to phase 2
const phaseModel& phase2_;
//- Blending model
const blendingMethod& blending_;
//- Model for region with no obvious dispersed phase
autoPtr<ModelType> model_;
//- Model for dispersed phase 1 in continuous phase 2
autoPtr<ModelType> model1In2_;
//- Model for dispersed phase 2 in continuous phase 1
autoPtr<ModelType> model2In1_;
//- If true set coefficients and forces to 0 at fixed-flux BCs
bool correctFixedFluxBCs_;
// Private Member Functions
//- Disallow default bitwise copy construct
BlendedInterfacialModel(const BlendedInterfacialModel<ModelType>&);
//- Disallow default bitwise assignment
void operator=(const BlendedInterfacialModel<ModelType>&);
//- Correct coeff/value on fixed flux boundary conditions
template<class GeometricField>
void correctFixedFluxBCs(GeometricField& field) const;
public:
// Constructors
//- Construct from two phases, blending method and three models
BlendedInterfacialModel
(
const phaseModel& phase1,
const phaseModel& phase2,
const blendingMethod& blending,
autoPtr<ModelType> model,
autoPtr<ModelType> model1In2,
autoPtr<ModelType> model2In1,
const bool correctFixedFluxBCs = true
);
//- Construct from the model table, dictionary and pairs
BlendedInterfacialModel
(
const phasePair::dictTable& modelTable,
const blendingMethod& blending,
const phasePair& pair,
const orderedPhasePair& pair1In2,
const orderedPhasePair& pair2In1,
const bool correctFixedFluxBCs = true
);
//- Destructor
~BlendedInterfacialModel();
// Member Functions
//- Return true if a model is specified for the supplied phase
bool hasModel(const phaseModel& phase) const;
//- Return the model for the supplied phase
const ModelType& model(const phaseModel& phase) const;
//- Return the sign of the explicit value for the supplied phase
scalar sign(const phaseModel& phase) const;
//- Return the blended force coefficient
tmp<volScalarField> K() const;
//- Return the blended force coefficient with a specified residual alpha
tmp<volScalarField> K(const scalar residualAlpha) const;
//- Return the face blended force coefficient
tmp<surfaceScalarField> Kf() const;
//- Return the blended force
template<class Type>
tmp<GeometricField<Type, fvPatchField, volMesh>> F() const;
//- Return the face blended force
tmp<surfaceScalarField> Ff() const;
//- Return the blended diffusivity
tmp<volScalarField> D() const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "BlendedInterfacialModel.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
|