Loading...
Searching...
No Matches
kinematicSingleLayer.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2025 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
30#include "fvm.H"
31#include "fvcDiv.H"
32#include "fvcLaplacian.H"
33#include "fvcSnGrad.H"
34#include "fvcReconstruct.H"
35#include "fvcVolumeIntegrate.H"
36#include "fvcFlux.H"
38#include "mappedWallPolyPatch.H"
39#include "mapDistribute.H"
40#include "filmThermoModel.H"
41
42// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
44namespace Foam
45{
46namespace regionModels
47{
49{
50
51// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
53defineTypeNameAndDebug(kinematicSingleLayer, 0);
54
55addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh);
56
57// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
58
60{
62 {
63 const dictionary& solution = this->solution().subDict("PISO");
64 solution.readEntry("momentumPredictor", momentumPredictor_);
65 solution.readIfPresent("nOuterCorr", nOuterCorr_);
66 solution.readEntry("nCorr", nCorr_);
67 solution.readEntry("nNonOrthCorr", nNonOrthCorr_);
68
69 return true;
70 }
71
72 return false;
73}
74
75
78 rho_ == filmThermo_->rho();
79 mu_ == filmThermo_->mu();
80 sigma_ == filmThermo_->sigma();
81}
82
83
85{
98 // Update fields from primary region via direct mapped
99 // (coupled) boundary conditions
100 UPrimary_.correctBoundaryConditions();
101 pPrimary_.correctBoundaryConditions();
102 rhoPrimary_.correctBoundaryConditions();
103 muPrimary_.correctBoundaryConditions();
104}
105
106
108{
110
111 auto& rhoSpPrimaryBf = rhoSpPrimary_.boundaryFieldRef();
112 auto& USpPrimaryBf = USpPrimary_.boundaryFieldRef();
113 auto& pSpPrimaryBf = pSpPrimary_.boundaryFieldRef();
114
115 // Convert accumulated source terms into per unit area per unit time
116 const scalar deltaT = time_.deltaTValue();
118 {
119 scalarField rpriMagSfdeltaT
120 (
121 (1.0/deltaT)
122 /primaryMesh().magSf().boundaryField()[patchi]
123 );
124
125 rhoSpPrimaryBf[patchi] *= rpriMagSfdeltaT;
126 USpPrimaryBf[patchi] *= rpriMagSfdeltaT;
127 pSpPrimaryBf[patchi] *= rpriMagSfdeltaT;
128 }
129
130 // Retrieve the source fields from the primary region via direct mapped
131 // (coupled) boundary conditions
132 // - fields require transfer of values for both patch AND to push the
133 // values into the first layer of internal cells
137
138 // update addedMassTotal counter
139 if (time().writeTime())
140 {
141 if (debug)
142 {
143 rhoSp_.write();
144 USp_.write();
145 pSp_.write();
146 }
147
148 scalar addedMassTotal = 0.0;
149 outputProperties().readIfPresent("addedMassTotal", addedMassTotal);
151 outputProperties().add("addedMassTotal", addedMassTotal, true);
152 addedMassTotal_ = 0.0;
153 }
154}
155
156
158{
160 (
163 (
164 pPrimary_ // pressure (mapped from primary region)
165 - pSp_ // accumulated particle impingement
166 - fvc::laplacian(sigma_, delta_) // surface tension
167 )
168 );
169}
170
171
173{
175 (
179 -rho_*gNormClipped() // hydrostatic effect only
180 )
181 );
182}
183
186{
188}
189
190
192{
194
195 // Update injection model - mass returned is mass available for injection
197
198 // Update transfer model - mass returned is mass available for transfer
200
201 // Update mass source field
203
204 turbulence_->correct();
205}
206
207
209{
210 const volScalarField deltaRho0(deltaRho_);
211
213
214 if (debug)
215 {
219 + dimensionedScalar("SMALL", dimMass*dimVolume, ROOTVSMALL);
220
221 const scalar sumLocalContErr =
222 (
224 ).value();
225
226 const scalar globalContErr =
227 (
229 ).value();
230
232
234 << "Surface film: " << type() << nl
235 << " time step continuity errors: sum local = "
236 << sumLocalContErr << ", global = " << globalContErr
237 << ", cumulative = " << cumulativeContErr_ << endl;
238 }
239}
240
241
243{
245
246 solve
247 (
250 ==
251 - rhoSp_
252 );
253}
254
255
257{
258 // Push boundary film velocity values into internal field
259 for (label i=0; i<intCoupledPatchIDs_.size(); i++)
260 {
261 label patchi = intCoupledPatchIDs_[i];
262 const polyPatch& pp = regionMesh().boundaryMesh()[patchi];
263 UIndirectList<vector>(Uw_, pp.faceCells()) =
264 U_.boundaryField()[patchi];
265 }
266 Uw_ -= nHat()*(Uw_ & nHat());
267 Uw_.correctBoundaryConditions();
268
269 Us_ = turbulence_->Us();
270}
271
272
274(
275 const volScalarField& pu,
276 const volScalarField& pp
277)
278{
280
281 // Momentum
283 (
285 + fvm::div(phi_, U_)
286 ==
287 - USp_
288 // - fvm::SuSp(rhoSp_, U_)
289 - rhoSp_*U_
291 + turbulence_->Su(U_)
292 );
293
294 fvVectorMatrix& UEqn = tUEqn.ref();
295
296 UEqn.relax();
297
299 {
300 solve
301 (
302 UEqn
303 ==
305 (
307 * (
308 regionMesh().magSf()
309 * (
310 fvc::snGrad(pu, "snGrad(p)")
311 + fvc::snGrad(pp, "snGrad(p)")*fvc::interpolate(delta_)
313 )
314 - fvc::flux(rho_*gTan())
315 )
316 )
317 );
318
319 // Remove any patch-normal components of velocity
320 U_ -= nHat()*(nHat() & U_);
321 U_.correctBoundaryConditions();
322 }
323
324 return tUEqn;
325}
326
327
329(
330 const volScalarField& pu,
331 const volScalarField& pp,
333)
334{
336
337 volScalarField rUA(1.0/UEqn.A());
338 U_ = rUA*UEqn.H();
339
342
343 surfaceScalarField phiAdd
344 (
345 "phiAdd",
346 regionMesh().magSf()
347 * (
348 fvc::snGrad(pu, "snGrad(p)")
349 + fvc::snGrad(pp, "snGrad(p)")*fvc::interpolate(delta_)
350 )
351 - fvc::flux(rho_*gTan())
352 );
353 constrainFilmField(phiAdd, 0.0);
354
356 (
357 "phid",
358 fvc::flux(U_*rho_) - deltarUAf*phiAdd*rhof
359 );
361
362 surfaceScalarField ddrhorUAppf
363 (
364 "deltaCoeff",
366 );
367
369
370 for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
371 {
372 // Film thickness equation
373 fvScalarMatrix deltaEqn
374 (
377 - fvm::laplacian(ddrhorUAppf, delta_)
378 ==
379 - rhoSp_
380 );
381
382 deltaEqn.solve();
383
384 if (nonOrth == nNonOrthCorr_)
385 {
386 phiAdd +=
389 * regionMesh().magSf();
390
391 phi_ == deltaEqn.flux();
392 }
393 }
394
395 // Bound film thickness by a minimum of zero
396 delta_.clamp_min(0);
397
398 // Update U field
399 U_ -= fvc::reconstruct(deltarUAf*phiAdd);
400
401 // Remove any patch-normal components of velocity
402 U_ -= nHat()*(nHat() & U_);
403
405
406 // Update film wall and surface velocities
408
409 // Continuity check
411}
412
413
414// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
415
416kinematicSingleLayer::kinematicSingleLayer
417(
418 const word& modelType,
419 const fvMesh& mesh,
420 const dimensionedVector& g,
421 const word& regionType,
422 const bool readFields
423)
424:
425 surfaceFilmRegionModel(modelType, mesh, g, regionType),
426
427 momentumPredictor_(solution().subDict("PISO").lookup("momentumPredictor")),
428 nOuterCorr_(solution().subDict("PISO").getOrDefault("nOuterCorr", 1)),
429 nCorr_(solution().subDict("PISO").get<label>("nCorr")),
430 nNonOrthCorr_
431 (
432 solution().subDict("PISO").get<label>("nNonOrthCorr")
433 ),
434
435 cumulativeContErr_(0.0),
436
437 deltaSmall_("deltaSmall", dimLength, SMALL),
438 deltaCoLimit_(solution().getOrDefault<scalar>("deltaCoLimit", 1e-4)),
439
440 rho_
441 (
443 (
444 "rhof",
445 time().timeName(),
446 regionMesh(),
447 IOobject::NO_READ,
448 IOobject::AUTO_WRITE
449 ),
450 regionMesh(),
452 fvPatchFieldBase::zeroGradientType()
453 ),
454 mu_
455 (
457 (
458 "muf",
459 time().timeName(),
460 regionMesh(),
461 IOobject::NO_READ,
462 IOobject::AUTO_WRITE
463 ),
464 regionMesh(),
466 fvPatchFieldBase::zeroGradientType()
467 ),
468 sigma_
469 (
471 (
472 "sigmaf",
473 time().timeName(),
474 regionMesh(),
475 IOobject::NO_READ,
476 IOobject::AUTO_WRITE
477 ),
478 regionMesh(),
480 fvPatchFieldBase::zeroGradientType()
481 ),
482
483 delta_
484 (
486 (
487 "deltaf",
488 time().timeName(),
489 regionMesh(),
490 IOobject::MUST_READ,
491 IOobject::AUTO_WRITE
492 ),
493 regionMesh()
494 ),
495 alpha_
496 (
498 (
499 "alpha",
500 time().timeName(),
501 regionMesh(),
502 IOobject::READ_IF_PRESENT,
503 IOobject::AUTO_WRITE
504 ),
505 regionMesh(),
507 fvPatchFieldBase::zeroGradientType()
508 ),
509 U_
510 (
512 (
513 "Uf",
514 time().timeName(),
515 regionMesh(),
516 IOobject::MUST_READ,
517 IOobject::AUTO_WRITE
518 ),
519 regionMesh()
520 ),
521 Us_
522 (
524 (
525 "Usf",
526 time().timeName(),
527 regionMesh(),
528 IOobject::NO_READ,
529 IOobject::NO_WRITE
530 ),
531 U_,
532 fvPatchFieldBase::zeroGradientType()
533 ),
534 Uw_
535 (
537 (
538 "Uwf",
539 time().timeName(),
540 regionMesh(),
541 IOobject::NO_READ,
542 IOobject::NO_WRITE
543 ),
544 U_,
545 fvPatchFieldBase::zeroGradientType()
546 ),
547 deltaRho_
548 (
550 (
551 delta_.name() + "*" + rho_.name(),
552 time().timeName(),
553 regionMesh(),
554 IOobject::NO_READ,
555 IOobject::NO_WRITE
556 ),
557 regionMesh(),
558 dimensionedScalar(delta_.dimensions()*rho_.dimensions(), Zero),
559 fvPatchFieldBase::zeroGradientType()
560 ),
561
562 phi_
563 (
565 (
566 "phi",
567 time().timeName(),
568 regionMesh(),
569 IOobject::NO_READ,
570 IOobject::AUTO_WRITE
571 ),
572 regionMesh(),
574 ),
575
576 primaryMassTrans_
577 (
579 (
580 "primaryMassTrans",
581 time().timeName(),
582 regionMesh(),
583 IOobject::NO_READ,
584 IOobject::NO_WRITE
585 ),
586 regionMesh(),
588 fvPatchFieldBase::zeroGradientType()
589 ),
590 cloudMassTrans_
591 (
593 (
594 "cloudMassTrans",
595 time().timeName(),
596 regionMesh(),
597 IOobject::NO_READ,
598 IOobject::NO_WRITE
599 ),
600 regionMesh(),
602 fvPatchFieldBase::zeroGradientType()
603 ),
604 cloudDiameterTrans_
605 (
607 (
608 "cloudDiameterTrans",
609 time().timeName(),
610 regionMesh(),
611 IOobject::NO_READ,
612 IOobject::NO_WRITE
613 ),
614 regionMesh(),
615 dimensionedScalar("minus1", dimLength, -1.0),
616 fvPatchFieldBase::zeroGradientType()
617 ),
618
619 USp_
620 (
622 (
623 "USpf",
624 time().timeName(),
625 regionMesh(),
626 IOobject::NO_READ,
627 IOobject::NO_WRITE
628 ),
629 regionMesh(),
631 this->mappedPushedFieldPatchTypes<vector>()
632 ),
633 pSp_
634 (
636 (
637 "pSpf",
638 time_.timeName(),
639 regionMesh(),
640 IOobject::NO_READ,
641 IOobject::NO_WRITE
642 ),
643 regionMesh(),
645 this->mappedPushedFieldPatchTypes<scalar>()
646 ),
647 rhoSp_
648 (
650 (
651 "rhoSpf",
652 time_.timeName(),
653 regionMesh(),
654 IOobject::NO_READ,
655 IOobject::NO_WRITE
656 ),
657 regionMesh(),
659 this->mappedPushedFieldPatchTypes<scalar>()
660 ),
661
662 USpPrimary_
663 (
665 (
666 USp_.name(), // must have same name as USp_ to enable mapping
667 time().timeName(),
668 primaryMesh(),
669 IOobject::NO_READ,
670 IOobject::NO_WRITE
671 ),
672 primaryMesh(),
673 dimensionedVector(USp_.dimensions(), Zero)
674 ),
675 pSpPrimary_
676 (
678 (
679 pSp_.name(), // must have same name as pSp_ to enable mapping
680 time().timeName(),
681 primaryMesh(),
682 IOobject::NO_READ,
683 IOobject::NO_WRITE
684 ),
685 primaryMesh(),
686 dimensionedScalar(pSp_.dimensions(), Zero)
687 ),
688 rhoSpPrimary_
689 (
691 (
692 rhoSp_.name(), // must have same name as rhoSp_ to enable mapping
693 time().timeName(),
694 primaryMesh(),
695 IOobject::NO_READ,
696 IOobject::NO_WRITE
697 ),
698 primaryMesh(),
699 dimensionedScalar(rhoSp_.dimensions(), Zero)
700 ),
701
702 UPrimary_
703 (
705 (
706 "U", // must have same name as U to enable mapping
707 time().timeName(),
708 regionMesh(),
709 IOobject::NO_READ,
710 IOobject::NO_WRITE
711 ),
712 regionMesh(),
714 this->mappedFieldAndInternalPatchTypes<vector>()
715 ),
716 pPrimary_
717 (
719 (
720 "p", // must have same name as p to enable mapping
721 time().timeName(),
722 regionMesh(),
723 IOobject::NO_READ,
724 IOobject::NO_WRITE
725 ),
726 regionMesh(),
728 this->mappedFieldAndInternalPatchTypes<scalar>()
729 ),
730 rhoPrimary_
731 (
733 (
734 "rho", // must have same name as rho to enable mapping
735 time().timeName(),
736 regionMesh(),
737 IOobject::NO_READ,
738 IOobject::NO_WRITE
739 ),
740 regionMesh(),
742 this->mappedFieldAndInternalPatchTypes<scalar>()
743 ),
744 muPrimary_
745 (
747 (
748 "thermo:mu", // must have same name as mu to enable mapping
749 time().timeName(),
750 regionMesh(),
751 IOobject::NO_READ,
752 IOobject::NO_WRITE
753 ),
754 regionMesh(),
756 this->mappedFieldAndInternalPatchTypes<scalar>()
757 ),
758
759 filmThermo_(filmThermoModel::New(*this, coeffs_)),
760
761 availableMass_(regionMesh().nCells(), Zero),
762
763 injection_(*this, coeffs_),
764
765 transfer_(*this, coeffs_),
766
767 turbulence_(filmTurbulenceModel::New(*this, coeffs_)),
768
769 forces_(*this, coeffs_),
770
771 addedMassTotal_(0.0)
772{
773 if (readFields)
774 {
776
777 correctAlpha();
778
780
782
784 (
786 (
787 "phi",
788 time().timeName(),
789 regionMesh(),
793 ),
795 );
796
797 phi_ == phi0;
798 }
799}
800
801
802// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
805{}
806
807
808// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
809
811(
812 const label patchi,
813 const label facei,
814 const scalar massSource,
815 const vector& momentumSource,
816 const scalar pressureSource,
817 const scalar energySource
818)
819{
821 << "\nSurface film: " << type() << ": adding to film source:" << nl
822 << " mass = " << massSource << nl
823 << " momentum = " << momentumSource << nl
824 << " pressure = " << pressureSource << endl;
825
826 rhoSpPrimary_.boundaryFieldRef()[patchi][facei] -= massSource;
827 USpPrimary_.boundaryFieldRef()[patchi][facei] -= momentumSource;
828 pSpPrimary_.boundaryFieldRef()[patchi][facei] -= pressureSource;
829
830 addedMassTotal_ += massSource;
831}
832
833
835{
837
839
841
843
845
847
848 correctAlpha();
849
850 // Reset transfer fields
855}
856
857
859{
861
862 // Update sub-models to provide updated source contributions
864
865 // Solve continuity for deltaRho_
867
868 // Implicit pressure source coefficient - constant
869 tmp<volScalarField> tpp(this->pp());
870
871 for (int oCorr=1; oCorr<=nOuterCorr_; oCorr++)
872 {
873 // Explicit pressure source contribution - varies with delta_
874 tmp<volScalarField> tpu(this->pu());
875
876 // Solve for momentum for U_
878 fvVectorMatrix& UEqn = tUEqn.ref();
879
880 // Film thickness correction loop
881 for (int corr=1; corr<=nCorr_; corr++)
882 {
883 // Solve thickness for delta_
884 solveThickness(tpu(), tpp(), UEqn);
885 }
887
888 // Update deltaRho_ with new delta_
890}
891
892
894{
896
897 // Reset source terms for next time integration
899}
900
901
903{
904 scalar CoNum = 0.0;
905
906 if (regionMesh().nInternalFaces() > 0)
907 {
908 const scalarField sumPhi
909 (
910 fvc::surfaceSum(mag(phi_))().primitiveField()
911 / (deltaRho_.primitiveField() + ROOTVSMALL)
912 );
913
914 forAll(delta_, i)
915 {
916 if ((delta_[i] > deltaCoLimit_) && (alpha_[i] > 0.5))
917 {
918 CoNum = max(CoNum, sumPhi[i]/(delta_[i]*magSf()[i]));
919 }
920 }
921
922 CoNum *= 0.5*time_.deltaTValue();
923 }
924
925 reduce(CoNum, maxOp<scalar>());
927 Info<< "Film max Courant number: " << CoNum << endl;
928
929 return CoNum;
930}
931
934{
935 return U_;
936}
937
940{
941 return Us_;
942}
943
946{
947 return Uw_;
948}
949
952{
953 return deltaRho_;
954}
955
958{
959 return phi_;
960}
961
964{
965 return rho_;
966}
967
968
970{
972 << "T field not available for " << type() << abort(FatalError);
973
974 return volScalarField::null();
975}
976
977
979{
981 << "Ts field not available for " << type() << abort(FatalError);
982
983 return volScalarField::null();
984}
985
986
988{
990 << "Tw field not available for " << type() << abort(FatalError);
991
992 return volScalarField::null();
993}
994
995
997{
999 << "hs field not available for " << type() << abort(FatalError);
1000
1001 return volScalarField::null();
1002}
1003
1004
1006{
1008 << "Cp field not available for " << type() << abort(FatalError);
1009
1010 return volScalarField::null();
1011}
1012
1013
1015{
1017 << "kappa field not available for " << type() << abort(FatalError);
1018
1019 return volScalarField::null();
1020}
1021
1026}
1027
1032}
1033
1036{
1037 return cloudDiameterTrans_;
1038}
1039
1040
1042{
1043 Info<< "\nSurface film: " << type() << endl;
1044
1045 const scalarField& deltaInternal = delta_;
1046 const vectorField& Uinternal = U_;
1047 scalar addedMassTotal = 0.0;
1048 outputProperties().readIfPresent("addedMassTotal", addedMassTotal);
1049 addedMassTotal += returnReduce(addedMassTotal_, sumOp<scalar>());
1050
1051 MinMax<scalar> limits_U = gMinMaxMag(Uinternal);
1052 MinMax<scalar> limits_delta = gMinMax(deltaInternal);
1053
1054 Info<< indent << "added mass = " << addedMassTotal << nl
1055 << indent << "current mass = "
1056 << gWeightedSum(magSf(), deltaRho_) << nl
1057 << indent << "min/max(mag(U)) = "
1058 << limits_U.min() << ", "
1059 << limits_U.max() << nl
1060 << indent << "min/max(delta) = "
1061 << limits_delta.min() << ", "
1062 << limits_delta.max() << nl
1063 << indent << "coverage = "
1064 << gWeightedAverage(magSf(), alpha_.primitiveField()) << nl;
1065
1066 injection_.info(Info);
1067 transfer_.info(Info);
1068}
1069
1070
1072{
1074 (
1077 primaryMesh(),
1079 );
1080}
1081
1082
1084(
1085 const label i
1086) const
1087{
1089 (
1090 IOobject::scopedName(typeName, "Srho(" + Foam::name(i) + ")"),
1101 (
1104 primaryMesh(),
1106 );
1107}
1108
1109
1110// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1111
1112} // End namespace surfaceFilmModels
1113} // End namespace regionModels
1114} // End namespace Foam
1115
1116// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
dimensionedScalar totalMass
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const uniformDimensionedVectorField & g
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
void clamp_min(const Type &lower)
Impose lower (floor) clamp on the field values (in-place).
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
void correctBoundaryConditions()
Correct boundary field.
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.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ 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 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
A min/max value pair with additional methods. In addition to conveniently storing values,...
Definition MinMax.H:106
const T & max() const noexcept
The max value.
Definition MinMax.H:217
const T & min() const noexcept
The min value.
Definition MinMax.H:207
scalar deltaTValue() const noexcept
Return time step value.
Definition TimeStateI.H:49
dimensionedScalar deltaT() const
Return time step.
Definition TimeStateI.H:61
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
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.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition dictionary.C:625
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition fvMatrix.C:1415
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Template invariant parts for fvPatchField.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Lookup type of boundary radiation properties.
Definition lookup.H:60
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
const Time & time_
Reference to the time database.
Definition regionModel.H:97
const dictionary & solution() const
Return the solution dictionary.
const Time & time() const noexcept
Return the reference to the time database.
dictionary coeffs_
Model coefficients dictionary.
const fvMesh & regionMesh() const
Return the region mesh database.
virtual void preEvolveRegion()
Pre-evolve region.
const fvMesh & primaryMesh() const noexcept
Return the reference to the primary mesh database.
const IOdictionary & outputProperties() const
Return const access to the output properties dictionary.
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
wordList mappedFieldAndInternalPatchTypes() const
Return boundary types for mapped field patches.
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
virtual const volVectorField & nHat() const
Return the patch normal vectors.
wordList mappedPushedFieldPatchTypes() const
Return boundary types for pushed mapped field patches.
static autoPtr< surfaceFilmModel > New(const fvMesh &mesh, const dimensionedVector &g, const word &regionType="surfaceFilm")
Return a reference to the selected surface film model.
tmp< fvVectorMatrix > correct(volVectorField &U)
Return (net) force system.
Definition forceList.C:73
virtual void correct(scalarField &availableMass, volScalarField &massToInject, volScalarField &diameterToInject)
Correct.
Kinematic form of single-cell layer surface film model.
tmp< volScalarField > gNormClipped() const
Return the gravity normal-to-patch component contribution.
virtual void solveThickness(const volScalarField &pu, const volScalarField &pp, fvVectorMatrix &UEqn)
Solve coupled velocity-thickness equations.
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
virtual const volVectorField & U() const
Return the film velocity [m/s].
scalar addedMassTotal_
Cumulative mass added via sources [kg].
virtual const volScalarField & hs() const
Return the film surface enthalpy [J/kg].
tmp< volVectorField > gTan() const
Return the gravity tangential component contributions.
virtual const volScalarField & rho() const
Return the film density [kg/m3].
virtual void updateSurfaceVelocities()
Update film surface velocities.
tmp< volScalarField > mass() const
Return the current film mass.
virtual const volScalarField & cloudDiameterTrans() const
Return the parcel diameters originating from film to cloud.
virtual const volScalarField & kappa() const
Return the film thermal conductivity [W/m/K].
virtual void addSources(const label patchi, const label facei, const scalar massSource, const vector &momentumSource, const scalar pressureSource, const scalar energySource=0)
External hook to add sources to the film.
virtual const volVectorField & Us() const
Return the film surface velocity [m/s].
virtual const surfaceScalarField & phi() const
Return the film flux [kg.m/s].
virtual const volScalarField & T() const
Return the film mean temperature [K].
virtual const volScalarField & deltaRho() const
Return the film thickness*density (helper field) [kg/m3].
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
virtual const volScalarField & Tw() const
Return the film wall temperature [K].
volScalarField deltaRho_
Film thickness*density (helper field) [kg/m2].
volScalarField cloudDiameterTrans_
Parcel diameters originating from film to cloud.
transferModelList transfer_
Transfer with the continuous phase.
virtual const volScalarField & Cp() const
Return the film specific heat capacity [J/kg/K].
virtual void resetPrimaryRegionSourceTerms()
Reset source term fields.
virtual void correctAlpha()
Correct film coverage field.
surfaceScalarField phi_
Mass flux (includes film thickness) [kg.m/s].
virtual tmp< fvVectorMatrix > solveMomentum(const volScalarField &pu, const volScalarField &pp)
Solve for film velocity.
scalar deltaCoLimit_
Film thickness above which Courant number calculation in valid.
virtual void correctThermoFields()
Correct the thermo fields.
virtual const volScalarField & cloudMassTrans() const
Return the film mass available for transfer to cloud.
virtual tmp< volScalarField > pp()
Implicit pressure source coefficient.
virtual void updateSubmodels()
Update the film sub-models.
virtual void transferPrimaryRegionSourceFields()
Transfer source fields from the primary region to the film region.
volScalarField primaryMassTrans_
Film mass available for transfer to the primary region.
void constrainFilmField(Type &field, const typename Type::cmptType &value)
Constrain a film region master/slave boundaries of a field to a.
virtual scalar CourantNumber() const
Courant number evaluation.
volScalarField cloudMassTrans_
Film mass available for transfer to cloud.
autoPtr< filmThermoModel > filmThermo_
Film thermo model.
virtual const volVectorField & Uw() const
Return the film wall velocity [m/s].
scalarField availableMass_
Available mass for transfer via sub-models.
virtual tmp< volScalarField::Internal > Srho() const
Return total mass source - Eulerian phase only.
virtual tmp< volScalarField > primaryMassTrans() const
Return mass transfer source - Eulerian phase only.
virtual bool read()
Read control parameters from dictionary.
virtual tmp< volScalarField > pu()
Explicit pressure source contribution.
autoPtr< filmTurbulenceModel > turbulence_
Turbulence model.
volScalarField alpha_
Film coverage indicator, 1 = covered, 0 = uncovered [].
virtual tmp< volScalarField::Internal > Sh() const
Return enthalpy source - Eulerian phase only.
const dimensionedVector & g() const
Return the acceleration due to gravity.
virtual bool read()
Read control parameters from dictionary.
virtual void correct(scalarField &availableMass, volScalarField &massToTransfer)
Correct kinematic transfers.
void setFluxRequired(const word &name) const
Set flux-required for given name (mutable).
Selector class for relaxation factors, solver type and solution.
Definition solution.H:95
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
fvVectorMatrix & UEqn
Definition UEqn.H:13
surfaceScalarField phid("phid", fvc::interpolate(psi) *(fvc::flux(HbyA)+MRF.zeroFilter(rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho))))
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
Calculate the laplacian of the given field.
Reconstruct volField from a face flux field.
Calculate the snGrad of the given volField.
Volume integrate volField creating a volField.
word timeName
Definition getTimeIndex.H:3
#define DebugInFunction
Report an information message using Foam::Info.
#define InfoInFunction
Report an information message using Foam::Info.
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
Namespace for handling debugging switches.
Definition debug.C:45
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< 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
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
Namespace for OpenFOAM.
const dimensionSet dimPressure
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
Type gWeightedAverage(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted average of a field, using the mag() of the weights.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimEnergy
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
fvMatrix< scalar > fvScalarMatrix
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Type gWeightedSum(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted sum (integral) of a field, using the mag() of the weights.
scalarMinMax gMinMaxMag(const FieldField< Field, Type > &f)
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fvMatrix< vector > fvVectorMatrix
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
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...
const dimensionSet dimVolume(pow3(dimLength))
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
const dimensionSet dimDensity
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
CEqn solve()
volScalarField & e
scalar sumLocalContErr
scalar globalContErr
scalar CoNum
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299