kinematicSingleLayer.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "kinematicSingleLayer.H"
27 
28 #include "fvcDdt.H"
29 #include "fvcDiv.H"
30 #include "fvcLaplacian.H"
31 #include "fvcSnGrad.H"
32 #include "fvcReconstruct.H"
33 #include "fvcVolumeIntegrate.H"
34 #include "fvcFlux.H"
35 
36 #include "fvmDdt.H"
37 #include "fvmDiv.H"
38 #include "fvmLaplacian.H"
39 #include "fvmSup.H"
40 #include "constrainHbyA.H"
41 
42 #include "mappedWallPolyPatch.H"
43 #include "distributionMap.H"
44 #include "filmViscosityModel.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 namespace regionModels
51 {
52 namespace surfaceFilmModels
53 {
54 
55 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
56 
57 defineTypeNameAndDebug(kinematicSingleLayer, 0);
58 
59 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
60 
62 {
64 }
65 
66 
68 {
70 
72  USpPrimary_ == Zero;
73  pSpPrimary_ == Zero;
74 }
75 
76 
78 {
80 
81  // Update fields from primary region via direct mapped
82  // (coupled) boundary conditions
87 }
88 
89 
91 {
93 
94  volScalarField::Boundary& rhoSpPrimaryBf =
96 
97  volVectorField::Boundary& USpPrimaryBf =
99 
100  volScalarField::Boundary& pSpPrimaryBf =
102 
103  // Convert accumulated source terms into per unit area per unit time
104  const scalar deltaT = time_.deltaTValue();
106  {
107  scalarField rpriMagSfdeltaT
108  (
109  (1/deltaT)
110  /primaryMesh().magSf().boundaryField()[patchi]
111  );
112 
113  rhoSpPrimaryBf[patchi] *= rpriMagSfdeltaT;
114  USpPrimaryBf[patchi] *= rpriMagSfdeltaT;
115  pSpPrimaryBf[patchi] *= rpriMagSfdeltaT;
116  }
117 
118  // Retrieve the source fields from the primary region
119  toRegion(rhoSp_, rhoSpPrimaryBf);
120  rhoSp_.field() /= VbyA();
121  toRegion(USp_, USpPrimaryBf);
122  USp_.field() /= VbyA();
123  toRegion(pSp_, pSpPrimaryBf);
124 
125  // update addedMassTotal counter
126  if (time().writeTime())
127  {
128  scalar addedMassTotal = 0;
129  outputProperties().readIfPresent("addedMassTotal", addedMassTotal);
130  addedMassTotal += returnReduce(addedMassTotal_, sumOp<scalar>());
131  outputProperties().add("addedMassTotal", addedMassTotal, true);
132  addedMassTotal_ = 0;
133  }
134 }
135 
136 
138 {
139  return -fvc::laplacian(sigma(), delta_);
140 }
141 
142 
144 {
146  (
148  (
149  "pSp",
150  regionMesh(),
152  zeroGradientFvPatchScalarField::typeName
153  )
154  );
155 
156  tpSp.ref().primitiveFieldRef() = pSp_;
157  tpSp.ref().correctBoundaryConditions();
158 
159  return volScalarField::New
160  (
161  IOobject::modelName("pe", typeName),
162  p_ // Pressure (mapped from primary region)
163  - tpSp // Accumulated particle impingement
164  );
165 }
166 
167 
169 {
170  return
172  (
173  max(nHat() & -g(), dimensionedScalar(g().dimensions(), 0))*VbyA()
174  )*fvc::interpolate(rho());
175 }
176 
177 
179 {
180  return
182  (
183  max(nHat() & -g(), dimensionedScalar(g().dimensions(), 0))*VbyA()
184  )*fvc::snGrad(rho());
185 }
186 
187 
189 {
191 }
192 
193 
195 {
197 
198  // Update ejection model - mass returned is mass available for ejection
200 
201  // Update transfer model - mass returned is mass available for transfer
203 
204  const volScalarField::Internal rVDt
205  (
206  1/(time().deltaT()*regionMesh().V())
207  );
208 
209  // Update mass source field
210  rhoSp_ += rVDt*(cloudMassTrans_() + primaryMassTrans_());
211  USp_ += rVDt*(cloudMassTrans_()*U_() + primaryMomentumTrans_());
212 
213  momentumTransport_->correct();
214 }
215 
216 
218 {
220 
222 
223  // Bound film volume fraction
224  alpha_.max(0);
225 
226  delta_ == alpha_*VbyA();
227 
228  // Update continuity error caused by the delta_ bounding
230 }
231 
232 
234 {
236 }
237 
238 
240 {
242 
243  if (totalMass.value() > small)
244  {
245  const volScalarField::Internal massErr
246  (
248  );
249 
250  const scalar sumLocalContErr =
251  (
253  ).value();
254 
255  const scalar globalContErr =
256  (
258  ).value();
259 
261 
262  Info<< "time step continuity errors: sum local = "
263  << sumLocalContErr
264  << ", global = " << globalContErr
265  << ", cumulative = " << cumulativeContErr_
266  << endl;
267  }
268 }
269 
270 
272 {
274  (
276  (
277  "Uw",
278  regionMesh(),
280  )
281  );
282 
284 
285  // Push boundary film velocity values into internal field
286  for (label i=0; i<intCoupledPatchIDs_.size(); i++)
287  {
288  const label patchi = intCoupledPatchIDs_[i];
289  const polyPatch& pp = regionMesh().boundaryMesh()[patchi];
292  }
293 
294  Uw -= nHat()*(Uw_ & nHat());
295 
296  return tUw;
297 }
298 
299 
301 (
302  const volScalarField& pc,
303  const volScalarField& pe
304 )
305 {
307 
308  // Evaluate viscosity from user-model
309  viscosity_->correct(thermo_->p(), thermo_->T());
310 
311  const volScalarField::Internal rVDt
312  (
313  1/(time().deltaT()*regionMesh().V())
314  );
315 
316  // Momentum equation
318  (
319  fvm::ddt(alpha_, rho(), U_) + fvm::div(phi_, U_)
321  ==
322  - USp_
323  + forces_.correct(U_)
324  + momentumTransport_->Su(U_)
325  );
326 
327  fvVectorMatrix& UEqn = tUEqn.ref();
328 
329  UEqn.relax();
330 
332  {
334 
335  solve
336  (
337  UEqn
338  ==
340  (
342  (
343  alphaf
344  *(
345  (
346  fvc::snGrad(pe + pc, "snGrad(p)")
347  + gGradRho()*alphaf
348  + rhog()*fvc::snGrad(alpha_)
349  )*regionMesh().magSf()
350  - fvc::interpolate(rho())*(g() & regionMesh().Sf())
351  ), 0
352  )
353  )
354  );
355 
356  // Remove any patch-normal components of velocity
357  U_ -= nHat()*(nHat() & U_);
358 
360  }
361 
362  return tUEqn;
363 }
364 
365 
367 (
368  const fvVectorMatrix& UEqn,
369  const volScalarField& pc,
370  const volScalarField& pe
371 )
372 {
374 
375  const volScalarField rAU(1/UEqn.A());
376  const volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U_, alpha_));
377 
380  const surfaceScalarField alpharAUf(fvc::interpolate(alpha_*rAU));
381  const surfaceScalarField rhogf(rhog());
382 
383  const surfaceScalarField phiu
384  (
385  "phiu",
386  (
388  (
389  (
390  fvc::snGrad(pe + pc, "snGrad(p)")
391  + gGradRho()*alphaf
392  )*regionMesh().magSf()
393  - rhof*(g() & regionMesh().Sf()),
394  0
395  )
396  )
397  );
398 
399  surfaceScalarField phid
400  (
401  "phid",
402  // constrainFilmField
403  // (
404  // rhof
405  // *constrainPhiHbyA(fvc::flux(HbyA) - alpharAUf*phiu, U_, alpha_),
406  // 0
407  // )
408  rhof*constrainPhiHbyA(fvc::flux(HbyA) - alpharAUf*phiu, U_, alpha_)
409  );
410 
411  const surfaceScalarField ddrhorAUrhogf
412  (
413  "alphaCoeff",
414  alphaf*rhof*alpharAUf*rhogf
415  );
416 
418 
419  while (pimple_.correctNonOrthogonal())
420  {
421  // Film thickness equation
422  fvScalarMatrix alphaEqn
423  (
424  fvm::ddt(rho(), alpha_)
425  + fvm::div(phid, alpha_)
426  - fvm::laplacian(ddrhorAUrhogf, alpha_)
427  ==
428  -rhoSp_
429  );
430 
431  alphaEqn.solve();
432 
434  {
435  phi_ == alphaEqn.flux();
436 
437  const surfaceScalarField phiGradAlpha
438  (
440  (
441  rhogf*fvc::snGrad(alpha_)*regionMesh().magSf(),
442  0
443  )
444  );
445 
446  phiU_ = constrainFilmField(phid/rhof - alpharAUf*phiGradAlpha, 0);
447 
448  // Update U field
449  U_ = HbyA - rAU*fvc::reconstruct(alphaf*(phiu + phiGradAlpha));
450 
451  // Remove any patch-normal components of velocity
452  U_ -= nHat()*(nHat() & U_);
453 
455  }
456  }
457 
458  // Bound film volume fraction
459  alpha_.max(0);
460 
461  delta_ == alpha_*VbyA();
462 
464 
465  // Continuity check
466  continuityCheck();
467 }
468 
469 
470 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
471 
473 (
474  const word& modelType,
475  const fvMesh& mesh,
476  const dimensionedVector& g,
477  const word& regionType,
478  const bool readFields
479 )
480 :
481  surfaceFilmRegionModel(modelType, mesh, g, regionType),
483  pimple_(regionMesh()),
484 
486 
487  deltaSmall_("deltaSmall", dimLength, small),
488  maxCo_(solution().lookupOrDefault<scalar>("maxCo", 0)),
489 
490  p_
491  (
492  IOobject
493  (
494  "p",
495  time().timeName(),
496  regionMesh()
497  ),
498  regionMesh(),
500  this->mappedFieldAndInternalPatchTypes<scalar>()
501  ),
502 
504 
505  mu_
506  (
507  IOobject
508  (
509  "mu",
510  time().timeName(),
511  regionMesh(),
514  ),
515  regionMesh(),
516  dimensionedScalar(dimPressure*dimTime, 0),
517  zeroGradientFvPatchScalarField::typeName
518  ),
519 
520  delta_
521  (
522  IOobject
523  (
524  "delta",
525  time().timeName(),
526  regionMesh(),
529  ),
530  regionMesh()
531  ),
532 
533  alpha_
534  (
535  IOobject
536  (
537  "alpha",
538  time().timeName(),
539  regionMesh(),
542  ),
543  delta_/VbyA(),
545  ),
546 
547  U_
548  (
549  IOobject
550  (
551  "U",
552  time().timeName(),
553  regionMesh(),
556  ),
557  regionMesh()
558  ),
559 
560  Uw_
561  (
562  IOobject
563  (
564  "Uw",
565  time().timeName(),
566  regionMesh(),
569  ),
570  U_
571  ),
572 
573  phi_
574  (
575  IOobject
576  (
577  "phi",
578  time().timeName(),
579  regionMesh(),
582  ),
583  regionMesh(),
584  dimensionedScalar(dimMass/dimTime, 0)
585  ),
586 
587  phiU_
588  (
589  IOobject
590  (
591  "phiU",
592  time().timeName(),
593  regionMesh()
594  ),
595  regionMesh(),
596  dimensionedScalar(dimVolume/dimTime, 0)
597  ),
598 
600  (
601  IOobject
602  (
603  "continuityErr",
604  time().timeName(),
605  regionMesh()
606  ),
607  regionMesh(),
609  ),
610 
611  coverage_
612  (
613  IOobject
614  (
615  "coverage",
616  time().timeName(),
617  regionMesh(),
620  ),
621  regionMesh(),
623  zeroGradientFvPatchScalarField::typeName
624  ),
625 
627  (
628  IOobject
629  (
630  "primaryMassTrans",
631  time().timeName(),
632  regionMesh(),
635  ),
636  regionMesh(),
638  zeroGradientFvPatchScalarField::typeName
639  ),
640 
642  (
643  IOobject
644  (
645  "cloudMassTrans",
646  time().timeName(),
647  regionMesh(),
650  ),
651  regionMesh(),
653  zeroGradientFvPatchScalarField::typeName
654  ),
655 
657  (
658  IOobject
659  (
660  "cloudDiameterTrans",
661  time().timeName(),
662  regionMesh(),
665  ),
666  regionMesh(),
668  zeroGradientFvPatchScalarField::typeName
669  ),
670 
672  (
673  IOobject
674  (
675  "primaryMomentumTrans",
676  time().timeName(),
677  regionMesh(),
680  ),
681  regionMesh(),
683  zeroGradientFvPatchVectorField::typeName
684  ),
685 
686  rhoSp_
687  (
688  IOobject
689  (
690  "rhoSp",
691  time_.timeName(),
692  regionMesh(),
695  ),
696  regionMesh(),
697  dimensionedScalar(dimDensity/dimTime, 0)
698  ),
699 
700  USp_
701  (
702  IOobject
703  (
704  "USp",
705  time().timeName(),
706  regionMesh(),
709  ),
710  regionMesh(),
712  ),
713 
714  pSp_
715  (
716  IOobject
717  (
718  "pSp",
719  time_.timeName(),
720  regionMesh(),
723  ),
724  regionMesh(),
726  ),
727 
729  (
730  IOobject
731  (
732  rhoSp_.name(),
733  time().timeName(),
734  primaryMesh(),
737  ),
738  primaryMesh(),
740  ),
741 
743  (
744  IOobject
745  (
746  USp_.name(),
747  time().timeName(),
748  primaryMesh(),
751  ),
752  primaryMesh(),
754  ),
755 
757  (
758  IOobject
759  (
760  pSp_.name(),
761  time().timeName(),
762  primaryMesh(),
765  ),
766  primaryMesh(),
768  ),
769 
770  UPrimary_
771  (
772  IOobject
773  (
774  "U", // Must have same name as U to enable mapping
775  time().timeName(),
776  regionMesh(),
779  ),
780  regionMesh(),
782  this->mappedFieldAndInternalPatchTypes<vector>()
783  ),
784 
786  (
787  IOobject
788  (
789  // Must have same name as rho to enable mapping
790  IOobject::groupName("thermo:rho", phaseName_),
791  time().timeName(),
792  regionMesh(),
795  ),
796  regionMesh(),
798  this->mappedFieldAndInternalPatchTypes<scalar>()
799  ),
800 
801  muPrimary_
802  (
803  IOobject
804  (
805  // Must have same name as rho to enable mapping
806  IOobject::groupName("thermo:mu", phaseName_),
807  time().timeName(),
808  regionMesh(),
811  ),
812  regionMesh(),
813  dimensionedScalar(dimPressure*dimTime, 0),
814  this->mappedFieldAndInternalPatchTypes<scalar>()
815  ),
816 
818 
819  sigma_(Function1<scalar>::New("sigma", coeffs())),
820 
821  availableMass_(regionMesh().nCells(), 0),
822 
823  ejection_(*this, coeffs_),
824 
825  transfer_(*this, coeffs_),
826 
828 
829  forces_(*this, coeffs_),
830 
831  addedMassTotal_(0)
832 {
833  alpha_ == delta_/VbyA();
834 
835  if (readFields)
836  {
838 
839  correctCoverage();
840 
842  (
843  IOobject
844  (
845  "phi",
846  time().timeName(),
847  regionMesh(),
850  false
851  ),
853  );
854 
855  phi_ == phi;
856  phiU_ = fvc::flux(U_);
857  }
858 }
859 
860 
861 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
862 
864 {}
865 
866 
867 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
868 
870 {
871  tmp<volScalarField> tsigma
872  (
874  (
875  type() + ":sigma",
876  regionMesh(),
877  dimensionedScalar(dimMass/sqr(dimTime), 0),
878  extrapolatedCalculatedFvPatchScalarField::typeName
879  )
880  );
881 
882  tsigma.ref().primitiveFieldRef() = sigma_->value(thermo_->T());
883 
884  tsigma.ref().correctBoundaryConditions();
885 
886  return tsigma;
887 }
888 
889 
891 (
892  const label patchi,
893  const label facei,
894  const scalar massSource,
895  const vector& momentumSource,
896  const scalar pressureSource,
897  const scalar energySource
898 )
899 {
901  << "\nSurface film: " << type() << ": adding to film source:" << nl
902  << " mass = " << massSource << nl
903  << " momentum = " << momentumSource << nl
904  << " pressure = " << pressureSource << endl;
905 
906  rhoSpPrimary_.boundaryFieldRef()[patchi][facei] -= massSource;
907  USpPrimary_.boundaryFieldRef()[patchi][facei] -= momentumSource;
908  pSpPrimary_.boundaryFieldRef()[patchi][facei] -= pressureSource;
909 
910  addedMassTotal_ += massSource;
911 }
912 
913 
915 {
917 
919 
921 
923 
924  // Reset transfer fields
925  availableMass_ = mass();
930 }
931 
932 
934 {
936 
937  // Update film coverage indicator
938  correctCoverage();
939 
940  // Predict delta_ from continuity
941  predictDelta();
942 
943  // Update sub-models to provide updated source contributions
944  updateSubmodels();
945 
946  // Predict delta_ from continuity with updated source
947  predictDelta();
948 
949  // Capillary pressure
950  const volScalarField pc(this->pc());
951 
952  while (pimple_.loop())
953  {
954  // External pressure
955  const volScalarField pe(this->pe());
956 
957  // Solve for momentum
958  const fvVectorMatrix UEqn(solveMomentum(pc, pe));
959 
960  // Film thickness correction loop
961  while (pimple_.correct())
962  {
963  solveAlpha(UEqn, pc, pe);
964  }
965  }
966 
967  // Reset source terms for next time integration
969 }
970 
971 
973 {
974  const scalarField sumPhi(fvc::surfaceSum(mag(phiU_))().primitiveField());
975 
976  return 0.5*gMax(sumPhi/regionMesh().V().field())*time_.deltaTValue();
977 }
978 
979 
981 {
982  if (maxCo_ > 0)
983  {
984  return maxCo_*time_.deltaTValue()/(CourantNumber() + small);
985  }
986  else
987  {
988  return great;
989  }
990 }
991 
992 
994 {
995  return primaryMassTrans_;
996 }
997 
998 
1000 {
1001  return cloudMassTrans_;
1002 }
1003 
1004 
1006 {
1007  return cloudDiameterTrans_;
1008 }
1009 
1010 
1012 {
1013  return primaryMomentumTrans_;
1014 }
1015 
1016 
1018 {
1019  Info<< "\nSurface film: " << type() << endl;
1020 
1021  const scalarField& deltaInternal = delta_;
1022  const vectorField& Uinternal = U_;
1023  scalar addedMassTotal = 0;
1024  outputProperties().readIfPresent("addedMassTotal", addedMassTotal);
1025  addedMassTotal += returnReduce(addedMassTotal_, sumOp<scalar>());
1026 
1027  Info<< indent << "added mass = " << addedMassTotal << nl
1028  << indent << "current mass = "
1029  << gSum((delta_*rho()*magSf())()) << nl
1030  << indent << "min/max(mag(U)) = " << gMin(mag(Uinternal)) << ", "
1031  << gMax(mag(Uinternal)) << nl
1032  << indent << "max Courant number = " << CourantNumber() << nl
1033  << indent << "min/max(delta) = " << gMin(deltaInternal) << ", "
1034  << gMax(deltaInternal) << nl
1035  << indent << "coverage = "
1036  << gSum(coverage_.primitiveField()*magSf())/gSum(magSf()) << nl;
1037 
1038  ejection_.info(Info);
1039  transfer_.info(Info);
1040 }
1041 
1042 
1044 {
1046  (
1048  (
1049  "thermoSingleLayer::Srho",
1050  primaryMesh(),
1051  dimensionedScalar(dimMass/dimVolume/dimTime, 0)
1052  )
1053  );
1054 
1055  scalarField& Srho = tSrho.ref();
1056  const scalarField& V = primaryMesh().V();
1057  const scalar dt = time_.deltaTValue();
1058 
1060  {
1061  const label filmPatchi = intCoupledPatchIDs()[i];
1062 
1063  scalarField patchMass =
1064  primaryMassTrans_.boundaryField()[filmPatchi];
1065 
1066  toPrimary(filmPatchi, patchMass);
1067 
1068  const label primaryPatchi = primaryPatchIDs()[i];
1069  const unallocLabelList& cells =
1070  primaryMesh().boundaryMesh()[primaryPatchi].faceCells();
1071 
1072  forAll(patchMass, j)
1073  {
1074  Srho[cells[j]] += patchMass[j]/(V[cells[j]]*dt);
1075  }
1076  }
1077 
1078  return tSrho;
1079 }
1080 
1081 
1084  const label i
1085 ) const
1086 {
1088  (
1089  IOobject::modelName("SY(" + Foam::name(i) + ")", typeName),
1090  primaryMesh(),
1091  dimensionedScalar(dimMass/dimVolume/dimTime, 0)
1092  );
1093 }
1094 
1095 
1097 {
1099  (
1101  (
1102  IOobject::modelName("SU", typeName),
1103  primaryMesh(),
1105  )
1106  );
1107 
1108  vectorField& SU = tSU.ref();
1109  const scalarField& V = primaryMesh().V();
1110  const scalar dt = time_.deltaTValue();
1111 
1113  {
1114  const label filmPatchi = intCoupledPatchIDs_[i];
1115 
1116  vectorField patchMomentum =
1117  primaryMomentumTrans_.boundaryField()[filmPatchi];
1118 
1119  toPrimary(filmPatchi, patchMomentum);
1120 
1121  const unallocLabelList& cells =
1122  primaryMesh().boundaryMesh()[primaryPatchIDs()[i]].faceCells();
1123 
1124  forAll(patchMomentum, j)
1125  {
1126  SU[cells[j]] += patchMomentum[j]/(V[cells[j]]*dt);
1127  }
1128  }
1129 
1130  return tSU;
1131 }
1132 
1133 
1135 {
1137  (
1138  IOobject::modelName("Sh", typeName),
1139  primaryMesh(),
1141  );
1142 }
1143 
1144 
1145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1146 
1147 } // End namespace surfaceFilmModels
1148 } // End namespace regionModels
1149 } // End namespace Foam
1150 
1151 // ************************************************************************* //
virtual bool read()
Read control parameters from dictionary.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
virtual bool read()
Read control parameters from dictionary.
tmp< volScalarField::Internal > mass() const
Return the current film mass.
Run-time selectable general function of one variable.
Definition: Function1.H:52
bool finalNonOrthogonalIter() const
Flag to indicate the last non-orthogonal iteration.
surfaceFilmRegionModel(const word &modelType, const fvMesh &mesh, const dimensionedVector &g, const word &regionType)
Construct from type name, mesh and gravity vector.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const volScalarField & VbyA() const
Return the cell layer volume/area [m].
virtual scalar maxDeltaT() const
Return the maximum time-step for stable operation.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const word & name() const
Return name.
Definition: IOobject.H:315
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
bool correctNonOrthogonal()
Non-orthogonal corrector loop.
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1672
static autoPtr< momentumTransportModel > New(surfaceFilmRegionModel &film, const dictionary &dict)
Return a reference to the selected ejection model.
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
volScalarField primaryMassTrans_
Film mass available for transfer to the primary region.
volVectorField primaryMomentumTrans_
Film momentum transfer.
autoPtr< Function1< scalar > > sigma_
Surface tension function.
Type gMin(const FieldField< Field, Type > &f)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
volScalarField::Internal continuityErr_
Current continuity error caused by delta_ bounding.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
const dimensionSet dimPressure
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static autoPtr< viscosityModel > New(surfaceFilmRegionModel &film, const dictionary &dict, volScalarField &mu)
Return a reference to the selected phase change model.
Reconstruct volField from a face flux field.
wordList types() const
Return a list of the patch field types.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
virtual void correctCoverage()
Correct film coverage field.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:831
tmp< surfaceScalarField > constrainPhiHbyA(const tmp< surfaceScalarField > &tphiHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:79
Calculate the matrix for the laplacian of the field.
static autoPtr< rhoThermo > New(const fvMesh &, const word &phaseName=word::null)
Standard selection based on fvMesh.
Definition: rhoThermo.C:92
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:487
volVectorField::Internal USp_
Momentum [kg/m/s^2].
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
virtual tmp< volScalarField::Internal > Srho() const
Return total mass source - Eulerian phase only.
Calculate the snGrad of the given volField.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
scalar globalContErr
scalarField availableMass_
Available mass for transfer via sub-models.
volScalarField cloudDiameterTrans_
Parcel diameters originating from film to cloud.
const dimensionSet dimless
volScalarField coverage_
Film coverage indicator, 1 = covered, 0 = uncovered [].
scalar addedMassTotal_
Cumulative mass added via sources [kg].
volVectorField USpPrimary_
Primary region tangential momentum source [kg m/s].
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
virtual tmp< volScalarField > primaryMassTrans() const
Return mass transfer source - Eulerian phase only.
virtual tmp< volScalarField::Internal > SYi(const label i) const
Return mass source for specie i - Eulerian phase only.
tmp< surfaceScalarField > rhog() const
Hydrostatic pressure coefficient.
virtual void evolveRegion()
Evolve the film equations.
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1153
const dimensionSet dimLength
tmp< fvVectorMatrix > correct(volVectorField &U)
Return (net) force system.
Definition: forceList.C:103
virtual tmp< fvVectorMatrix > solveMomentum(const volScalarField &pc, const volScalarField &pe)
Solve for film velocity.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
virtual void transferPrimaryRegionSourceFields()
Transfer source fields from the primary region to the film region.
Calculate the first temporal derivative.
tmp< volVectorField::Internal > Uw() const
Return the film wall velocity [m/s].
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet dimTime
const dimensionSet & dimensions() const
Return dimensions.
void toPrimary(const label regionPatchi, List< Type > &regionField) const
Convert a local region field to the primary region.
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionModelI.H:31
virtual void correct(scalarField &availableMass, volScalarField &massToTransfer, volVectorField &momentumToTransfer)
Correct kinematic transfers.
Type gSum(const FieldField< Field, Type > &f)
const cellShapeList & cells
virtual tmp< volVectorField > primaryMomentumTrans() const
Return momentum transfer source - Eulerian phase only.
virtual const volScalarField & cloudMassTrans() const
Return the film mass available for transfer to cloud.
virtual void preEvolveRegion()
Pre-evolve region.
Definition: regionModel.C:400
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:340
virtual void info(Ostream &os)
Provide some info.
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
A class for handling words, derived from string.
Definition: word.H:59
virtual const volScalarField & cloudDiameterTrans() const
Return the parcel diameters originating from film to cloud.
Calculate the face-flux of the given field.
static word groupName(Name name, const word &group)
#define DebugInFunction
Report an information message using Foam::Info.
Calculate the laplacian of the given field.
const IOdictionary & outputProperties() const
Return const access to the output properties dictionary.
Definition: regionModelI.H:104
scalar sumLocalContErr
const dimensionSet dimDensity
const dictionary & solution() const
Return the solution dictionary.
Definition: regionModelI.H:97
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calculate the matrix for the first temporal derivative.
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
bool loop()
Pimple loop.
Definition: pimpleControl.C:83
const Type & value() const
Return const reference to value.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
const volVectorField & nHat() const
Return the patch normal vectors.
static const word null
An empty word.
Definition: word.H:77
const surfaceScalarField & phi() const
Return the film flux [kg m/s].
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:55
Volume integrate volField creating a volField.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
word timeName
Definition: getTimeIndex.H:3
virtual void solveAlpha(const fvVectorMatrix &UEqn, const volScalarField &pc, const volScalarField &pe)
Solve for film volume fraction and thickness.
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
static const zero Zero
Definition: zero.H:97
volScalarField pSpPrimary_
Primary region normal momentum source (pressure) [kg m/s].
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
autoPtr< momentumTransportModel > momentumTransport_
Momentum transport model.
Calculate the divergence of the given field.
const dimensionedVector & g() const
Return the acceleration due to gravity.
dimensionedScalar totalMass
Definition: continuityErrs.H:4
bool correct()
Piso loop.
Definition: pisoControl.C:80
const labelList & primaryPatchIDs() const
Return the list of patch IDs on the primary region coupled.
Definition: regionModelI.H:166
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
virtual void info(Ostream &os)
Provide some info.
static const char nl
Definition: Ostream.H:260
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
const dimensionSet dimVelocity
scalar maxCo_
Optional maximum Courant number for stable film solution.
Type gMax(const FieldField< Field, Type > &f)
virtual void resetPrimaryRegionSourceTerms()
Reset source term fields.
const dimensionSet dimEnergy
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevTau(U)==fvModels.source(rho, U))
const Field< Type > & field() const
const dimensionSet dimMass
tmp< volScalarField > sigma() const
Return the surface tension coefficient [kg/s^2].
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
transferModelList transfer_
Transfer with the continuous phase.
kinematicSingleLayer(const word &modelType, const fvMesh &mesh, const dimensionedVector &g, const word &regionType, const bool readFields=true)
Construct from components.
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.
dictionary coeffs_
Model coefficients dictionary.
Definition: regionModel.H:94
virtual void updateSubmodels()
Update the film sub-models.
Calculate the matrix for the divergence of the given field and flux.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:809
label patchi
const volScalarField::Internal & magSf() const
Return the face area magnitudes [m^2].
virtual void predictDelta()
Predict delta_ from the continuity equation.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
tmp< volScalarField > rAU
const Time & time_
Reference to the time database.
Definition: regionModel.H:82
void max(const dimensioned< Type > &)
volVectorField & HbyA
Definition: pEqn.H:13
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
volScalarField alpha_
Film volume fraction in the cell layer [].
A List with indirect addressing.
Definition: fvMatrix.H:106
bool momentumPredictor() const
Flag to indicate to solve for momentum.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void correct(scalarField &availableMass, volScalarField &massToEject, volScalarField &diameterToEject)
Correct.
const dimensionSet dimVolume
void correctBoundaryConditions()
Correct boundary field.
tmp< surfaceScalarField > gGradRho() const
Hydrostatic pressure coefficient gradient.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
messageStream Info
volScalarField rhoSpPrimary_
Primary region mass source [kg].
void toRegion(const label regionPatchi, List< Type > &primaryFieldField) const
Convert a primary region field to the local region.
dimensioned< scalar > mag(const dimensioned< Type > &)
surfaceScalarField phiU_
Film velocity flux [m^3/s].
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
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.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:927
const volScalarField & rho() const
Return the film density [kg/m^3].
volVectorField::Internal Uw_
Velocity - wall [m/s].
volScalarField cloudMassTrans_
Film mass available for transfer to cloud.
static word modelName(Name name, const word &model)
Return the name of the object within the given model.
rDeltaTY field()
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
Definition: regionModel.H:106
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
Definition: regionModelI.H:173
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
const dictionary & coeffs() const
Return the model coefficients dictionary.
Definition: regionModelI.H:90
tmp< Type > constrainFilmField(const tmp< Type > &field, const typename Type::cmptType &value)
Constrain a film region master/slave boundaries of a field to a.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
defineTypeNameAndDebug(kinematicSingleLayer, 0)
virtual tmp< volScalarField::Internal > Sh() const
Return enthalpy source - Eulerian phase only.
autoPtr< viscosityModel > viscosity_
Viscosity model.
virtual tmp< volVectorField::Internal > SU() const
Return momentum source - Eulerian phase only.
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:46
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.
fvVectorMatrix & UEqn
Definition: UEqn.H:13