BlendedInterfacialModel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2014-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::BlendedInterfacialModel
26 
27 Description
28 
29 SourceFiles
30  BlendedInterfacialModel.C
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef BlendedInterfacialModel_H
35 #define BlendedInterfacialModel_H
36 
37 #include "blendingMethod.H"
38 #include "phasePair.H"
39 #include "orderedPhasePair.H"
40 
41 #include "geometricZeroField.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 /*---------------------------------------------------------------------------*\
49  Class BlendedInterfacialModel Declaration
50 \*---------------------------------------------------------------------------*/
51 
52 template<class ModelType>
54 {
55  // Private data
56 
57  //- Reference to phase 1
58  const phaseModel& phase1_;
59 
60  //- Reference to phase 2
61  const phaseModel& phase2_;
62 
63  //- Blending model
65 
66  //- Model for region with no obvious dispersed phase
68 
69  //- Model for dispersed phase 1 in continuous phase 2
71 
72  //- Model for dispersed phase 2 in continuous phase 1
74 
75  //- If true set coefficients and forces to 0 at fixed-flux BCs
77 
78 
79  // Private Member Functions
80 
81  //- Disallow default bitwise copy construct
83 
84  //- Disallow default bitwise assignment
86 
87  //- Correct coeff/value on fixed flux boundary conditions
88  template<class GeometricField>
89  void correctFixedFluxBCs(GeometricField& field) const;
90 
91 
92 public:
93 
94  // Constructors
95 
96  //- Construct from two phases, blending method and three models
98  (
99  const phaseModel& phase1,
100  const phaseModel& phase2,
101  const blendingMethod& blending,
103  autoPtr<ModelType> model1In2,
104  autoPtr<ModelType> model2In1,
105  const bool correctFixedFluxBCs = true
106  );
107 
108 
109  //- Construct from the model table, dictionary and pairs
111  (
112  const phasePair::dictTable& modelTable,
113  const blendingMethod& blending,
114  const phasePair& pair,
115  const orderedPhasePair& pair1In2,
116  const orderedPhasePair& pair2In1,
117  const bool correctFixedFluxBCs = true
118  );
119 
120 
121  //- Destructor
123 
124 
125  // Member Functions
126 
127  //- Return true if a model is specified for the supplied phase
128  bool hasModel(const phaseModel& phase) const;
129 
130  //- Return the model for the supplied phase
131  const ModelType& model(const phaseModel& phase) const;
132 
133  //- Return the sign of the explicit value for the supplied phase
134  scalar sign(const phaseModel& phase) const;
135 
136  //- Return the blended force coefficient
137  tmp<volScalarField> K() const;
138 
139  //- Return the blended force coefficient with a specified residual alpha
140  tmp<volScalarField> K(const scalar residualAlpha) const;
141 
142  //- Return the face blended force coefficient
143  tmp<surfaceScalarField> Kf() const;
144 
145  //- Return the blended force
146  template<class Type>
148 
149  //- Return the face blended force
150  tmp<surfaceScalarField> Ff() const;
151 
152  //- Return the blended diffusivity
153  tmp<volScalarField> D() const;
154 };
155 
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 } // End namespace Foam
160 
161 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162 
163 #ifdef NoRepository
164  #include "BlendedInterfacialModel.C"
165 #endif
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 #endif
170 
171 // ************************************************************************* //
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Foam::BlendedInterfacialModel::model_
autoPtr< ModelType > model_
Model for region with no obvious dispersed phase.
Definition: BlendedInterfacialModel.H:66
Foam::phasePair
Definition: phasePair.H:49
Foam::BlendedInterfacialModel::sign
scalar sign(const phaseModel &phase) const
Return the sign of the explicit value for the supplied phase.
Foam::BlendedInterfacialModel::BlendedInterfacialModel
BlendedInterfacialModel(const BlendedInterfacialModel< ModelType > &)
Disallow default bitwise copy construct.
Foam::phase
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::BlendedInterfacialModel::Ff
tmp< surfaceScalarField > Ff() const
Return the face blended force.
Foam::BlendedInterfacialModel::~BlendedInterfacialModel
~BlendedInterfacialModel()
Destructor.
Foam::BlendedInterfacialModel::phase1_
const phaseModel & phase1_
Reference to phase 1.
Definition: BlendedInterfacialModel.H:57
Foam::BlendedInterfacialModel::F
tmp< GeometricField< Type, fvPatchField, volMesh > > F() const
Return the blended force.
Foam::BlendedInterfacialModel::correctFixedFluxBCs_
bool correctFixedFluxBCs_
If true set coefficients and forces to 0 at fixed-flux BCs.
Definition: BlendedInterfacialModel.H:75
Foam::BlendedInterfacialModel::hasModel
bool hasModel(const phaseModel &phase) const
Return true if a model is specified for the supplied phase.
Foam::BlendedInterfacialModel::blending_
const blendingMethod & blending_
Blending model.
Definition: BlendedInterfacialModel.H:63
Foam::BlendedInterfacialModel::model1In2_
autoPtr< ModelType > model1In2_
Model for dispersed phase 1 in continuous phase 2.
Definition: BlendedInterfacialModel.H:69
Foam::BlendedInterfacialModel::phase2_
const phaseModel & phase2_
Reference to phase 2.
Definition: BlendedInterfacialModel.H:60
Foam::orderedPhasePair
Definition: orderedPhasePair.H:47
Foam::BlendedInterfacialModel::model
const ModelType & model(const phaseModel &phase) const
Return the model for the supplied phase.
Foam::BlendedInterfacialModel::D
tmp< volScalarField > D() const
Return the blended diffusivity.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::HashTable
An STL-conforming hash table.
Definition: HashTable.H:61
Foam::blendingMethod
Definition: blendingMethod.H:49
Foam::autoPtr< ModelType >
Foam::BlendedInterfacialModel::operator=
void operator=(const BlendedInterfacialModel< ModelType > &)
Disallow default bitwise assignment.
Foam::BlendedInterfacialModel
Definition: BlendedInterfacialModel.H:52
Foam::BlendedInterfacialModel::model2In1_
autoPtr< ModelType > model2In1_
Model for dispersed phase 2 in continuous phase 1.
Definition: BlendedInterfacialModel.H:72
Foam::BlendedInterfacialModel::Kf
tmp< surfaceScalarField > Kf() const
Return the face blended force coefficient.
phase2
phaseModel & phase2
Definition: createFields.H:13
Foam::BlendedInterfacialModel::K
tmp< volScalarField > K() const
Return the blended force coefficient.
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
geometricZeroField.H
phase1
phaseModel & phase1
Definition: createFields.H:12
Foam::BlendedInterfacialModel::correctFixedFluxBCs
void correctFixedFluxBCs(GeometricField &field) const
Correct coeff/value on fixed flux boundary conditions.