BlendedInterfacialModel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2014-2019 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 #include "geometricZeroField.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 /*---------------------------------------------------------------------------*\
48  Class blendedInterfacialModel Declaration
49 \*---------------------------------------------------------------------------*/
50 
52 {
53  public:
54 
55  //- Convenience function to interpolate blending values. Needs to be
56  // specialised, so can't sit in the templated class.
57  template<class GeoField>
59 };
60 
61 
62 /*---------------------------------------------------------------------------*\
63  Class BlendedInterfacialModel Declaration
64 \*---------------------------------------------------------------------------*/
65 
66 template<class ModelType>
68 :
69  public regIOobject
70 {
71  // Private Data
72 
73  //- Reference to phase 1
74  const phaseModel& phase1_;
75 
76  //- Reference to phase 2
77  const phaseModel& phase2_;
78 
79  //- Blending model
80  const blendingMethod& blending_;
81 
82  //- Model for region with no obvious dispersed phase
83  autoPtr<ModelType> model_;
84 
85  //- Model for dispersed phase 1 in continuous phase 2
86  autoPtr<ModelType> model1In2_;
87 
88  //- Model for dispersed phase 2 in continuous phase 1
89  autoPtr<ModelType> model2In1_;
90 
91  //- If true set coefficients and forces to 0 at fixed-flux BCs
92  bool correctFixedFluxBCs_;
93 
94 
95  // Private Member Functions
96 
97  //- Correct coeff/value on fixed flux boundary conditions
98  template<class GeoField>
99  void correctFixedFluxBCs(GeoField& field) const;
100 
101  //- Return the blended coeff/value
102  template
103  <
104  class Type,
105  template<class> class PatchField,
106  class GeoMesh,
107  class ... Args
108  >
110  (
112  (ModelType::*method)(Args ...) const,
113  const word& name,
114  const dimensionSet& dims,
115  const bool subtract,
116  Args ... args
117  ) const;
118 
119 
120 public:
121 
122  //- Runtime type information
123  TypeName("BlendedInterfacialModel");
124 
125 
126  // Constructors
127 
128  //- Construct from two phases, blending method and three models
130  (
131  const phaseModel& phase1,
132  const phaseModel& phase2,
133  const blendingMethod& blending,
134  autoPtr<ModelType> model,
135  autoPtr<ModelType> model1In2,
136  autoPtr<ModelType> model2In1,
137  const bool correctFixedFluxBCs = true
138  );
139 
140 
141  //- Construct from the model table, dictionary and pairs
143  (
144  const phasePair::dictTable& modelTable,
145  const blendingMethod& blending,
146  const phasePair& pair,
147  const orderedPhasePair& pair1In2,
148  const orderedPhasePair& pair2In1,
149  const bool correctFixedFluxBCs = true
150  );
151 
152  //- Disallow default bitwise copy construction
154  (
156  ) = delete;
157 
158 
159  //- Destructor
161 
162 
163  // Member Functions
164 
165  //- Return true if a model is specified for the supplied phase
166  bool hasModel(const phaseModel& phase) const;
167 
168  //- Return the model for the supplied phase
169  const ModelType& model(const phaseModel& phase) const;
170 
171  //- Return the sign of the explicit value for the supplied phase
172  scalar sign(const phaseModel& phase) const;
173 
174  //- Return the blended force coefficient
175  tmp<volScalarField> K() const;
176 
177  //- Return the blended force coefficient with a specified residual alpha
178  tmp<volScalarField> K(const scalar residualAlpha) const;
179 
180  //- Return the face blended force coefficient
181  tmp<surfaceScalarField> Kf() const;
182 
183  //- Return the blended force
184  template<class Type>
186 
187  //- Return the face blended force
188  tmp<surfaceScalarField> Ff() const;
189 
190  //- Return the blended diffusivity
191  tmp<volScalarField> D() const;
192 
193  //- Return the blended mass transfer rate
194  tmp<volScalarField> dmdt() const;
195 
196  //- Dummy write for regIOobject
197  bool writeData(Ostream& os) const;
198 
199 
200  // Member Operators
201 
202  //- Disallow default bitwise assignment
203  void operator=(const BlendedInterfacialModel<ModelType>&) = delete;
204 };
205 
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 #define defineBlendedInterfacialModelTypeNameAndDebug(ModelType, DebugSwitch) \
210  \
211  defineTemplateTypeNameAndDebugWithName \
212  ( \
213  BlendedInterfacialModel<ModelType>, \
214  ( \
215  word(BlendedInterfacialModel<ModelType>::typeName_()) + "<" \
216  + ModelType::typeName_() + ">" \
217  ).c_str(), \
218  DebugSwitch \
219  );
220 
221 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 
223 } // End namespace Foam
224 
225 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 
227 #ifdef NoRepository
228  #include "BlendedInterfacialModel.C"
229 #endif
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 #endif
234 
235 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
Generic GeometricField class.
CGAL::Exact_predicates_exact_constructions_kernel K
void subtract(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Dimension set for the base types.
Definition: dimensionSet.H:120
A class for handling words, derived from string.
Definition: word.H:59
static tmp< GeoField > interpolate(tmp< volScalarField > f)
Convenience function to interpolate blending values. Needs to be.
volVectorField F(fluid.F())
An STL-conforming hash table.
Definition: HashTable.H:61
phaseModel & phase1
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:70
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const bool writeData(readBool(pdfDictionary.lookup("writeData")))
labelList f(nPoints)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:52
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
phaseModel & phase2
A class for managing temporary objects.
Definition: PtrList.H:53
surfaceScalarField Ff(fluid.Ff())
Foam::argList args(argc, argv)
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Namespace for OpenFOAM.