Loading...
Searching...
No Matches
multiphaseSystem.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) 2013-2019 OpenFOAM Foundation
9 Copyright (C) 2020-2022 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
29#include "multiphaseSystem.H"
30#include "alphaContactAngleFvPatchScalarField.H"
31
32#include "MULES.H"
33#include "subCycle.H"
34#include "UniformField.H"
35
36#include "fvcDdt.H"
37#include "fvcDiv.H"
38#include "fvcSnGrad.H"
39#include "fvcFlux.H"
40#include "fvcMeshPhi.H"
41#include "fvcSup.H"
42
43#include "fvmDdt.H"
44#include "fvmLaplacian.H"
45#include "fvmSup.H"
48
49// * * * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * //
50
51namespace Foam
52{
55}
56
57
58// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
59
60void Foam::multiphaseSystem::calcAlphas()
61{
62 scalar level = 0.0;
63 alphas_ == 0.0;
64
65 forAll(phases(), i)
66 {
67 alphas_ += level*phases()[i];
68 level += 1.0;
69 }
70}
71
72
73void Foam::multiphaseSystem::solveAlphas()
74{
76 {
77 phases()[phasei].correctBoundaryConditions();
78 }
79
80 // Calculate the void fraction
81 volScalarField alphaVoid
82 (
84 (
85 "alphaVoid",
86 mesh_.time().timeName(),
87 mesh_
88 ),
89 mesh_,
91 );
92 forAll(stationaryPhases(), stationaryPhasei)
93 {
94 alphaVoid -= stationaryPhases()[stationaryPhasei];
95 }
96
97 // Generate face-alphas
100 {
103 (
104 phasei,
106 (
107 IOobject::groupName("alphaf", phase.name()),
108 upwind<scalar>(mesh_, phi_).interpolate(phase)
109 )
110 );
111 }
112
113 // Create correction fluxes
114 PtrList<surfaceScalarField> alphaPhiCorrs(phases().size());
115 forAll(stationaryPhases(), stationaryPhasei)
116 {
117 phaseModel& phase = stationaryPhases()[stationaryPhasei];
118
119 alphaPhiCorrs.set
120 (
121 phase.index(),
123 (
124 IOobject::groupName("alphaPhiCorr", phase.name()),
125 - upwind<scalar>(mesh_, phi_).flux(phase)
126 )
127 );
128 }
129 forAll(movingPhases(), movingPhasei)
130 {
131 phaseModel& phase = movingPhases()[movingPhasei];
133
134 alphaPhiCorrs.set
135 (
136 phase.index(),
138 (
139 IOobject::groupName("alphaPhiCorr", phase.name()),
140 fvc::flux(phi_, phase, "div(phi," + alpha.name() + ')')
141 )
142 );
143
144 surfaceScalarField& alphaPhiCorr = alphaPhiCorrs[phase.index()];
145
147 {
150
151 if (&phase2 == &phase) continue;
152
154
155 cAlphaTable::const_iterator cAlpha
156 (
157 cAlphas_.find(phasePairKey(phase.name(), phase2.name()))
158 );
159
160 if (cAlpha != cAlphas_.end())
161 {
163 (
164 (mag(phi_) + mag(phir))/mesh_.magSf()
165 );
166
167 phir += min(cAlpha()*phic, max(phic))*nHatf(phase, phase2);
168 }
169
170 word phirScheme
171 (
172 "div(phir," + alpha2.name() + ',' + alpha.name() + ')'
173 );
174
175 alphaPhiCorr += fvc::flux
176 (
177 -fvc::flux(-phir, phase2, phirScheme),
178 phase,
179 phirScheme
180 );
181 }
182
183 phase.correctInflowOutflow(alphaPhiCorr);
184
186 (
188 phase,
189 phi_,
190 alphaPhiCorr,
191 zeroField(),
192 zeroField(),
193 min(alphaVoid.primitiveField(), phase.alphaMax())(),
194 zeroField(),
195 true
196 );
197 }
198
199 // Limit the flux sums, fixing those of the stationary phases
200 labelHashSet fixedAlphaPhiCorrs;
201 forAll(stationaryPhases(), stationaryPhasei)
202 {
203 fixedAlphaPhiCorrs.insert(stationaryPhases()[stationaryPhasei].index());
204 }
205 MULES::limitSum(alphafs, alphaPhiCorrs, fixedAlphaPhiCorrs);
206
207 // Solve for the moving phase alphas
208 forAll(movingPhases(), movingPhasei)
209 {
210 phaseModel& phase = movingPhases()[movingPhasei];
212
213 surfaceScalarField& alphaPhi = alphaPhiCorrs[phase.index()];
214 alphaPhi += upwind<scalar>(mesh_, phi_).flux(phase);
215 phase.correctInflowOutflow(alphaPhi);
216
218 (
220 (
221 "Sp",
222 mesh_.time().timeName(),
223 mesh_
224 ),
225 mesh_,
227 );
228
230 (
231 "Su",
232 min(alpha, scalar(1))*fvc::div(fvc::absolute(phi_, phase.U()))
233 );
234
235 if (phase.divU().valid())
236 {
237 const scalarField& dgdt = phase.divU()();
238
239 forAll(dgdt, celli)
240 {
241 if (dgdt[celli] > 0.0)
242 {
243 Sp[celli] -= dgdt[celli];
244 Su[celli] += dgdt[celli];
245 }
246 else if (dgdt[celli] < 0.0)
247 {
248 Sp[celli] +=
249 dgdt[celli]
250 *(1 - alpha[celli])/max(alpha[celli], 1e-4);
251 }
252 }
253 }
254
255 forAll(phases(), phasej)
256 {
257 const phaseModel& phase2 = phases()[phasej];
259
260 if (&phase2 == &phase) continue;
261
262 if (phase2.divU().valid())
263 {
264 const scalarField& dgdt2 = phase2.divU()();
265
266 forAll(dgdt2, celli)
267 {
268 if (dgdt2[celli] < 0.0)
269 {
270 Sp[celli] +=
271 dgdt2[celli]
272 *(1 - alpha2[celli])/max(alpha2[celli], 1e-4);
273
274 Su[celli] -=
275 dgdt2[celli]
276 *alpha[celli]/max(alpha2[celli], 1e-4);
277 }
278 else if (dgdt2[celli] > 0.0)
279 {
280 Sp[celli] -= dgdt2[celli];
281 }
282 }
283 }
284 }
285
287 (
289 alpha,
290 alphaPhi,
291 Sp,
292 Su
293 );
294
295 phase.alphaPhiRef() = alphaPhi;
296 }
297
298 // Report the phase fractions and the phase fraction sum
300 {
302
303 Info<< phase.name() << " fraction, min, max = "
304 << phase.weightedAverage(mesh_.V()).value()
305 << ' ' << min(phase).value()
306 << ' ' << max(phase).value()
307 << endl;
308 }
309
310 volScalarField sumAlphaMoving
311 (
313 (
314 "sumAlphaMoving",
315 mesh_.time().timeName(),
316 mesh_
317 ),
318 mesh_,
320 );
321 forAll(movingPhases(), movingPhasei)
322 {
323 sumAlphaMoving += movingPhases()[movingPhasei];
324 }
325
326 Info<< "Phase-sum volume fraction, min, max = "
327 << (sumAlphaMoving + 1 - alphaVoid)().weightedAverage(mesh_.V()).value()
328 << ' ' << min(sumAlphaMoving + 1 - alphaVoid).value()
329 << ' ' << max(sumAlphaMoving + 1 - alphaVoid).value()
330 << endl;
331
332 // Correct the sum of the phase fractions to avoid drift
333 forAll(movingPhases(), movingPhasei)
334 {
335 movingPhases()[movingPhasei] *= alphaVoid/sumAlphaMoving;
336 }
337}
338
339
340Foam::tmp<Foam::surfaceVectorField> Foam::multiphaseSystem::nHatfv
341(
342 const volScalarField& alpha1,
344) const
345{
346 /*
347 // Cell gradient of alpha
348 volVectorField gradAlpha =
349 alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2);
350
351 // Interpolated face-gradient of alpha
352 surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
353 */
354
355 surfaceVectorField gradAlphaf
356 (
359 );
360
361 // Face unit interface normal
362 return gradAlphaf/(mag(gradAlphaf) + deltaN_);
363}
364
365
366Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseSystem::nHatf
367(
368 const volScalarField& alpha1,
370) const
371{
372 // Face unit interface normal flux
373 return nHatfv(alpha1, alpha2) & mesh_.Sf();
374}
375
376
377// Correction for the boundary condition on the unit normal nHat on
378// walls to produce the correct contact angle.
379
380// The dynamic contact angle is calculated from the component of the
381// velocity on the direction of the interface, parallel to the wall.
382
383void Foam::multiphaseSystem::correctContactAngle
384(
385 const phaseModel& phase1,
386 const phaseModel& phase2,
388) const
389{
390 const volScalarField::Boundary& gbf
392
393 const fvBoundaryMesh& boundary = mesh_.boundary();
394
395 forAll(boundary, patchi)
396 {
397 if
398 (
400 (
401 gbf[patchi]
402 )
403 )
404 {
405 const auto& acap =
406 refCast
407 <
408 const reactingMultiphaseEuler::
409 alphaContactAngleFvPatchScalarField
410 >(gbf[patchi]);
411
412 vectorField& nHatPatch = nHatb[patchi];
413
414 vectorField AfHatPatch
415 (
416 mesh_.Sf().boundaryField()[patchi]
417 /mesh_.magSf().boundaryField()[patchi]
418 );
419
420 reactingMultiphaseEuler::alphaContactAngleFvPatchScalarField::
421 thetaPropsTable::const_iterator tp =
422 acap.thetaProps()
423 .find(phasePairKey(phase1.name(), phase2.name()));
424
425 if (tp == acap.thetaProps().end())
426 {
428 << "Cannot find interface "
430 << "\n in table of theta properties for patch "
431 << acap.patch().name()
432 << exit(FatalError);
433 }
434
435 bool matched = (tp.key().first() == phase1.name());
436
437 scalar theta0 = degToRad(tp().theta0(matched));
438 scalarField theta(boundary[patchi].size(), theta0);
439
440 scalar uTheta = tp().uTheta();
441
442 // Calculate the dynamic contact angle if required
443 if (uTheta > SMALL)
444 {
445 const scalar thetaA = degToRad(tp().thetaA(matched));
446 const scalar thetaR = degToRad(tp().thetaR(matched));
447
448 // Calculated the component of the velocity parallel to the wall
449 vectorField Uwall
450 (
451 phase1.U()().boundaryField()[patchi].patchInternalField()
452 - phase1.U()().boundaryField()[patchi]
453 );
454 Uwall -= (AfHatPatch & Uwall)*AfHatPatch;
455
456 // Find the direction of the interface parallel to the wall
457 vectorField nWall
458 (
459 nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch
460 );
461
462 // Normalise nWall
463 nWall /= (mag(nWall) + SMALL);
464
465 // Calculate Uwall resolved normal to the interface parallel to
466 // the interface
467 scalarField uwall(nWall & Uwall);
468
469 theta += (thetaA - thetaR)*tanh(uwall/uTheta);
470 }
471
472
473 // Reset nHatPatch to correspond to the contact angle
474
475 scalarField a12(nHatPatch & AfHatPatch);
476
477 scalarField b1(cos(theta));
478
479 scalarField b2(nHatPatch.size());
480
481 forAll(b2, facei)
482 {
483 b2[facei] = cos(acos(a12[facei]) - theta[facei]);
484 }
485
486 scalarField det(1 - a12*a12);
487
488 scalarField a((b1 - a12*b2)/det);
489 scalarField b((b2 - a12*b1)/det);
490
491 nHatPatch = a*AfHatPatch + b*nHatPatch;
492
493 nHatPatch /= (mag(nHatPatch) + deltaN_.value());
494 }
495 }
496}
497
498
499Foam::tmp<Foam::volScalarField> Foam::multiphaseSystem::K
500(
501 const phaseModel& phase1,
502 const phaseModel& phase2
503) const
504{
505 tmp<surfaceVectorField> tnHatfv = nHatfv(phase1, phase2);
506
507 correctContactAngle(phase1, phase2, tnHatfv.ref().boundaryFieldRef());
508
509 // Simple expression for curvature
510 return -fvc::div(tnHatfv & mesh_.Sf());
511}
512
513
514// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
515
517(
518 const fvMesh& mesh
519)
520:
522
523 alphas_
524 (
526 (
527 "alphas",
528 mesh_.time().timeName(),
529 mesh,
530 IOobject::NO_READ,
531 IOobject::AUTO_WRITE
532 ),
533 mesh,
535 ),
536
537 cAlphas_(lookup("interfaceCompression")),
538
539 deltaN_
540 (
541 "deltaN",
542 1e-8/pow(average(mesh_.V()), 1.0/3.0)
543 )
544{
546 {
547 volScalarField& alphai = phases()[phasei];
548 mesh_.setFluxRequired(alphai.name());
549 }
550}
551
552
553// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
554
555Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseSystem::surfaceTension
556(
557 const phaseModel& phase1
558) const
559{
560 auto tSurfaceTension = surfaceScalarField::New
561 (
562 "surfaceTension",
564 mesh_,
565 dimensionedScalar(dimensionSet(1, -2, -2, 0, 0), Zero)
566 );
567
568 tSurfaceTension.ref().setOriented();
569
570 forAll(phases(), phasej)
571 {
572 const phaseModel& phase2 = phases()[phasej];
573
574 if (&phase2 != &phase1)
575 {
576 phasePairKey key12(phase1.name(), phase2.name());
577
578 cAlphaTable::const_iterator cAlpha(cAlphas_.find(key12));
579
580 if (cAlpha != cAlphas_.end())
581 {
582 tSurfaceTension.ref() +=
584 *(
587 );
588 }
589 }
590 }
591
592 tSurfaceTension->setOriented();
593
594 return tSurfaceTension;
595}
596
597
598Foam::tmp<Foam::volScalarField>
600{
601 auto tnearInt = volScalarField::New
602 (
603 "nearInterface",
605 mesh_,
607 );
608
610 {
611 tnearInt.ref() = max
612 (
613 tnearInt(),
614 pos0(phases()[phasei] - 0.01)*pos0(0.99 - phases()[phasei])
615 );
616 }
617
618 return tnearInt;
619}
620
621
623{
624 const Time& runTime = mesh_.time();
625
626 const dictionary& alphaControls = mesh_.solverDict("alpha");
627 label nAlphaSubCycles(alphaControls.get<label>("nAlphaSubCycles"));
628
629 bool LTS = fv::localEulerDdt::enabled(mesh_);
630
631 if (nAlphaSubCycles > 1)
632 {
633 tmp<volScalarField> trSubDeltaT;
634
635 if (LTS)
636 {
637 trSubDeltaT =
639 }
640
641 List<volScalarField*> alphaPtrs(phases().size());
642 PtrList<surfaceScalarField> alphaPhiSums(phases().size());
643
645 {
648
649 alphaPtrs[phasei] = &alpha;
650
651 alphaPhiSums.set
652 (
653 phasei,
655 (
657 (
658 "phiSum" + alpha.name(),
660 mesh_
661 ),
662 mesh_,
663 dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero)
664 )
665 );
666 }
667
668 for
669 (
670 subCycleTime alphaSubCycle
671 (
672 const_cast<Time&>(runTime),
674 );
675 !(++alphaSubCycle).end();
676 )
677 {
678 solveAlphas();
679
681 {
682 alphaPhiSums[phasei] += phases()[phasei].alphaPhi();
683 }
684 }
685
687 {
689 if (phase.stationary()) continue;
690
691 phase.alphaPhiRef() = alphaPhiSums[phasei]/nAlphaSubCycles;
692 }
693 }
694 else
695 {
696 solveAlphas();
697 }
698
700 {
702 if (phase.stationary()) continue;
703
704 phase.alphaRhoPhiRef() =
705 fvc::interpolate(phase.rho())*phase.alphaPhi();
706
707 phase.clamp_range(SMALL, 1 - SMALL);
708 }
709
710 calcAlphas();
711}
712
713
714// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
MULES: Multidimensional universal limiter for explicit solution.
const volScalarField & alpha1
phaseModel & phase1
const volScalarField & alpha2
phaseModel & phase2
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weights, const label comm=UPstream::worldComm) const
Return the global weighted average.
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvsPatchField< scalar >::calculatedType())
DimensionedField< scalar, volMesh > Internal
void clamp_range(const dimensioned< MinMax< Type > > &range)
Clamp field values (in-place) to the specified range.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
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 groupName(StringType base, const word &group)
Create dot-delimited name.group string.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
bool set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition List.H:469
autoPtr< T > set(const label i, const word &key, T *ptr)
Set element to pointer provided and return old element.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static tmp< volScalarField > localRSubDeltaT(const fvMesh &mesh, const label nAlphaSubCycles)
Calculate and return the reciprocal of the local sub-cycling.
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
Incompressible multi-phase mixture with built in solution for the phase fractions with interface comp...
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
multiphaseSystem(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
const PtrDictionary< phaseModel > & phases() const
Return the phases.
tmp< surfaceScalarField > surfaceTension(const phaseModel &phase) const
void solve()
Solve for the mixture phase-fractions.
const Time & time() const noexcept
Return time registry.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
virtual tmp< volScalarField > divU() const =0
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi)).
const volVectorField & U() const
Definition phaseModel.H:198
const surfaceScalarField & phi() const
Definition phaseModel.H:218
const word & name() const
Definition phaseModel.H:166
An ordered or unorder pair of phase names. Typically specified as follows.
Class to represent a system of phases and model interfacial transfers between them.
Definition phaseSystem.H:72
const fvMesh & mesh_
Reference to the mesh.
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
const fvMesh & mesh() const
Return the mesh.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
const word & name() const
Definition phase.H:114
const dimensionedScalar & rho() const
Return const-access to phase1 density.
Definition phase.H:151
Lookup type of boundary radiation properties.
Definition lookup.H:60
void setFluxRequired(const word &name) const
Set flux-required for given name (mutable).
A class for managing sub-cycling times.
bool end() const
Return true if the number of sub-cycles has been reached.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
Return the face-interpolate of the given cell field.
A class for managing temporary objects.
Definition tmp.H:75
Upwind differencing scheme class.
Definition upwind.H:55
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition zeroField.H:50
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
faceListList boundary
dynamicFvMesh & mesh
engineTime & runTime
bool LTS
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Calculate the snGrad of the given volField.
Calculate the field for explicit evaluation of implicit and explicit sources.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the laplacian of the field.
Calculate the finiteVolume matrix for implicit and explicit sources.
zeroField Su
Definition alphaSuSp.H:1
zeroField Sp
Definition alphaSuSp.H:2
word timeName
Definition getTimeIndex.H:3
label phasei
Definition pEqn.H:27
PtrList< surfaceScalarField > alphafs(phases.size())
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin, const bool returnCorr)
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition MULES.C:27
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.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvcSnGrad.C:40
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition fvcFlux.C:27
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition fvcMeshPhi.C:183
Namespace for OpenFOAM.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
dimensionedScalar det(const dimensionedSphericalTensor &dt)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
void Sp(fvMatrix< typename Expr::value_type > &m, const Expr2 &mult, const Expr &expression)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
volScalarField & alpha
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
volScalarField & b
volScalarField & e
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
multiphaseSystem::phaseModelList & phases
const dictionary & alphaControls
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
Unit conversion functions.