Loading...
Searching...
No Matches
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-2018 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::BlendedInterfacialModel
28
29Description
30
31SourceFiles
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#include "geometricZeroField.H"
43
44// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45
46namespace Foam
47{
49/*---------------------------------------------------------------------------*\
50 Class blendedInterfacialModel Declaration
51\*---------------------------------------------------------------------------*/
52
54{
55 public:
56
57 //- Convenience function to interpolate blending values. Needs to be
58 // specialised, so can't sit in the templated class.
59 template<class GeoField>
61};
62
63
64/*---------------------------------------------------------------------------*\
65 Class BlendedInterfacialModel Declaration
66\*---------------------------------------------------------------------------*/
68template<class ModelType>
69class BlendedInterfacialModel
70:
71 public regIOobject
72{
73 // Private data
74
75 //- Reference to phase 1
76 const phaseModel& phase1_;
77
78 //- Reference to phase 2
79 const phaseModel& phase2_;
80
81 //- Blending model
82 const blendingMethod& blending_;
83
84 //- Model for region with no obvious dispersed phase
85 autoPtr<ModelType> model_;
86
87 //- Model for dispersed phase 1 in continuous phase 2
88 autoPtr<ModelType> model1In2_;
89
90 //- Model for dispersed phase 2 in continuous phase 1
91 autoPtr<ModelType> model2In1_;
92
93 //- If true set coefficients and forces to 0 at fixed-flux BCs
94 bool correctFixedFluxBCs_;
95
96
97 // Private Member Functions
98
99 //- Disallow default bitwise copy construct
100 BlendedInterfacialModel(const BlendedInterfacialModel<ModelType>&);
101
102 //- Disallow default bitwise assignment
103 void operator=(const BlendedInterfacialModel<ModelType>&);
104
105 //- Correct coeff/value on fixed flux boundary conditions
106 template<class GeoField>
107 void correctFixedFluxBCs(GeoField& field) const;
108
109 //- Return the blended coeff/value
110 template
111 <
112 class Type,
113 template<class> class PatchField,
114 class GeoMesh,
115 class ... Args
116 >
118 (
120 (ModelType::*method)(Args ...) const,
121 const word& name,
122 const dimensionSet& dims,
123 const bool subtract,
124 Args ... args
125 ) const;
126
127
128public:
129
130 //- Runtime type information
131 TypeName("BlendedInterfacialModel");
132
133
134 // Constructors
135
136 //- Construct from two phases, blending method and three models
137 BlendedInterfacialModel
138 (
139 const phaseModel& phase1,
140 const phaseModel& phase2,
141 const blendingMethod& blending,
143 autoPtr<ModelType> model1In2,
144 autoPtr<ModelType> model2In1,
145 const bool correctFixedFluxBCs = true
146 );
147
148
149 //- Construct from the model table, dictionary and pairs
150 BlendedInterfacialModel
151 (
152 const phasePair::dictTable& modelTable,
153 const blendingMethod& blending,
154 const phasePair& pair,
155 const orderedPhasePair& pair1In2,
156 const orderedPhasePair& pair2In1,
157 const bool correctFixedFluxBCs = true
158 );
159
160
161 //- Destructor
163
164
165 // Member Functions
166
167 //- Return true if a model is specified for the supplied phase
168 bool hasModel(const phaseModel& phase) const;
169
170 //- Return the model for the supplied phase
171 const ModelType& model(const phaseModel& phase) const;
172
173 //- Return the sign of the explicit value for the supplied phase
174 scalar sign(const phaseModel& phase) const;
175
176 //- Return the blended force coefficient
177 tmp<volScalarField> K() const;
178
179 //- Return the blended force coefficient with a specified residual alpha
180 tmp<volScalarField> K(const scalar residualAlpha) const;
181
182 //- Return the face blended force coefficient
184
185 //- Return the blended force
186 template<class Type>
188
189 //- Return the face blended force
191
192 //- Return the blended diffusivity
193 tmp<volScalarField> D() const;
194
195 //- Return the blended mass transfer rate
197
198 //- Dummy write for regIOobject
199 bool writeData(Ostream& os) const;
200};
201
202
203// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204
205#define defineBlendedInterfacialModelTypeNameAndDebug(ModelType, DebugSwitch) \
206 \
207 defineTemplateTypeNameAndDebugWithName \
208 ( \
209 BlendedInterfacialModel<ModelType>, \
210 ( \
211 word(BlendedInterfacialModel<ModelType>::typeName_()) + "<" \
212 + ModelType::typeName_() + ">" \
213 ).c_str(), \
214 DebugSwitch \
215 );
216
217// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218
219} // End namespace Foam
220
221// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222
223#ifdef NoRepository
225#endif
226
227// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228
229#endif
230
231// ************************************************************************* //
phaseModel & phase1
phaseModel & phase2
tmp< surfaceScalarField > Kf() const
Return the face blended force coefficient.
tmp< surfaceScalarField > Ff() const
Return the face blended force.
tmp< volScalarField > D() const
Return the blended diffusivity.
bool hasModel(const phaseModel &phase) const
Return true if a model is specified for the supplied phase.
tmp< GeometricField< Type, fvPatchField, volMesh > > F() const
Return the blended force.
bool writeData(Ostream &os) const
Dummy write for regIOobject.
scalar sign(const phaseModel &phase) const
Return the sign of the explicit value for the supplied phase.
tmp< volScalarField > K() const
Return the blended force coefficient.
TypeName("BlendedInterfacialModel")
Runtime type information.
const ModelType & model(const phaseModel &phase) const
Return the model for the supplied phase.
const modelType & phaseModel(const phaseModel &phase) const
Return the model for the supplied phase.
tmp< volScalarField > dmdt() const
Return the blended mass transfer rate.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition GeoMesh.H:46
Generic GeometricField class.
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static tmp< GeoField > interpolate(tmp< volScalarField > f)
Convenience function to interpolate blending values. Needs to be.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition phasePair.H:52
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Dictionary hash table.
Definition phasePair.H:59
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
regIOobject(const IOobject &io, const bool isTimeObject=false)
Construct from IOobject. The optional flag adds special handling if the object is the top-level regIO...
Definition regIOobject.C:43
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
rDeltaTY field()
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
void subtract(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
labelList f(nPoints)
Foam::argList args(argc, argv)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68