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-2020 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{
48namespace surfaceFilmModels
49{
50
51// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
52
54
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
77{
78 rho_ == filmThermo_->rho();
79 mu_ == filmThermo_->mu();
80 sigma_ == filmThermo_->sigma();
81}
82
83
85{
87
91}
92
93
95{
97
98 // Update fields from primary region via direct mapped
99 // (coupled) boundary conditions
104}
105
106
108{
110
111 volScalarField::Boundary& rhoSpPrimaryBf =
113
114 volVectorField::Boundary& USpPrimaryBf =
116
117 volScalarField::Boundary& pSpPrimaryBf =
119
120 // Convert accumulated source terms into per unit area per unit time
121 const scalar deltaT = time_.deltaTValue();
123 {
124 scalarField rpriMagSfdeltaT
125 (
126 (1.0/deltaT)
127 /primaryMesh().magSf().boundaryField()[patchi]
128 );
129
130 rhoSpPrimaryBf[patchi] *= rpriMagSfdeltaT;
131 USpPrimaryBf[patchi] *= rpriMagSfdeltaT;
132 pSpPrimaryBf[patchi] *= rpriMagSfdeltaT;
133 }
134
135 // Retrieve the source fields from the primary region via direct mapped
136 // (coupled) boundary conditions
137 // - fields require transfer of values for both patch AND to push the
138 // values into the first layer of internal cells
142
143 // update addedMassTotal counter
144 if (time().writeTime())
145 {
146 if (debug)
147 {
148 rhoSp_.write();
149 USp_.write();
150 pSp_.write();
151 }
152
153 scalar addedMassTotal = 0.0;
154 outputProperties().readIfPresent("addedMassTotal", addedMassTotal);
155 addedMassTotal += returnReduce(addedMassTotal_, sumOp<scalar>());
156 outputProperties().add("addedMassTotal", addedMassTotal, true);
157 addedMassTotal_ = 0.0;
158 }
159}
160
161
163{
165 (
167 (
169 (
170 typeName + ":pu",
171 time_.timeName(),
172 regionMesh(),
175 ),
176 pPrimary_ // pressure (mapped from primary region)
177 - pSp_ // accumulated particle impingement
178 - fvc::laplacian(sigma_, delta_) // surface tension
179 )
180 );
181}
182
183
185{
187 (
189 (
191 (
192 typeName + ":pp",
193 time_.timeName(),
194 regionMesh(),
197 ),
198 -rho_*gNormClipped() // hydrostatic effect only
199 )
200 );
201}
202
203
205{
207}
208
209
211{
213
214 // Update injection model - mass returned is mass available for injection
216
217 // Update transfer model - mass returned is mass available for transfer
219
220 // Update mass source field
222
223 turbulence_->correct();
224}
225
226
228{
229 const volScalarField deltaRho0(deltaRho_);
230
232
233 if (debug)
234 {
238 + dimensionedScalar("SMALL", dimMass*dimVolume, ROOTVSMALL);
239
240 const scalar sumLocalContErr =
241 (
243 ).value();
244
245 const scalar globalContErr =
246 (
248 ).value();
249
251
253 << "Surface film: " << type() << nl
254 << " time step continuity errors: sum local = "
255 << sumLocalContErr << ", global = " << globalContErr
256 << ", cumulative = " << cumulativeContErr_ << endl;
257 }
258}
259
260
262{
264
265 solve
266 (
268 + fvc::div(phi_)
269 ==
270 - rhoSp_
271 );
272}
273
274
276{
277 // Push boundary film velocity values into internal field
278 for (label i=0; i<intCoupledPatchIDs_.size(); i++)
279 {
280 label patchi = intCoupledPatchIDs_[i];
281 const polyPatch& pp = regionMesh().boundaryMesh()[patchi];
282 UIndirectList<vector>(Uw_, pp.faceCells()) =
283 U_.boundaryField()[patchi];
284 }
285 Uw_ -= nHat()*(Uw_ & nHat());
287
288 Us_ = turbulence_->Us();
289}
290
291
293(
294 const volScalarField& pu,
295 const volScalarField& pp
296)
297{
299
300 // Momentum
302 (
304 + fvm::div(phi_, U_)
305 ==
306 - USp_
307 // - fvm::SuSp(rhoSp_, U_)
308 - rhoSp_*U_
310 + turbulence_->Su(U_)
311 );
312
313 fvVectorMatrix& UEqn = tUEqn.ref();
314
315 UEqn.relax();
316
318 {
319 solve
320 (
321 UEqn
322 ==
324 (
326 * (
327 regionMesh().magSf()
328 * (
329 fvc::snGrad(pu, "snGrad(p)")
330 + fvc::snGrad(pp, "snGrad(p)")*fvc::interpolate(delta_)
332 )
333 - fvc::flux(rho_*gTan())
334 )
335 )
336 );
337
338 // Remove any patch-normal components of velocity
339 U_ -= nHat()*(nHat() & U_);
341 }
342
343 return tUEqn;
344}
345
346
348(
349 const volScalarField& pu,
350 const volScalarField& pp,
352)
353{
355
356 volScalarField rUA(1.0/UEqn.A());
357 U_ = rUA*UEqn.H();
358
361
362 surfaceScalarField phiAdd
363 (
364 "phiAdd",
365 regionMesh().magSf()
366 * (
367 fvc::snGrad(pu, "snGrad(p)")
368 + fvc::snGrad(pp, "snGrad(p)")*fvc::interpolate(delta_)
369 )
370 - fvc::flux(rho_*gTan())
371 );
372 constrainFilmField(phiAdd, 0.0);
373
375 (
376 "phid",
377 fvc::flux(U_*rho_) - deltarUAf*phiAdd*rhof
378 );
380
381 surfaceScalarField ddrhorUAppf
382 (
383 "deltaCoeff",
385 );
386
388
389 for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
390 {
391 // Film thickness equation
392 fvScalarMatrix deltaEqn
393 (
396 - fvm::laplacian(ddrhorUAppf, delta_)
397 ==
398 - rhoSp_
399 );
400
401 deltaEqn.solve();
402
403 if (nonOrth == nNonOrthCorr_)
404 {
405 phiAdd +=
408 * regionMesh().magSf();
409
410 phi_ == deltaEqn.flux();
411 }
412 }
413
414 // Bound film thickness by a minimum of zero
415 delta_.max(0.0);
416
417 // Update U field
418 U_ -= fvc::reconstruct(deltarUAf*phiAdd);
419
420 // Remove any patch-normal components of velocity
421 U_ -= nHat()*(nHat() & U_);
422
424
425 // Update film wall and surface velocities
427
428 // Continuity check
430}
431
432
433// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
434
436(
437 const word& modelType,
438 const fvMesh& mesh,
439 const dimensionedVector& g,
440 const word& regionType,
441 const bool readFields
442)
443:
444 surfaceFilmRegionModel(modelType, mesh, g, regionType),
445
446 momentumPredictor_(solution().subDict("PISO").lookup("momentumPredictor")),
447 nOuterCorr_(solution().subDict("PISO").getOrDefault("nOuterCorr", 1)),
448 nCorr_(solution().subDict("PISO").get<label>("nCorr")),
449 nNonOrthCorr_
450 (
451 solution().subDict("PISO").get<label>("nNonOrthCorr")
452 ),
453
454 cumulativeContErr_(0.0),
455
456 deltaSmall_("deltaSmall", dimLength, SMALL),
457 deltaCoLimit_(solution().getOrDefault<scalar>("deltaCoLimit", 1e-4)),
458
459 rho_
460 (
462 (
463 "rhof",
464 time().timeName(),
465 regionMesh(),
466 IOobject::NO_READ,
467 IOobject::AUTO_WRITE
468 ),
469 regionMesh(),
471 zeroGradientFvPatchScalarField::typeName
472 ),
473 mu_
474 (
476 (
477 "muf",
478 time().timeName(),
479 regionMesh(),
480 IOobject::NO_READ,
481 IOobject::AUTO_WRITE
482 ),
483 regionMesh(),
485 zeroGradientFvPatchScalarField::typeName
486 ),
487 sigma_
488 (
490 (
491 "sigmaf",
492 time().timeName(),
493 regionMesh(),
494 IOobject::NO_READ,
495 IOobject::AUTO_WRITE
496 ),
497 regionMesh(),
499 zeroGradientFvPatchScalarField::typeName
500 ),
501
502 delta_
503 (
505 (
506 "deltaf",
507 time().timeName(),
508 regionMesh(),
509 IOobject::MUST_READ,
510 IOobject::AUTO_WRITE
511 ),
512 regionMesh()
513 ),
514 alpha_
515 (
517 (
518 "alpha",
519 time().timeName(),
520 regionMesh(),
521 IOobject::READ_IF_PRESENT,
522 IOobject::AUTO_WRITE
523 ),
524 regionMesh(),
526 zeroGradientFvPatchScalarField::typeName
527 ),
528 U_
529 (
531 (
532 "Uf",
533 time().timeName(),
534 regionMesh(),
535 IOobject::MUST_READ,
536 IOobject::AUTO_WRITE
537 ),
538 regionMesh()
539 ),
540 Us_
541 (
543 (
544 "Usf",
545 time().timeName(),
546 regionMesh(),
547 IOobject::NO_READ,
548 IOobject::NO_WRITE
549 ),
550 U_,
551 zeroGradientFvPatchScalarField::typeName
552 ),
553 Uw_
554 (
556 (
557 "Uwf",
558 time().timeName(),
559 regionMesh(),
560 IOobject::NO_READ,
561 IOobject::NO_WRITE
562 ),
563 U_,
564 zeroGradientFvPatchScalarField::typeName
565 ),
566 deltaRho_
567 (
569 (
570 delta_.name() + "*" + rho_.name(),
571 time().timeName(),
572 regionMesh(),
573 IOobject::NO_READ,
574 IOobject::NO_WRITE
575 ),
576 regionMesh(),
577 dimensionedScalar(delta_.dimensions()*rho_.dimensions(), Zero),
578 zeroGradientFvPatchScalarField::typeName
579 ),
580
581 phi_
582 (
584 (
585 "phi",
586 time().timeName(),
587 regionMesh(),
588 IOobject::NO_READ,
589 IOobject::AUTO_WRITE
590 ),
591 regionMesh(),
593 ),
594
595 primaryMassTrans_
596 (
598 (
599 "primaryMassTrans",
600 time().timeName(),
601 regionMesh(),
602 IOobject::NO_READ,
603 IOobject::NO_WRITE
604 ),
605 regionMesh(),
607 zeroGradientFvPatchScalarField::typeName
608 ),
609 cloudMassTrans_
610 (
612 (
613 "cloudMassTrans",
614 time().timeName(),
615 regionMesh(),
616 IOobject::NO_READ,
617 IOobject::NO_WRITE
618 ),
619 regionMesh(),
621 zeroGradientFvPatchScalarField::typeName
622 ),
623 cloudDiameterTrans_
624 (
626 (
627 "cloudDiameterTrans",
628 time().timeName(),
629 regionMesh(),
630 IOobject::NO_READ,
631 IOobject::NO_WRITE
632 ),
633 regionMesh(),
634 dimensionedScalar("minus1", dimLength, -1.0),
635 zeroGradientFvPatchScalarField::typeName
636 ),
637
638 USp_
639 (
641 (
642 "USpf",
643 time().timeName(),
644 regionMesh(),
645 IOobject::NO_READ,
646 IOobject::NO_WRITE
647 ),
648 regionMesh(),
650 this->mappedPushedFieldPatchTypes<vector>()
651 ),
652 pSp_
653 (
655 (
656 "pSpf",
657 time_.timeName(),
658 regionMesh(),
659 IOobject::NO_READ,
660 IOobject::NO_WRITE
661 ),
662 regionMesh(),
664 this->mappedPushedFieldPatchTypes<scalar>()
665 ),
666 rhoSp_
667 (
669 (
670 "rhoSpf",
671 time_.timeName(),
672 regionMesh(),
673 IOobject::NO_READ,
674 IOobject::NO_WRITE
675 ),
676 regionMesh(),
678 this->mappedPushedFieldPatchTypes<scalar>()
679 ),
680
681 USpPrimary_
682 (
684 (
685 USp_.name(), // must have same name as USp_ to enable mapping
686 time().timeName(),
687 primaryMesh(),
688 IOobject::NO_READ,
689 IOobject::NO_WRITE
690 ),
691 primaryMesh(),
692 dimensionedVector(USp_.dimensions(), Zero)
693 ),
694 pSpPrimary_
695 (
697 (
698 pSp_.name(), // must have same name as pSp_ to enable mapping
699 time().timeName(),
700 primaryMesh(),
701 IOobject::NO_READ,
702 IOobject::NO_WRITE
703 ),
704 primaryMesh(),
705 dimensionedScalar(pSp_.dimensions(), Zero)
706 ),
707 rhoSpPrimary_
708 (
710 (
711 rhoSp_.name(), // must have same name as rhoSp_ to enable mapping
712 time().timeName(),
713 primaryMesh(),
714 IOobject::NO_READ,
715 IOobject::NO_WRITE
716 ),
717 primaryMesh(),
718 dimensionedScalar(rhoSp_.dimensions(), Zero)
719 ),
720
721 UPrimary_
722 (
724 (
725 "U", // must have same name as U to enable mapping
726 time().timeName(),
727 regionMesh(),
728 IOobject::NO_READ,
729 IOobject::NO_WRITE
730 ),
731 regionMesh(),
733 this->mappedFieldAndInternalPatchTypes<vector>()
734 ),
735 pPrimary_
736 (
738 (
739 "p", // must have same name as p to enable mapping
740 time().timeName(),
741 regionMesh(),
742 IOobject::NO_READ,
743 IOobject::NO_WRITE
744 ),
745 regionMesh(),
747 this->mappedFieldAndInternalPatchTypes<scalar>()
748 ),
749 rhoPrimary_
750 (
752 (
753 "rho", // must have same name as rho to enable mapping
754 time().timeName(),
755 regionMesh(),
756 IOobject::NO_READ,
757 IOobject::NO_WRITE
758 ),
759 regionMesh(),
761 this->mappedFieldAndInternalPatchTypes<scalar>()
762 ),
763 muPrimary_
764 (
766 (
767 "thermo:mu", // must have same name as mu to enable mapping
768 time().timeName(),
769 regionMesh(),
770 IOobject::NO_READ,
771 IOobject::NO_WRITE
772 ),
773 regionMesh(),
775 this->mappedFieldAndInternalPatchTypes<scalar>()
776 ),
777
778 filmThermo_(filmThermoModel::New(*this, coeffs_)),
779
780 availableMass_(regionMesh().nCells(), Zero),
781
782 injection_(*this, coeffs_),
783
784 transfer_(*this, coeffs_),
785
786 turbulence_(filmTurbulenceModel::New(*this, coeffs_)),
787
788 forces_(*this, coeffs_),
789
790 addedMassTotal_(0.0)
791{
792 if (readFields)
793 {
795
796 correctAlpha();
797
799
801
803 (
805 (
806 "phi",
807 time().timeName(),
808 regionMesh(),
811 false
812 ),
814 );
815
816 phi_ == phi0;
817 }
818}
819
820
821// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
822
824{}
825
826
827// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
828
830(
831 const label patchi,
832 const label facei,
833 const scalar massSource,
834 const vector& momentumSource,
835 const scalar pressureSource,
836 const scalar energySource
837)
838{
840 << "\nSurface film: " << type() << ": adding to film source:" << nl
841 << " mass = " << massSource << nl
842 << " momentum = " << momentumSource << nl
843 << " pressure = " << pressureSource << endl;
844
845 rhoSpPrimary_.boundaryFieldRef()[patchi][facei] -= massSource;
846 USpPrimary_.boundaryFieldRef()[patchi][facei] -= momentumSource;
847 pSpPrimary_.boundaryFieldRef()[patchi][facei] -= pressureSource;
848
849 addedMassTotal_ += massSource;
850}
851
852
854{
856
858
860
862
864
866
867 correctAlpha();
868
869 // Reset transfer fields
874}
875
876
878{
880
881 // Update sub-models to provide updated source contributions
883
884 // Solve continuity for deltaRho_
886
887 // Implicit pressure source coefficient - constant
888 tmp<volScalarField> tpp(this->pp());
889
890 for (int oCorr=1; oCorr<=nOuterCorr_; oCorr++)
891 {
892 // Explicit pressure source contribution - varies with delta_
893 tmp<volScalarField> tpu(this->pu());
894
895 // Solve for momentum for U_
897 fvVectorMatrix& UEqn = tUEqn.ref();
898
899 // Film thickness correction loop
900 for (int corr=1; corr<=nCorr_; corr++)
901 {
902 // Solve thickness for delta_
903 solveThickness(tpu(), tpp(), UEqn);
904 }
905 }
906
907 // Update deltaRho_ with new delta_
909}
910
911
913{
915
916 // Reset source terms for next time integration
918}
919
920
922{
923 scalar CoNum = 0.0;
924
925 if (regionMesh().nInternalFaces() > 0)
926 {
927 const scalarField sumPhi
928 (
929 fvc::surfaceSum(mag(phi_))().primitiveField()
930 / (deltaRho_.primitiveField() + ROOTVSMALL)
931 );
932
933 forAll(delta_, i)
934 {
935 if ((delta_[i] > deltaCoLimit_) && (alpha_[i] > 0.5))
936 {
937 CoNum = max(CoNum, sumPhi[i]/(delta_[i]*magSf()[i]));
938 }
939 }
940
941 CoNum *= 0.5*time_.deltaTValue();
942 }
943
945
946 Info<< "Film max Courant number: " << CoNum << endl;
947
948 return CoNum;
949}
950
951
953{
954 return U_;
955}
956
957
959{
960 return Us_;
961}
962
963
965{
966 return Uw_;
967}
968
969
971{
972 return deltaRho_;
973}
974
975
977{
978 return phi_;
979}
980
981
983{
984 return rho_;
985}
986
987
989{
991 << "T field not available for " << type() << abort(FatalError);
992
993 return volScalarField::null();
994}
995
996
998{
1000 << "Ts field not available for " << type() << abort(FatalError);
1001
1002 return volScalarField::null();
1003}
1004
1005
1007{
1009 << "Tw field not available for " << type() << abort(FatalError);
1010
1011 return volScalarField::null();
1012}
1013
1014
1016{
1018 << "hs field not available for " << type() << abort(FatalError);
1019
1020 return volScalarField::null();
1021}
1022
1023
1025{
1027 << "Cp field not available for " << type() << abort(FatalError);
1028
1029 return volScalarField::null();
1030}
1031
1032
1034{
1036 << "kappa field not available for " << type() << abort(FatalError);
1037
1038 return volScalarField::null();
1039}
1040
1041
1043{
1044 return primaryMassTrans_;
1045}
1046
1047
1049{
1050 return cloudMassTrans_;
1051}
1052
1053
1055{
1056 return cloudDiameterTrans_;
1057}
1058
1059
1061{
1062 Info<< "\nSurface film: " << type() << endl;
1063
1064 const scalarField& deltaInternal = delta_;
1065 const vectorField& Uinternal = U_;
1066 scalar addedMassTotal = 0.0;
1067 outputProperties().readIfPresent("addedMassTotal", addedMassTotal);
1068 addedMassTotal += returnReduce(addedMassTotal_, sumOp<scalar>());
1069
1070 Info<< indent << "added mass = " << addedMassTotal << nl
1071 << indent << "current mass = "
1072 << gSum((deltaRho_*magSf())()) << nl
1073 << indent << "min/max(mag(U)) = " << gMin(mag(Uinternal)) << ", "
1074 << gMax(mag(Uinternal)) << nl
1075 << indent << "min/max(delta) = " << gMin(deltaInternal) << ", "
1076 << gMax(deltaInternal) << nl
1077 << indent << "coverage = "
1078 << gSum(alpha_.primitiveField()*magSf())/gSum(magSf()) << nl;
1079
1082}
1083
1084
1086{
1088 (
1089 IOobject
1090 (
1091 typeName + ":Srho",
1092 time().timeName(),
1093 primaryMesh(),
1096 false
1097 ),
1098 primaryMesh(),
1100 );
1101}
1102
1103
1105(
1106 const label i
1107) const
1108{
1110 (
1111 IOobject
1112 (
1113 typeName + ":Srho(" + Foam::name(i) + ")",
1114 time().timeName(),
1115 primaryMesh(),
1118 false
1119 ),
1120 primaryMesh(),
1122 );
1123}
1124
1125
1127{
1129 (
1130 IOobject
1131 (
1132 typeName + ":Sh",
1133 time().timeName(),
1134 primaryMesh(),
1137 false
1138 ),
1139 primaryMesh(),
1141 );
1142}
1143
1144
1145// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1146
1147} // End namespace surfaceFilmModels
1148} // End namespace regionModels
1149} // End namespace Foam
1150
1151// ************************************************************************* //
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
Definition: continuityErrs.H:4
const uniformDimensionedVectorField & g
const dimensionSet & dimensions() const
Return dimensions.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static const GeometricField< scalar, fvPatchField, volMesh > & null()
Return a null geometric field.
void correctBoundaryConditions()
Correct boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:43
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:55
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
A List with indirect addressing. Like IndirectList but does not store addressing.
Definition: IndirectList.H:79
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:640
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1092
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:1443
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:1303
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:1332
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
Lookup type of boundary radiation properties.
Definition: lookup.H:66
virtual bool write(const bool valid=true) const
Write using setting from DB.
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
tmp< faVectorMatrix > correct(areaVectorField &U)
Return (net) force system.
Definition: forceList.C:82
virtual void correct(scalarField &availableMass, volScalarField &massToInject, volScalarField &diameterToInject)
Correct.
const Time & time_
Reference to the time database.
Definition: regionModel.H:90
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:40
const dictionary & solution() const
Return the solution dictionary.
Definition: regionModelI.H:99
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionModelI.H:34
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
const IOdictionary & outputProperties() const
Return const access to the output properties dictionary.
Definition: regionModelI.H:106
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
Definition: regionModel.H:114
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
virtual const volVectorField & nHat() const
Return the patch normal vectors.
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 void solveContinuity()
Solve continuity equation.
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.
label nNonOrthCorr_
Number of non-orthogonal correctors.
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.
virtual bool read()
Read control parameters from dictionary.
virtual void correct(scalarField &availableMass, volScalarField &massToTransfer)
Correct kinematic transfers.
virtual void info(Ostream &os)
Provide some info.
void setFluxRequired(const word &name) const
Get flux-required for given name, or default.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:66
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
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:453
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)"))
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:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:34
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)
Definition: fvcLaplacian.C:47
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Namespace for OpenFOAM.
const dimensionSet dimPressure
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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:47
Type gSum(const FieldField< Field, Type > &f)
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimEnergy
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:342
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
Type gMin(const FieldField< Field, Type > &f)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
const dimensionSet dimDensity
Type gMax(const FieldField< Field, Type > &f)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
CEqn solve()
volScalarField & e
Definition: createFields.H:11
scalar sumLocalContErr
scalar globalContErr
scalar CoNum
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333