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-2023 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  Struct blendingParameter Declaration
50 \*---------------------------------------------------------------------------*/
51 
52 struct blendingParameter
53 {
54  //- Is this parameter valid?
55  bool valid;
56 
57  //- The parameter value. Could be templated if needed.
58  scalar value;
59 };
60 
61 
62 /*---------------------------------------------------------------------------*\
63  Class blendingMethod Declaration
64 \*---------------------------------------------------------------------------*/
65 
66 class blendingMethod
67 {
68 protected:
69 
70  // Protected Data
71 
72  //- Interface
74 
75 
76  // Protected Static Member Functions
77 
78  //- Read a parameter and check it lies within specified bounds
80  (
81  const word& name,
82  const dictionary& dict,
83  const Pair<scalar>& bounds,
84  const bool allowNone
85  );
86 
87  //- Read a parameter for each phase in the interface
89  (
90  const word& name,
91  const dictionary& dict,
92  const phaseInterface& interface,
93  const Pair<scalar>& bounds,
94  const bool allowNone
95  );
96 
97 
98  // Protected Member Functions
99 
100  //- Return a constant field with the given value
102  (
103  const UPtrList<const volScalarField>& alphas,
104  const scalar k
105  ) const;
106 
107  //- Get the volume fraction of the given set
109  (
110  const UPtrList<const volScalarField>& alphas,
111  const label set,
112  const bool protect
113  ) const;
114 
115  //- Get a blending parameter averaged for the given set
117  (
118  const UPtrList<const volScalarField>& alphas,
119  const label set,
120  const Pair<blendingParameter>& parameters
121  ) const;
122 
123  //- Return the coordinate of the blending function
125  (
126  const UPtrList<const volScalarField>& alphas,
127  const label phaseSet,
128  const label systemSet
129  ) const;
130 
131  //- Evaluate the blending function for sets in which all phases can be
132  // continuous
134  (
135  const UPtrList<const volScalarField>& alphas,
136  const label phaseSet,
137  const label systemSet
138  ) const = 0;
139 
140  //- Evaluate the blending function. Filters out phases that cannot
141  // become continuous from the sets, then calls fContinuous
142  virtual tmp<volScalarField> f
143  (
144  const UPtrList<const volScalarField>& alphas,
145  const label phaseSet,
146  const label systemSet
147  ) const;
148 
149 
150 public:
151 
152  //- Runtime type information
153  TypeName("blendingMethod");
154 
155 
156  // Declare runtime construction
158  (
159  autoPtr,
161  dictionary,
162  (
163  const dictionary& dict,
164  const phaseInterface& interface
165  ),
166  (dict, interface)
167  );
168 
169 
170  // Constructors
171 
172  //- Construct from a dictionary and an interface
174  (
175  const dictionary& dict,
176  const phaseInterface& interface
177  );
178 
179 
180  // Selector
181 
183  (
184  const word& modelTypeName,
185  const dictionary& dict,
186  const phaseInterface& interface
187  );
188 
189 
190  //- Destructor
191  virtual ~blendingMethod();
192 
193 
194  // Member Functions
195 
196  //- Return whether or not a phase can be considered continuous
197  virtual bool canBeContinuous(const label index) const = 0;
198 
199  //- Return whether or not this interface can segregate
200  virtual bool canSegregate() const = 0;
201 
202  //- Return the coefficient for models in which phase 1 is dispersed in
203  // phase 2
205  (
206  const UPtrList<const volScalarField>& alphas
207  ) const;
208 
209  //- Return the coefficient for models in which phase 2 is dispersed in
210  // phase 1
212  (
213  const UPtrList<const volScalarField>& alphas
214  ) const;
215 
216  //- Return the coefficient for when the interface is displaced by a
217  // third phase
219  (
220  const UPtrList<const volScalarField>& alphas
221  ) const;
222 };
223 
224 
225 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 
227 } // End namespace Foam
228 
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 
231 #endif
232 
233 // ************************************************************************* //
label k
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:66
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Abstract base class for functions that are used to combine interfacial sub-models according to the vo...
virtual ~blendingMethod()
Destructor.
tmp< volScalarField > x(const UPtrList< const volScalarField > &alphas, const label phaseSet, const label systemSet) const
Return the coordinate of the blending function.
tmp< volScalarField > f1DispersedIn2(const UPtrList< const volScalarField > &alphas) const
Return the coefficient for models in which phase 1 is dispersed in.
virtual bool canBeContinuous(const label index) const =0
Return whether or not a phase can be considered continuous.
TypeName("blendingMethod")
Runtime type information.
tmp< volScalarField > parameter(const UPtrList< const volScalarField > &alphas, const label set, const Pair< blendingParameter > &parameters) const
Get a blending parameter averaged for the given set.
blendingMethod(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
static autoPtr< blendingMethod > New(const word &modelTypeName, const dictionary &dict, const phaseInterface &interface)
static blendingParameter readParameter(const word &name, const dictionary &dict, const Pair< scalar > &bounds, const bool allowNone)
Read a parameter and check it lies within specified bounds.
tmp< volScalarField > f2DispersedIn1(const UPtrList< const volScalarField > &alphas) const
Return the coefficient for models in which phase 2 is dispersed in.
const phaseInterface interface_
Interface.
tmp< volScalarField > alpha(const UPtrList< const volScalarField > &alphas, const label set, const bool protect) const
Get the volume fraction of the given set.
virtual bool canSegregate() const =0
Return whether or not this interface can segregate.
tmp< volScalarField > constant(const UPtrList< const volScalarField > &alphas, const scalar k) const
Return a constant field with the given value.
declareRunTimeSelectionTable(autoPtr, blendingMethod, dictionary,(const dictionary &dict, const phaseInterface &interface),(dict, interface))
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.
static Pair< blendingParameter > 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.
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.
tmp< volScalarField > fDisplaced(const UPtrList< const volScalarField > &alphas) const
Return the coefficient for when the interface is displaced by a.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Class to represent an interface between phases. Derivations can further specify the configuration of ...
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Macros to ease declaration of run-time selection tables.
dictionary dict
bool valid
Is this parameter valid?
scalar value
The parameter value. Could be templated if needed.