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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2014-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::BlendedInterfacialModel
28 
29 Description
30 
31 SourceFiles
32  BlendedInterfacialModel.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef BlendedInterfacialModel_H
37 #define BlendedInterfacialModel_H
38 
39 #include "blendingMethod.H"
40 #include "phasePair.H"
41 #include "orderedPhasePair.H"
42 
43 #include "geometricZeroField.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class BlendedInterfacialModel Declaration
52 \*---------------------------------------------------------------------------*/
53 
54 template<class modelType>
55 class BlendedInterfacialModel
56 {
57  // Private data
58 
59  //- Unordered phase pair
60  const phasePair& pair_;
61 
62  //- Ordered phase pair for dispersed phase 1 in continuous phase 2
63  const orderedPhasePair& pair1In2_;
64 
65  //- Ordered phase pair for dispersed phase 2 in continuous phase 1
66  const orderedPhasePair& pair2In1_;
67 
68  //- Model for region with no obvious dispersed phase
69  autoPtr<modelType> model_;
70 
71  //- Model for dispersed phase 1 in continuous phase 2
72  autoPtr<modelType> model1In2_;
73 
74  //- Model for dispersed phase 2 in continuous phase 1
75  autoPtr<modelType> model2In1_;
76 
77  //- Blending model
78  const blendingMethod& blending_;
79 
80  //- If true set coefficients and forces to 0 at fixed-flux BCs
81  bool correctFixedFluxBCs_;
82 
83 
84  // Private Member Functions
85 
86  //- No copy construct
87  BlendedInterfacialModel
88  (
89  const BlendedInterfacialModel&
90  ) = delete;
91 
92  //- No copy assignment
93  void operator=(const BlendedInterfacialModel<modelType>&) = delete;
94 
95  //- Correct coeff/value on fixed flux boundary conditions
96  template<class GeometricField>
97  void correctFixedFluxBCs(GeometricField& field) const;
98 
99 
100 public:
101 
102  // Constructors
103 
104  //- Construct from the model table, dictionary and pairs
105  BlendedInterfacialModel
106  (
107  const phasePair::dictTable& modelTable,
108  const blendingMethod& blending,
109  const phasePair& pair,
110  const orderedPhasePair& pair1In2,
111  const orderedPhasePair& pair2In1,
112  const bool correctFixedFluxBCs = true
113  );
114 
115 
116  //- Destructor
118 
119 
120  // Member Functions
121 
122  //- Return true if a model is specified for the supplied phase
123  bool hasModel(const phaseModel& phase) const;
124 
125  //- Return the model for the supplied phase
126  const modelType& phaseModel(const phaseModel& phase) const;
127 
128  //- Return the blended force coefficient
129  tmp<volScalarField> K() const;
130 
131  //- Return the face blended force coefficient
132  tmp<surfaceScalarField> Kf() const;
133 
134  //- Return the blended force
135  template<class Type>
136  tmp<GeometricField<Type, fvPatchField, volMesh>> F() const;
137 
138  //- Return the face blended force
139  tmp<surfaceScalarField> Ff() const;
140 
141  //- Return the blended diffusivity
142  tmp<volScalarField> D() const;
143 };
144 
145 
146 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
147 
148 } // End namespace Foam
149 
150 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
151 
152 #ifdef NoRepository
153  #include "BlendedInterfacialModel.C"
154 #endif
155 
156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 
158 #endif
159 
160 // ************************************************************************* //
Foam::BlendedInterfacialModel::phaseModel
const modelType & phaseModel(const phaseModel &phase) const
Return the model for the supplied phase.
Definition: BlendedInterfacialModel.C:463
Foam::BlendedInterfacialModel::Kf
tmp< surfaceScalarField > Kf() const
Return the face blended force coefficient.
Definition: BlendedInterfacialModel.C:341
Foam::BlendedInterfacialModel::~BlendedInterfacialModel
~BlendedInterfacialModel()
Destructor.
Definition: BlendedInterfacialModel.C:290
Foam::BlendedInterfacialModel::F
tmp< GeometricField< Type, fvPatchField, volMesh > > F() const
Return the blended force.
Foam::BlendedInterfacialModel::hasModel
bool hasModel(const phaseModel &phase) const
Return true if a model is specified for the supplied phase.
field
rDeltaTY field()
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::BlendedInterfacialModel::K
tmp< volScalarField > K() const
Return the blended force coefficient.
Definition: BlendedInterfacialModel.C:321
Foam::BlendedInterfacialModel::Ff
tmp< surfaceScalarField > Ff() const
Return the face blended force.
Definition: BlendedInterfacialModel.C:358
Foam::phasePair::dictTable
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Dictionary hash table.
Definition: phasePair.H:61
Foam::BlendedInterfacialModel::D
tmp< volScalarField > D() const
Return the blended diffusivity.
Definition: BlendedInterfacialModel.C:366
geometricZeroField.H