blendingMethod.C
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 \*---------------------------------------------------------------------------*/
25 
26 #include "blendingMethod.H"
27 #include "phaseSystem.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
35 }
36 
37 
38 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
39 
41 (
42  const word& name,
43  const dictionary& dict,
44  const Pair<scalar>& bounds,
45  const bool allowNone
46 )
47 {
48  if (dict.found(name) || allowNone)
49  {
50  const token t(dict.lookup(name));
51 
52  if (allowNone && t.isWord() && t.wordToken() == "none")
53  {
54  return {false, NaN};
55  }
56 
57  if (t.isNumber())
58  {
59  forAll(bounds, i)
60  {
61  const label s = i == 0 ? -1 : +1;
62 
63  if (s*t.number() > s*bounds[i])
64  {
66  << "Blending parameter " << name << " is "
67  << (i == 0 ? "less" : "greater") << " than "
68  << bounds[i] << exit(FatalError);
69  }
70  }
71 
72  return {true, t.number()};
73  }
74 
76  << "wrong token type - expected Scalar or the word 'none', found "
77  << t.info() << exit(FatalIOError);
78  }
79 
80  return {false, NaN};
81 }
82 
83 
85 (
86  const word& name,
87  const dictionary& dict,
88  const phaseInterface& interface,
89  const Pair<scalar>& bounds,
90  const bool allowNone
91 )
92 {
93  const word name1 = IOobject::groupName(name, interface.phase1().name());
94  const word name2 = IOobject::groupName(name, interface.phase2().name());
95  return
97  (
98  readParameter(name1, dict, bounds, allowNone),
99  readParameter(name2, dict, bounds, allowNone)
100  );
101 }
102 
103 
104 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
105 
107 (
108  const UPtrList<const volScalarField>& alphas,
109  const scalar k
110 ) const
111 {
112  return
114  (
115  name(k),
116  alphas.first().mesh(),
118  );
119 }
120 
121 
123 (
124  const UPtrList<const volScalarField>& alphas,
125  const label set,
126  const bool protect
127 ) const
128 {
129  tmp<volScalarField> talpha = constant(alphas, 0);
130 
131  forAllConstIter(phaseInterface, interface_, iter)
132  {
133  if (0b01 << iter.index() & set)
134  {
135  talpha.ref() +=
136  protect
137  ? max(iter().residualAlpha(), alphas[iter().index()])()
138  : alphas[iter().index()];
139  }
140  }
141 
142  return talpha;
143 }
144 
145 
147 (
148  const UPtrList<const volScalarField>& alphas,
149  const label set,
150  const Pair<blendingParameter>& parameters
151 ) const
152 {
153  tmp<volScalarField> talphaParameter = constant(alphas, 0);
154 
155  forAllConstIter(phaseInterface, interface_, iter)
156  {
157  if (0b01 << iter.index() & set)
158  {
159  talphaParameter.ref() +=
160  max(iter().residualAlpha(), alphas[iter().index()])
161  *parameters[iter.index()].value;
162  }
163  }
164 
165  return talphaParameter/alpha(alphas, set, true);
166 }
167 
168 
170 (
171  const UPtrList<const volScalarField>& alphas,
172  const label phaseSet,
173  const label systemSet
174 ) const
175 {
176  return
177  systemSet == 0b00
178  ? alpha(alphas, phaseSet, false)
179  : alpha(alphas, phaseSet, false)/alpha(alphas, systemSet, true);
180 }
181 
182 
184 (
185  const UPtrList<const volScalarField>& alphas,
186  const label phaseSet,
187  const label systemSet
188 ) const
189 {
190  label canBeContinuousPhaseSet = 0b00;
191  label canBeContinuousSystemSet = 0b00;
192  forAllConstIter(phaseInterface, interface_, iter)
193  {
194  if (canBeContinuous(iter.index()))
195  {
196  canBeContinuousPhaseSet += 0b01 << iter.index() & phaseSet;
197  canBeContinuousSystemSet += 0b01 << iter.index() & systemSet;
198  }
199  }
200 
201  if (canBeContinuousPhaseSet == 0)
202  {
203  return constant(alphas, 0);
204  }
205 
206  if (canBeContinuousPhaseSet == canBeContinuousSystemSet)
207  {
208  return constant(alphas, 1);
209  }
210 
211  return
212  fContinuous
213  (
214  alphas,
215  canBeContinuousPhaseSet,
216  canBeContinuousSystemSet
217  );
218 }
219 
220 
221 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
222 
224 (
225  const dictionary& dict,
226  const phaseInterface& interface
227 )
228 :
229  interface_(interface)
230 {}
231 
232 
233 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
234 
236 {}
237 
238 
239 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
240 
242 (
243  const UPtrList<const volScalarField>& alphas
244 ) const
245 {
246  return f(alphas, 0b10, 0b11);
247 }
248 
249 
251 (
252  const UPtrList<const volScalarField>& alphas
253 ) const
254 {
255  return f(alphas, 0b01, 0b11);
256 }
257 
258 
260 (
261  const UPtrList<const volScalarField>& alphas
262 ) const
263 {
264  return 1 - f(alphas, 0b11, 0b00);
265 }
266 
267 
268 // ************************************************************************* //
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
static word groupName(Name name, const word &group)
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
T & first()
Return reference to the first element of the list.
Definition: UPtrListI.H:43
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.
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 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.
tmp< volScalarField > alpha(const UPtrList< const volScalarField > &alphas, const label set, const bool protect) const
Get the volume fraction of the given set.
tmp< volScalarField > constant(const UPtrList< const volScalarField > &alphas, const scalar k) const
Return a constant field with the given value.
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.
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 ...
const phaseModel & phase1() const
Return phase 1.
const phaseModel & phase2() const
Return phase 2.
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:109
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A token holds items read from Istream.
Definition: token.H:73
bool isNumber() const
Definition: tokenI.H:498
InfoProxy< token > info() const
Return info proxy.
Definition: token.H:391
bool isWord() const
Definition: tokenI.H:261
const word & wordToken() const
Definition: tokenI.H:266
scalar number() const
Definition: tokenI.H:503
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static Type NaN()
Return a primitive with all components set to NaN.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
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
const dimensionSet dimless
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
error FatalError
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
labelList f(nPoints)
dictionary dict