Loading...
Searching...
No Matches
BlendedInterfacialModel.C
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 Copyright (C) 2020-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31#include "surfaceInterpolate.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37
38// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
39
40template<>
44 return f;
45}
46
47
48template<>
51{
52 return fvc::interpolate(f);
53}
54
55// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56
57} // End namespace Foam
58
59// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
60
61template<class ModelType>
62template<class GeoField>
63void Foam::BlendedInterfacialModel<ModelType>::correctFixedFluxBCs
64(
65 GeoField& field
66) const
67{
68 typename GeoField::Boundary& fieldBf = field.boundaryFieldRef();
69
70 forAll(phase1_.phi()().boundaryField(), patchi)
71 {
72 if
73 (
74 isA<fixedValueFvsPatchScalarField>
75 (
76 phase1_.phi()().boundaryField()[patchi]
77 )
78 )
79 {
80 fieldBf[patchi] = Zero;
81 }
82 }
83}
84
85
86template<class ModelType>
87template
89 class Type,
90 template<class> class PatchField,
91 class GeoMesh,
92 class ... Args
93>
95Foam::BlendedInterfacialModel<ModelType>::evaluate
96(
98 (ModelType::*method)(Args ...) const,
99 const word& name,
100 const dimensionSet& dims,
101 const bool subtract,
102 Args ... args
103) const
104{
105 typedef GeometricField<scalar, PatchField, GeoMesh> scalarGeoField;
107
108 tmp<scalarGeoField> f1, f2;
109
110 if (model_ || model1In2_)
111 {
112 f1 =
114 (
115 blending_.f1(phase1_, phase2_)
116 );
117 }
118
119 if (model_ || model2In1_)
120 {
121 f2 =
123 (
124 blending_.f2(phase1_, phase2_)
125 );
126 }
127
128 auto x = typeGeoField::New
129 (
130 IOobject::scopedName(ModelType::typeName, name),
132 phase1_.mesh(),
134 );
135
136 if (model_)
137 {
138 if (subtract)
139 {
141 << "Cannot treat an interfacial model with no distinction "
142 << "between continuous and dispersed phases as signed"
143 << exit(FatalError);
144 }
145
146 x.ref() += (model_().*method)(args ...)*(scalar(1) - f1() - f2());
147 }
148
149 if (model1In2_)
150 {
151 x.ref() += (model1In2_().*method)(args ...)*f1;
152 }
153
154 if (model2In1_)
155 {
156 tmp<typeGeoField> dx = (model2In1_().*method)(args ...)*f2;
157
158 if (subtract)
159 {
160 x.ref() -= dx;
161 }
162 else
163 {
164 x.ref() += dx;
165 }
166 }
167
168 if
169 (
170 correctFixedFluxBCs_
171 && (model_ || model1In2_ || model2In1_)
172 )
173 {
174 correctFixedFluxBCs(x.ref());
175 }
176
177 return x;
178}
179
180
181// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
182
183template<class ModelType>
184Foam::BlendedInterfacialModel<ModelType>::BlendedInterfacialModel
185(
186 const phaseModel& phase1,
187 const phaseModel& phase2,
188 const blendingMethod& blending,
189 autoPtr<ModelType> model,
191 autoPtr<ModelType> model2In1,
192 const bool correctFixedFluxBCs
193)
194:
196 (
198 (
200 phase1.mesh().time().timeName(),
201 phase1.mesh()
202 )
203 ),
204 phase1_(phase1),
205 phase2_(phase2),
206 blending_(blending),
207 model_(model),
208 model1In2_(model1In2),
209 model2In1_(model2In1),
210 correctFixedFluxBCs_(correctFixedFluxBCs)
211{}
212
214template<class ModelType>
215Foam::BlendedInterfacialModel<ModelType>::BlendedInterfacialModel
216(
217 const phasePair::dictTable& modelTable,
218 const blendingMethod& blending,
219 const phasePair& pair,
220 const orderedPhasePair& pair1In2,
221 const orderedPhasePair& pair2In1,
222 const bool correctFixedFluxBCs
224:
226 (
228 (
230 pair.phase1().mesh().time().timeName(),
231 pair.phase1().mesh()
232 )
233 ),
234 phase1_(pair.phase1()),
235 phase2_(pair.phase2()),
236 blending_(blending),
237 correctFixedFluxBCs_(correctFixedFluxBCs)
238{
239 if (modelTable.found(pair))
240 {
241 model_.set
242 (
243 ModelType::New
245 modelTable[pair],
246 pair
247 ).ptr()
248 );
250
251 if (modelTable.found(pair1In2))
252 {
253 model1In2_.set
254 (
255 ModelType::New
256 (
257 modelTable[pair1In2],
258 pair1In2
259 ).ptr()
260 );
261 }
262
263 if (modelTable.found(pair2In1))
264 {
265 model2In1_.set
266 (
267 ModelType::New
268 (
269 modelTable[pair2In1],
270 pair2In1
271 ).ptr()
272 );
274}
275
276
277// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
278
279template<class ModelType>
281{}
282
283
284// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
285
286template<class ModelType>
288(
289 const class phaseModel& phase
290) const
291{
292 return
293 &phase == &(phase1_)
294 ? bool(model1In2_)
295 : bool(model2In1_);
296}
297
298
299template<class ModelType>
301(
302 const class phaseModel& phase
303) const
305 return &phase == &(phase1_) ? model1In2_ : model2In1_;
306}
307
308
309template<class ModelType>
312{
313 tmp<volScalarField> (ModelType::*k)() const = &ModelType::K;
315 return evaluate(k, "K", ModelType::dimK, false);
316}
317
318
319template<class ModelType>
321Foam::BlendedInterfacialModel<ModelType>::K(const scalar residualAlpha) const
322{
323 tmp<volScalarField> (ModelType::*k)(const scalar) const = &ModelType::K;
325 return evaluate(k, "K", ModelType::dimK, false, residualAlpha);
326}
327
328
329template<class ModelType>
332{
333 return evaluate(&ModelType::Kf, "Kf", ModelType::dimK, false);
334}
335
336
337template<class ModelType>
338template<class Type>
342 return evaluate(&ModelType::F, "F", ModelType::dimF, true);
343}
344
345
346template<class ModelType>
350 return evaluate(&ModelType::Ff, "Ff", ModelType::dimF*dimArea, true);
351}
352
353
354template<class ModelType>
358 return evaluate(&ModelType::D, "D", ModelType::dimD, false);
359}
360
361
362template<class ModelType>
365{
366 return evaluate(&ModelType::dmdt, "dmdt", ModelType::dimDmdt, false);
367}
368
369
370template<class ModelType>
372{
373 return os.good();
374}
375
376
377// ************************************************************************* //
label k
if(patchID !=-1)
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.
Foam::tmp< Foam::GeometricField< Type, PatchField, GeoMesh > > evaluate(tmp< GeometricField< Type, PatchField, GeoMesh > >(Foam::dragModel::*method)(Args ...) const, const word &name, const dimensionSet &dims, const bool subtract, Args ... args) const
tmp< volScalarField > K() const
Return the blended force coefficient.
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.
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
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...
Generic dimensioned Type class.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
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()
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
word timeName
Definition getTimeIndex.H:3
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Namespace for OpenFOAM.
void subtract(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
const dimensionSet dimArea(sqr(dimLength))
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
labelList f(nPoints)
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299