All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
blendingMethod.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-2022 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::blendingMethod
26 
27 Description
28  Abstract base class for functions that are used to combine interfacial
29  sub-models according to the volume fractions of the phases that they apply
30  to.
31 
32 SourceFiles
33  blendingMethod.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef blendingMethod_H
38 #define blendingMethod_H
39 
40 #include "phaseInterface.H"
41 #include "runTimeSelectionTables.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 /*---------------------------------------------------------------------------*\
49  Class blendingMethod Declaration
50 \*---------------------------------------------------------------------------*/
51 
52 class blendingMethod
53 {
54 protected:
55 
56  // Protected Data
57 
58  //- Interface
60 
61 
62  // Protected Static Member Functions
63 
64  //- Read a parameter and check it lies within specified bounds
65  static scalar readParameter
66  (
67  const word& name,
68  const dictionary& dict,
69  const Pair<scalar>& bounds,
70  const bool allowNone
71  );
72 
73  //- Read a parameter for each phase in the interface
75  (
76  const word& name,
77  const dictionary& dict,
78  const phaseInterface& interface,
79  const Pair<scalar>& bounds,
80  const bool allowNone
81  );
82 
83  //- Test if the given scalar is a valid parameter
84  static bool isParameter(const scalar parameter);
85 
86 
87  // Protected Member Functions
88 
89  //- Return a constant field with the given value
91  (
92  const UPtrList<const volScalarField>& alphas,
93  const scalar k
94  ) const;
95 
96  //- Get the volume fraction of the given set
98  (
99  const UPtrList<const volScalarField>& alphas,
100  const label set,
101  const bool protect
102  ) const;
103 
104  //- Get a blending parameter averaged for the given set
106  (
107  const UPtrList<const volScalarField>& alphas,
108  const label set,
109  const Pair<scalar>& parameters
110  ) const;
111 
112  //- Return the coordinate of the blending function
114  (
115  const UPtrList<const volScalarField>& alphas,
116  const label phaseSet,
117  const label systemSet
118  ) const;
119 
120  //- Evaluate the blending function for sets in which all phases can be
121  // continuous
123  (
124  const UPtrList<const volScalarField>& alphas,
125  const label phaseSet,
126  const label systemSet
127  ) const = 0;
128 
129  //- Evaluate the blending function. Filters out phases that cannot
130  // become continuous from the sets, then calls fContinuous
131  virtual tmp<volScalarField> f
132  (
133  const UPtrList<const volScalarField>& alphas,
134  const label phaseSet,
135  const label systemSet
136  ) const;
137 
138 
139 public:
140 
141  //- Runtime type information
142  TypeName("blendingMethod");
143 
144 
145  // Declare runtime construction
147  (
148  autoPtr,
150  dictionary,
151  (
152  const dictionary& dict,
153  const phaseInterface& interface
154  ),
155  (dict, interface)
156  );
157 
158 
159  // Constructors
160 
161  //- Construct from a dictionary and an interface
163  (
164  const dictionary& dict,
165  const phaseInterface& interface
166  );
167 
168 
169  // Selector
170 
172  (
173  const word& modelTypeName,
174  const dictionary& dict,
175  const phaseInterface& interface
176  );
177 
178 
179  //- Destructor
180  virtual ~blendingMethod();
181 
182 
183  // Member Functions
184 
185  //- Return whether or not a phase can be considered continuous
186  virtual bool canBeContinuous(const label index) const = 0;
187 
188  //- Return whether or not this interface can segregate
189  virtual bool canSegregate() const = 0;
190 
191  //- Return the coefficient for models in which phase 1 is dispersed in
192  // phase 2
194  (
195  const UPtrList<const volScalarField>& alphas
196  ) const;
197 
198  //- Return the coefficient for models in which phase 2 is dispersed in
199  // phase 1
201  (
202  const UPtrList<const volScalarField>& alphas
203  ) const;
204 
205  //- Return the coefficient for when the interface is displaced by a
206  // third phase
208  (
209  const UPtrList<const volScalarField>& alphas
210  ) const;
211 };
212 
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
216 } // End namespace Foam
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 #endif
221 
222 // ************************************************************************* //
dictionary dict
static Pair< scalar > readParameters(const word &name, const dictionary &dict, const phaseInterface &interface, const Pair< scalar > &bounds, const bool allowNone)
Read a parameter for each phase in the interface.
tmp< volScalarField > alpha(const UPtrList< const volScalarField > &alphas, const label set, const bool protect) const
Get the volume fraction of the given set.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
static autoPtr< blendingMethod > New(const word &modelTypeName, const dictionary &dict, const phaseInterface &interface)
tmp< volScalarField > parameter(const UPtrList< const volScalarField > &alphas, const label set, const Pair< scalar > &parameters) const
Get a blending parameter averaged for the given set.
tmp< volScalarField > f2DispersedIn1(const UPtrList< const volScalarField > &alphas) const
Return the coefficient for models in which phase 2 is dispersed in.
label k
Boltzmann constant.
virtual tmp< volScalarField > fContinuous(const UPtrList< const volScalarField > &alphas, const label phaseSet, const label systemSet) const =0
Evaluate the blending function for sets in which all phases can be.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
virtual bool canBeContinuous(const label index) const =0
Return whether or not a phase can be considered continuous.
tmp< volScalarField > constant(const UPtrList< const volScalarField > &alphas, const scalar k) const
Return a constant field with the given value.
A class for handling words, derived from string.
Definition: word.H:59
Abstract base class for functions that are used to combine interfacial sub-models according to the vo...
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:54
static scalar readParameter(const word &name, const dictionary &dict, const Pair< scalar > &bounds, const bool allowNone)
Read a parameter and check it lies within specified bounds.
static bool isParameter(const scalar parameter)
Test if the given scalar is a valid parameter.
virtual ~blendingMethod()
Destructor.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual tmp< volScalarField > f(const UPtrList< const volScalarField > &alphas, const label phaseSet, const label systemSet) const
Evaluate the blending function. Filters out phases that cannot.
tmp< volScalarField > f1DispersedIn2(const UPtrList< const volScalarField > &alphas) const
Return the coefficient for models in which phase 1 is dispersed in.
blendingMethod(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
tmp< volScalarField > x(const UPtrList< const volScalarField > &alphas, const label phaseSet, const label systemSet) const
Return the coordinate of the blending function.
TypeName("blendingMethod")
Runtime type information.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
declareRunTimeSelectionTable(autoPtr, blendingMethod, dictionary,(const dictionary &dict, const phaseInterface &interface),(dict, interface))
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual bool canSegregate() const =0
Return whether or not this interface can segregate.
Namespace for OpenFOAM.
const phaseInterface interface_
Interface.
tmp< volScalarField > fDisplaced(const UPtrList< const volScalarField > &alphas) const
Return the coefficient for when the interface is displaced by a.