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 | Copyright (C) 2011-2015 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 #include "fvm.H"
28 #include "fvcDiv.H"
29 #include "fvcLaplacian.H"
30 #include "fvcSnGrad.H"
31 #include "fvcReconstruct.H"
32 #include "fvcVolumeIntegrate.H"
34 #include "mappedWallPolyPatch.H"
35 #include "mapDistribute.H"
36 #include "filmThermoModel.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace regionModels
43 {
44 namespace surfaceFilmModels
45 {
46 
47 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
48 
50 
52 
53 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
54 
56 {
58  {
59  const dictionary& solution = this->solution().subDict("PISO");
60  solution.lookup("momentumPredictor") >> momentumPredictor_;
61  solution.readIfPresent("nOuterCorr", nOuterCorr_);
62  solution.lookup("nCorr") >> nCorr_;
63  solution.lookup("nNonOrthCorr") >> nNonOrthCorr_;
64 
65  return true;
66  }
67  else
68  {
69  return false;
70  }
71 }
72 
73 
75 {
76  rho_ == filmThermo_->rho();
77  mu_ == filmThermo_->mu();
78  sigma_ == filmThermo_->sigma();
79 }
80 
81 
83 {
84  if (debug)
85  {
86  Info<< "kinematicSingleLayer::resetPrimaryRegionSourceTerms()" << endl;
87  }
88 
89  rhoSpPrimary_ == dimensionedScalar("zero", rhoSp_.dimensions(), 0.0);
90  USpPrimary_ == dimensionedVector("zero", USp_.dimensions(), vector::zero);
91  pSpPrimary_ == dimensionedScalar("zero", pSp_.dimensions(), 0.0);
92 }
93 
94 
96 {
97  if (debug)
98  {
99  Info<< "kinematicSingleLayer::"
100  << "transferPrimaryRegionThermoFields()" << endl;
101  }
102 
103  // Update fields from primary region via direct mapped
104  // (coupled) boundary conditions
109 }
110 
111 
113 {
114  if (debug)
115  {
116  Info<< "kinematicSingleLayer::"
117  << "transferPrimaryRegionSourceFields()" << endl;
118  }
119 
120  // Convert accummulated source terms into per unit area per unit time
121  const scalar deltaT = time_.deltaTValue();
123  {
124  const scalarField& priMagSf =
125  primaryMesh().magSf().boundaryField()[patchI];
126 
127  rhoSpPrimary_.boundaryField()[patchI] /= priMagSf*deltaT;
128  USpPrimary_.boundaryField()[patchI] /= priMagSf*deltaT;
129  pSpPrimary_.boundaryField()[patchI] /= priMagSf*deltaT;
130  }
131 
132  // Retrieve the source fields from the primary region via direct mapped
133  // (coupled) boundary conditions
134  // - fields require transfer of values for both patch AND to push the
135  // values into the first layer of internal cells
139 
140  // update addedMassTotal counter
141  if (time().outputTime())
142  {
143  if (debug)
144  {
145  rhoSp_.write();
146  USp_.write();
147  pSp_.write();
148  }
149 
150  scalar addedMassTotal = 0.0;
151  outputProperties().readIfPresent("addedMassTotal", addedMassTotal);
152  addedMassTotal += returnReduce(addedMassTotal_, sumOp<scalar>());
153  outputProperties().add("addedMassTotal", addedMassTotal, true);
154  addedMassTotal_ = 0.0;
155  }
156 }
157 
158 
160 {
161  return tmp<volScalarField>
162  (
163  new volScalarField
164  (
165  IOobject
166  (
167  typeName + ":pu",
168  time_.timeName(),
169  regionMesh(),
172  ),
173  pPrimary_ // pressure (mapped from primary region)
174  - pSp_ // accumulated particle impingement
175  - fvc::laplacian(sigma_, delta_) // surface tension
176  )
177  );
178 }
179 
180 
182 {
183  return tmp<volScalarField>
184  (
185  new volScalarField
186  (
187  IOobject
188  (
189  typeName + ":pp",
190  time_.timeName(),
191  regionMesh(),
194  ),
195  -rho_*gNormClipped() // hydrostatic effect only
196  )
197  );
198 }
199 
200 
202 {
204 }
205 
206 
208 {
209  if (debug)
210  {
211  Info<< "kinematicSingleLayer::updateSubmodels()" << endl;
212  }
213 
214  // Update injection model - mass returned is mass available for injection
216 
217  // Update source fields
218  const dimensionedScalar deltaT = time().deltaT();
219  rhoSp_ += cloudMassTrans_/magSf()/deltaT;
220 
221  turbulence_->correct();
222 }
223 
224 
226 {
227  const volScalarField deltaRho0(deltaRho_);
228 
229  solveContinuity();
230 
231  if (debug)
232  {
236  + dimensionedScalar("SMALL", dimMass*dimVolume, ROOTVSMALL);
237 
238  const scalar sumLocalContErr =
239  (
241  ).value();
242 
243  const scalar globalContErr =
244  (
246  ).value();
247 
249 
250  Info<< "Surface film: " << type() << nl
251  << " time step continuity errors: sum local = "
252  << sumLocalContErr << ", global = " << globalContErr
253  << ", cumulative = " << cumulativeContErr_ << endl;
254  }
255 }
256 
257 
259 {
260  if (debug)
261  {
262  Info<< "kinematicSingleLayer::solveContinuity()" << endl;
263  }
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 {
298  if (debug)
299  {
300  Info<< "kinematicSingleLayer::solveMomentum()" << endl;
301  }
302 
303  // Momentum
304  tmp<fvVectorMatrix> tUEqn
305  (
306  fvm::ddt(deltaRho_, U_)
307  + fvm::div(phi_, U_)
308  ==
309  - USp_
310 // - fvm::SuSp(rhoSp_, U_)
311  - rhoSp_*U_
312  + forces_.correct(U_)
313  + turbulence_->Su(U_)
314  );
315 
316  fvVectorMatrix& UEqn = tUEqn();
317 
318  UEqn.relax();
319 
320  if (momentumPredictor_)
321  {
322  solve
323  (
324  UEqn
325  ==
327  (
328  - fvc::interpolate(delta_)
329  * (
330  regionMesh().magSf()
331  * (
332  fvc::snGrad(pu, "snGrad(p)")
333  + fvc::snGrad(pp, "snGrad(p)")*fvc::interpolate(delta_)
334  + fvc::snGrad(delta_)*fvc::interpolate(pp)
335  )
336  - (fvc::interpolate(rho_*gTan()) & regionMesh().Sf())
337  )
338  )
339  );
340 
341  // Remove any patch-normal components of velocity
342  U_ -= nHat()*(nHat() & U_);
343  U_.correctBoundaryConditions();
344  }
345 
346  return tUEqn;
347 }
348 
349 
351 (
352  const volScalarField& pu,
353  const volScalarField& pp,
354  const fvVectorMatrix& UEqn
355 )
356 {
357  if (debug)
358  {
359  Info<< "kinematicSingleLayer::solveThickness()" << endl;
360  }
361 
362  volScalarField rUA(1.0/UEqn.A());
363  U_ = rUA*UEqn.H();
364 
365  surfaceScalarField deltarUAf(fvc::interpolate(delta_*rUA));
367 
368  surfaceScalarField phiAdd
369  (
370  "phiAdd",
371  regionMesh().magSf()
372  * (
373  fvc::snGrad(pu, "snGrad(p)")
374  + fvc::snGrad(pp, "snGrad(p)")*fvc::interpolate(delta_)
375  )
376  - (fvc::interpolate(rho_*gTan()) & regionMesh().Sf())
377  );
378  constrainFilmField(phiAdd, 0.0);
379 
381  (
382  "phid",
383  (fvc::interpolate(U_*rho_) & regionMesh().Sf())
384  - deltarUAf*phiAdd*rhof
385  );
386  constrainFilmField(phid, 0.0);
387 
388  surfaceScalarField ddrhorUAppf
389  (
390  "deltaCoeff",
391  fvc::interpolate(delta_)*deltarUAf*rhof*fvc::interpolate(pp)
392  );
393 
394  regionMesh().setFluxRequired(delta_.name());
395 
396  for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
397  {
398  // Film thickness equation
399  fvScalarMatrix deltaEqn
400  (
401  fvm::ddt(rho_, delta_)
402  + fvm::div(phid, delta_)
403  - fvm::laplacian(ddrhorUAppf, delta_)
404  ==
405  - rhoSp_
406  );
407 
408  deltaEqn.solve();
409 
410  if (nonOrth == nNonOrthCorr_)
411  {
412  phiAdd +=
413  fvc::interpolate(pp)
414  * fvc::snGrad(delta_)
415  * regionMesh().magSf();
416 
417  phi_ == deltaEqn.flux();
418  }
419  }
420 
421  // Bound film thickness by a minimum of zero
422  delta_.max(0.0);
423 
424  // Update U field
425  U_ -= fvc::reconstruct(deltarUAf*phiAdd);
426 
427  // Remove any patch-normal components of velocity
428  U_ -= nHat()*(nHat() & U_);
429 
430  U_.correctBoundaryConditions();
431 
432  // Update film wall and surface velocities
433  updateSurfaceVelocities();
434 
435  // Continuity check
436  continuityCheck();
437 }
438 
439 
440 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
441 
443 (
444  const word& modelType,
445  const fvMesh& mesh,
446  const dimensionedVector& g,
447  const word& regionType,
448  const bool readFields
449 )
450 :
451  surfaceFilmModel(modelType, mesh, g, regionType),
452 
453  momentumPredictor_(solution().subDict("PISO").lookup("momentumPredictor")),
454  nOuterCorr_(solution().subDict("PISO").lookupOrDefault("nOuterCorr", 1)),
455  nCorr_(readLabel(solution().subDict("PISO").lookup("nCorr"))),
456  nNonOrthCorr_
457  (
458  readLabel(solution().subDict("PISO").lookup("nNonOrthCorr"))
459  ),
460 
461  cumulativeContErr_(0.0),
462 
463  deltaSmall_("deltaSmall", dimLength, SMALL),
464  deltaCoLimit_(solution().lookupOrDefault("deltaCoLimit", 1e-4)),
465 
466  rho_
467  (
468  IOobject
469  (
470  "rhof",
471  time().timeName(),
472  regionMesh(),
475  ),
476  regionMesh(),
477  dimensionedScalar("zero", dimDensity, 0.0),
478  zeroGradientFvPatchScalarField::typeName
479  ),
480  mu_
481  (
482  IOobject
483  (
484  "muf",
485  time().timeName(),
486  regionMesh(),
489  ),
490  regionMesh(),
491  dimensionedScalar("zero", dimPressure*dimTime, 0.0),
492  zeroGradientFvPatchScalarField::typeName
493  ),
494  sigma_
495  (
496  IOobject
497  (
498  "sigmaf",
499  time().timeName(),
500  regionMesh(),
503  ),
504  regionMesh(),
505  dimensionedScalar("zero", dimMass/sqr(dimTime), 0.0),
506  zeroGradientFvPatchScalarField::typeName
507  ),
508 
509  delta_
510  (
511  IOobject
512  (
513  "deltaf",
514  time().timeName(),
515  regionMesh(),
518  ),
519  regionMesh()
520  ),
521  alpha_
522  (
523  IOobject
524  (
525  "alpha",
526  time().timeName(),
527  regionMesh(),
530  ),
531  regionMesh(),
532  dimensionedScalar("zero", dimless, 0.0),
533  zeroGradientFvPatchScalarField::typeName
534  ),
535  U_
536  (
537  IOobject
538  (
539  "Uf",
540  time().timeName(),
541  regionMesh(),
544  ),
545  regionMesh()
546  ),
547  Us_
548  (
549  IOobject
550  (
551  "Usf",
552  time().timeName(),
553  regionMesh(),
556  ),
557  U_,
558  zeroGradientFvPatchScalarField::typeName
559  ),
560  Uw_
561  (
562  IOobject
563  (
564  "Uwf",
565  time().timeName(),
566  regionMesh(),
569  ),
570  U_,
571  zeroGradientFvPatchScalarField::typeName
572  ),
573  deltaRho_
574  (
575  IOobject
576  (
577  delta_.name() + "*" + rho_.name(),
578  time().timeName(),
579  regionMesh(),
582  ),
583  regionMesh(),
584  dimensionedScalar("zero", delta_.dimensions()*rho_.dimensions(), 0.0),
585  zeroGradientFvPatchScalarField::typeName
586  ),
587 
588  phi_
589  (
590  IOobject
591  (
592  "phi",
593  time().timeName(),
594  regionMesh(),
597  ),
598  regionMesh(),
600  ),
601 
602  primaryMassTrans_
603  (
604  IOobject
605  (
606  "primaryMassTrans",
607  time().timeName(),
608  regionMesh(),
611  ),
612  regionMesh(),
613  dimensionedScalar("zero", dimMass, 0.0),
614  zeroGradientFvPatchScalarField::typeName
615  ),
616  cloudMassTrans_
617  (
618  IOobject
619  (
620  "cloudMassTrans",
621  time().timeName(),
622  regionMesh(),
625  ),
626  regionMesh(),
627  dimensionedScalar("zero", dimMass, 0.0),
628  zeroGradientFvPatchScalarField::typeName
629  ),
630  cloudDiameterTrans_
631  (
632  IOobject
633  (
634  "cloudDiameterTrans",
635  time().timeName(),
636  regionMesh(),
639  ),
640  regionMesh(),
641  dimensionedScalar("zero", dimLength, -1.0),
642  zeroGradientFvPatchScalarField::typeName
643  ),
644 
645  USp_
646  (
647  IOobject
648  (
649  "USpf",
650  time().timeName(),
651  regionMesh(),
654  ),
655  regionMesh(),
657  (
659  ),
660  this->mappedPushedFieldPatchTypes<vector>()
661  ),
662  pSp_
663  (
664  IOobject
665  (
666  "pSpf",
667  time_.timeName(),
668  regionMesh(),
671  ),
672  regionMesh(),
673  dimensionedScalar("zero", dimPressure, 0.0),
674  this->mappedPushedFieldPatchTypes<scalar>()
675  ),
676  rhoSp_
677  (
678  IOobject
679  (
680  "rhoSpf",
681  time_.timeName(),
682  regionMesh(),
685  ),
686  regionMesh(),
688  this->mappedPushedFieldPatchTypes<scalar>()
689  ),
690 
691  USpPrimary_
692  (
693  IOobject
694  (
695  USp_.name(), // must have same name as USp_ to enable mapping
696  time().timeName(),
697  primaryMesh(),
700  ),
701  primaryMesh(),
702  dimensionedVector("zero", USp_.dimensions(), vector::zero)
703  ),
704  pSpPrimary_
705  (
706  IOobject
707  (
708  pSp_.name(), // must have same name as pSp_ to enable mapping
709  time().timeName(),
710  primaryMesh(),
713  ),
714  primaryMesh(),
715  dimensionedScalar("zero", pSp_.dimensions(), 0.0)
716  ),
717  rhoSpPrimary_
718  (
719  IOobject
720  (
721  rhoSp_.name(), // must have same name as rhoSp_ to enable mapping
722  time().timeName(),
723  primaryMesh(),
726  ),
727  primaryMesh(),
728  dimensionedScalar("zero", rhoSp_.dimensions(), 0.0)
729  ),
730 
731  UPrimary_
732  (
733  IOobject
734  (
735  "U", // must have same name as U to enable mapping
736  time().timeName(),
737  regionMesh(),
740  ),
741  regionMesh(),
743  this->mappedFieldAndInternalPatchTypes<vector>()
744  ),
745  pPrimary_
746  (
747  IOobject
748  (
749  "p", // must have same name as p to enable mapping
750  time().timeName(),
751  regionMesh(),
754  ),
755  regionMesh(),
756  dimensionedScalar("zero", dimPressure, 0.0),
757  this->mappedFieldAndInternalPatchTypes<scalar>()
758  ),
759  rhoPrimary_
760  (
761  IOobject
762  (
763  "rho", // must have same name as rho to enable mapping
764  time().timeName(),
765  regionMesh(),
768  ),
769  regionMesh(),
770  dimensionedScalar("zero", dimDensity, 0.0),
771  this->mappedFieldAndInternalPatchTypes<scalar>()
772  ),
773  muPrimary_
774  (
775  IOobject
776  (
777  "thermo:mu", // must have same name as mu to enable mapping
778  time().timeName(),
779  regionMesh(),
782  ),
783  regionMesh(),
784  dimensionedScalar("zero", dimPressure*dimTime, 0.0),
785  this->mappedFieldAndInternalPatchTypes<scalar>()
786  ),
787 
788  filmThermo_(filmThermoModel::New(*this, coeffs_)),
789 
790  availableMass_(regionMesh().nCells(), 0.0),
791 
792  injection_(*this, coeffs_),
793 
794  turbulence_(filmTurbulenceModel::New(*this, coeffs_)),
795 
796  forces_(*this, coeffs_),
797 
798  addedMassTotal_(0.0)
799 {
800  if (readFields)
801  {
802  transferPrimaryRegionThermoFields();
803 
804  correctAlpha();
805 
806  correctThermoFields();
807 
808  deltaRho_ == delta_*rho_;
809 
811  (
812  IOobject
813  (
814  "phi",
815  time().timeName(),
816  regionMesh(),
819  false
820  ),
821  fvc::interpolate(deltaRho_*U_) & regionMesh().Sf()
822  );
823 
824  phi_ == phi0;
825  }
826 }
827 
828 
829 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
830 
832 {}
833 
834 
835 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
836 
838 (
839  const label patchI,
840  const label faceI,
841  const scalar massSource,
842  const vector& momentumSource,
843  const scalar pressureSource,
844  const scalar energySource
845 )
846 {
847  if (debug)
848  {
849  Info<< "\nSurface film: " << type() << ": adding to film source:" << nl
850  << " mass = " << massSource << nl
851  << " momentum = " << momentumSource << nl
852  << " pressure = " << pressureSource << endl;
853  }
854 
855  rhoSpPrimary_.boundaryField()[patchI][faceI] -= massSource;
856  USpPrimary_.boundaryField()[patchI][faceI] -= momentumSource;
857  pSpPrimary_.boundaryField()[patchI][faceI] -= pressureSource;
858 
859  addedMassTotal_ += massSource;
860 }
861 
862 
864 {
865  if (debug)
866  {
867  Info<< "kinematicSingleLayer::preEvolveRegion()" << endl;
868  }
869 
871 
873 
875 
877 
879 
880  correctAlpha();
881 
882  // Reset transfer fields
883 // availableMass_ = mass();
885  cloudMassTrans_ == dimensionedScalar("zero", dimMass, 0.0);
887 }
888 
889 
891 {
892  if (debug)
893  {
894  Info<< "kinematicSingleLayer::evolveRegion()" << endl;
895  }
896 
897  // Update sub-models to provide updated source contributions
898  updateSubmodels();
899 
900  // Solve continuity for deltaRho_
901  solveContinuity();
902 
903  // Implicit pressure source coefficient - constant
904  tmp<volScalarField> tpp(this->pp());
905 
906  for (int oCorr=1; oCorr<=nOuterCorr_; oCorr++)
907  {
908  // Explicit pressure source contribution - varies with delta_
909  tmp<volScalarField> tpu(this->pu());
910 
911  // Solve for momentum for U_
912  tmp<fvVectorMatrix> UEqn = solveMomentum(tpu(), tpp());
913 
914  // Film thickness correction loop
915  for (int corr=1; corr<=nCorr_; corr++)
916  {
917  // Solve thickness for delta_
918  solveThickness(tpu(), tpp(), UEqn());
919  }
920  }
921 
922  // Update deltaRho_ with new delta_
923  deltaRho_ == delta_*rho_;
924 }
925 
926 
928 {
929  if (debug)
930  {
931  Info<< "kinematicSingleLayer::postEvolveRegion()" << endl;
932  }
933 
934  // Reset source terms for next time integration
936 }
937 
938 
940 {
941  scalar CoNum = 0.0;
942 
943  if (regionMesh().nInternalFaces() > 0)
944  {
945  const scalarField sumPhi
946  (
948  / (deltaRho_.internalField() + ROOTVSMALL)
949  );
950 
951  forAll(delta_, i)
952  {
953  if ((delta_[i] > deltaCoLimit_) && (alpha_[i] > 0.5))
954  {
955  CoNum = max(CoNum, sumPhi[i]/(delta_[i]*magSf()[i]));
956  }
957  }
958 
959  CoNum *= 0.5*time_.deltaTValue();
960  }
961 
963 
964  Info<< "Film max Courant number: " << CoNum << endl;
965 
966  return CoNum;
967 }
968 
969 
971 {
972  return U_;
973 }
974 
975 
977 {
978  return Us_;
979 }
980 
981 
983 {
984  return Uw_;
985 }
986 
987 
989 {
990  return deltaRho_;
991 }
992 
993 
995 {
996  return phi_;
997 }
998 
999 
1001 {
1002  return rho_;
1003 }
1004 
1005 
1007 {
1009  << "T field not available for " << type() << abort(FatalError);
1010 
1011  return volScalarField::null();
1012 }
1013 
1014 
1016 {
1018  << "Ts field not available for " << type() << abort(FatalError);
1019 
1020  return volScalarField::null();
1021 }
1022 
1023 
1025 {
1027  << "Tw field not available for " << type() << abort(FatalError);
1028 
1029  return volScalarField::null();
1030 }
1031 
1032 
1034 {
1036  << "Cp field not available for " << type() << abort(FatalError);
1037 
1038  return volScalarField::null();
1039 }
1040 
1041 
1043 {
1045  << "kappa field not available for " << type() << abort(FatalError);
1046 
1047  return volScalarField::null();
1048 }
1049 
1050 
1052 {
1053  return tmp<volScalarField>
1054  (
1055  new volScalarField
1056  (
1057  IOobject
1058  (
1059  typeName + ":primaryMassTrans",
1060  time().timeName(),
1061  primaryMesh(),
1064  false
1065  ),
1066  primaryMesh(),
1068  )
1069  );
1070 }
1071 
1072 
1074 {
1075  return cloudMassTrans_;
1076 }
1077 
1078 
1080 {
1081  return cloudDiameterTrans_;
1082 }
1083 
1084 
1086 {
1087  Info<< "\nSurface film: " << type() << endl;
1088 
1089  const scalarField& deltaInternal = delta_.internalField();
1090  const vectorField& Uinternal = U_.internalField();
1091  scalar addedMassTotal = 0.0;
1092  outputProperties().readIfPresent("addedMassTotal", addedMassTotal);
1093  addedMassTotal += returnReduce(addedMassTotal_, sumOp<scalar>());
1094 
1095  Info<< indent << "added mass = " << addedMassTotal << nl
1096  << indent << "current mass = "
1097  << gSum((deltaRho_*magSf())()) << nl
1098  << indent << "min/max(mag(U)) = " << gMin(mag(Uinternal)) << ", "
1099  << gMax(mag(Uinternal)) << nl
1100  << indent << "min/max(delta) = " << gMin(deltaInternal) << ", "
1101  << gMax(deltaInternal) << nl
1102  << indent << "coverage = "
1103  << gSum(alpha_.internalField()*magSf())/gSum(magSf()) << nl;
1104 
1105  injection_.info(Info);
1106 }
1107 
1108 
1110 {
1112  (
1114  (
1115  IOobject
1116  (
1117  typeName + ":Srho",
1118  time().timeName(),
1119  primaryMesh(),
1122  false
1123  ),
1124  primaryMesh(),
1126  )
1127  );
1128 }
1129 
1130 
1133  const label i
1134 ) const
1135 {
1137  (
1139  (
1140  IOobject
1141  (
1142  typeName + ":Srho(" + Foam::name(i) + ")",
1143  time().timeName(),
1144  primaryMesh(),
1147  false
1148  ),
1149  primaryMesh(),
1151  )
1152  );
1153 }
1154 
1155 
1157 {
1159  (
1161  (
1162  IOobject
1163  (
1164  typeName + ":Sh",
1165  time().timeName(),
1166  primaryMesh(),
1169  false
1170  ),
1171  primaryMesh(),
1173  )
1174  );
1175 }
1176 
1177 
1178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1179 
1180 } // End namespace surfaceFilmModels
1181 } // End namespace regionModels
1182 } // End namespace Foam
1183 
1184 // ************************************************************************* //
totalMass
dimensionedScalar totalMass
Definition: continuityErrs.H:4
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::T
virtual const volScalarField & T() const
Return the film mean temperature [K].
Definition: kinematicSingleLayer.C:1006
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::read
virtual bool read()
Read control parameters from dictionary.
Definition: kinematicSingleLayer.C:55
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::correctThermoFields
virtual void correctThermoFields()
Correct the thermo fields.
Definition: kinematicSingleLayer.C:74
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rhoSpPrimary_
volScalarField rhoSpPrimary_
Mass / [kg/m2/s].
Definition: kinematicSingleLayer.H:182
Foam::fvc::reconstruct
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcReconstruct.C:54
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::deltaCoLimit_
scalar deltaCoLimit_
Film thickness above which Courant number calculation in valid.
Definition: kinematicSingleLayer.H:103
Foam::dimPressure
const dimensionSet dimPressure
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Foam::fvc::surfaceSum
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcSurfaceIntegrate.C:138
sumLocalContErr
scalar sumLocalContErr
Definition: compressibleContinuityErrs.H:32
Foam::maxOp
Definition: ops.H:172
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::~kinematicSingleLayer
virtual ~kinematicSingleLayer()
Destructor.
Definition: kinematicSingleLayer.C:831
filmThermoModel.H
mappedWallPolyPatch.H
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::transferPrimaryRegionThermoFields
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
Definition: kinematicSingleLayer.C:95
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Us_
volVectorField Us_
Velocity - surface / [m/s].
Definition: kinematicSingleLayer.H:132
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Srho
virtual tmp< DimensionedField< scalar, volMesh > > Srho() const
Return total mass source - Eulerian phase only.
Definition: kinematicSingleLayer.C:1109
Foam::solution
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::USpPrimary_
volVectorField USpPrimary_
Momementum / [kg/m/s2].
Definition: kinematicSingleLayer.H:176
phid
surfaceScalarField phid("phid", fvc::interpolate(psi) *((mesh.Sf() &fvc::interpolate(HbyA))+rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho)))
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::USp_
volVectorField USp_
Momementum / [kg/m/s2].
Definition: kinematicSingleLayer.H:163
Foam::TimeState::deltaT
dimensionedScalar deltaT() const
Return time step.
Definition: TimeState.C:79
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Us
virtual const volVectorField & Us() const
Return the film surface velocity [m/s].
Definition: kinematicSingleLayer.C:976
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::injection_
injectionModelList injection_
Cloud injection.
Definition: kinematicSingleLayer.H:210
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve film hook.
Definition: kinematicSingleLayer.C:863
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::regionModels::regionModel::outputProperties
const IOdictionary & outputProperties() const
Return const access to the output properties dictionary.
Definition: regionModelI.H:110
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::dimDensity
const dimensionSet dimDensity
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Ts
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
Definition: kinematicSingleLayer.C:1015
Foam::dictionary::readIfPresent
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Definition: dictionaryTemplates.C:94
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::turbulence_
autoPtr< filmTurbulenceModel > turbulence_
Turbulence model.
Definition: kinematicSingleLayer.H:213
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::filmThermo_
autoPtr< filmThermoModel > filmThermo_
Film thermo model.
Definition: kinematicSingleLayer.H:204
Foam::fvc::domainIntegrate
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcVolumeIntegrate.C:86
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::U
virtual const volVectorField & U() const
Return the film velocity [m/s].
Definition: kinematicSingleLayer.C:970
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::deltaSmall_
const dimensionedScalar deltaSmall_
Small delta.
Definition: kinematicSingleLayer.H:100
fvcSnGrad.H
Calculate the snGrad of the given volField.
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::evolveRegion
virtual void evolveRegion()
Evolve the film equations.
Definition: kinematicSingleLayer.C:890
Foam::fvc::interpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Definition: surfaceInterpolate.C:110
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::solveContinuity
virtual void solveContinuity()
Solve continuity equation.
Definition: kinematicSingleLayer.C:258
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::regionModels::surfaceFilmModels::injectionModelList::info
virtual void info(Ostream &os)
Provide some info.
Definition: injectionModelList.C:128
fvcDiv.H
Calculate the divergence of the given field.
Foam::regionModels::regionModel::time
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::solveThickness
virtual void solveThickness(const volScalarField &pu, const volScalarField &pp, const fvVectorMatrix &UEqn)
Solve coupled velocity-thickness equations.
Definition: kinematicSingleLayer.C:351
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rho_
volScalarField rho_
Density / [kg/m3].
Definition: kinematicSingleLayer.H:111
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::sigma_
volScalarField sigma_
Surface tension / [m/s2].
Definition: kinematicSingleLayer.H:117
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
g
const dimensionedVector & g
Definition: setRegionFluidFields.H:33
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::continuityCheck
virtual void continuityCheck()
Continuity check.
Definition: kinematicSingleLayer.C:225
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::readFields
This function object reads fields from the time directories and adds them to the mesh database for fu...
Definition: readFields.H:104
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::availableMass_
scalarField availableMass_
Available mass for transfer via sub-models.
Definition: kinematicSingleLayer.H:207
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pPrimary_
volScalarField pPrimary_
Pressure / [Pa].
Definition: kinematicSingleLayer.H:192
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::muPrimary_
volScalarField muPrimary_
Viscosity / [Pa.s].
Definition: kinematicSingleLayer.H:198
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:48
Foam::fvMesh::magSf
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Definition: fvMeshGeometry.C:358
Foam::regionModels::regionModel::primaryMesh
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionModelI.H:31
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::gNormClipped
tmp< volScalarField > gNormClipped() const
Return the gravity normal-to-patch component contribution.
Definition: kinematicSingleLayerI.H:235
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::cloudDiameterTrans_
volScalarField cloudDiameterTrans_
Parcel diameters originating from film to cloud.
Definition: kinematicSingleLayer.H:153
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::deltaRho_
volScalarField deltaRho_
Film thickness*density (helper field) / [kg/m2].
Definition: kinematicSingleLayer.H:138
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::cloudMassTrans
virtual const volScalarField & cloudMassTrans() const
Return the film mass available for transfer to cloud.
Definition: kinematicSingleLayer.C:1073
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::addedMassTotal_
scalar addedMassTotal_
Cumulative mass added via sources [kg].
Definition: kinematicSingleLayer.H:222
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::transferPrimaryRegionSourceFields
virtual void transferPrimaryRegionSourceFields()
Transfer source fields from the primary region to the film region.
Definition: kinematicSingleLayer.C:112
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::delta_
volScalarField delta_
Film thickness / [m].
Definition: kinematicSingleLayer.H:123
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::kappa
virtual const volScalarField & kappa() const
Return the film thermal conductivity [W/m/K].
Definition: kinematicSingleLayer.C:1042
Foam::regionModels::singleLayerRegion::magSf
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
Definition: singleLayerRegion.C:220
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pp
virtual tmp< volScalarField > pp()
Implicit pressure source coefficient.
Definition: kinematicSingleLayer.C:181
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::cloudDiameterTrans
virtual const volScalarField & cloudDiameterTrans() const
Return the parcel diameters originating from film to cloud.
Definition: kinematicSingleLayer.C:1079
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::info
virtual void info()
Provide some feedback.
Definition: kinematicSingleLayer.C:1085
Foam::regionModels::regionModel::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve region.
Definition: regionModel.C:543
UEqn
tmp< fvVectorMatrix > UEqn(fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::constant::electromagnetic::phi0
const dimensionedScalar phi0
Magnetic flux quantum: default SI units: [Wb].
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::TimeState::deltaTValue
scalar deltaTValue() const
Return time step value.
Definition: TimeState.H:100
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Uw_
volVectorField Uw_
Velocity - wall / [m/s].
Definition: kinematicSingleLayer.H:135
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::primaryMassTrans
virtual tmp< volScalarField > primaryMassTrans() const
Return mass transfer source - Eulerian phase only.
Definition: kinematicSingleLayer.C:1051
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::resetPrimaryRegionSourceTerms
virtual void resetPrimaryRegionSourceTerms()
Reset source term fields.
Definition: kinematicSingleLayer.C:82
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::solveMomentum
virtual tmp< fvVectorMatrix > solveMomentum(const volScalarField &pu, const volScalarField &pp)
Solve for film velocity.
Definition: kinematicSingleLayer.C:293
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::phi
virtual const surfaceScalarField & phi() const
Return the film flux [kg.m/s].
Definition: kinematicSingleLayer.C:994
CoNum
scalar CoNum
Definition: compressibleCourantNo.H:29
Foam::regionModels::regionModel::solution
const dictionary & solution() const
Return the solution dictionary.
Definition: regionModelI.H:103
Foam::regionModels::surfaceFilmModels::filmTurbulenceModel::New
static autoPtr< filmTurbulenceModel > New(surfaceFilmModel &owner, const dictionary &dict)
Return a reference to the selected injection model.
Definition: filmTurbulenceModelNew.C:40
Foam::regionModels::surfaceFilmModels::surfaceFilmModel::read
virtual bool read()
Read control parameters from dictionary.
Definition: surfaceFilmModel.C:45
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::correctAlpha
virtual void correctAlpha()
Correct film coverage field.
Definition: kinematicSingleLayer.C:201
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::momentumPredictor_
Switch momentumPredictor_
Momentum predictor.
Definition: kinematicSingleLayer.H:85
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::addSources
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.
Definition: kinematicSingleLayer.C:838
fvm.H
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
kinematicSingleLayer.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::U_
volVectorField U_
Velocity - mean / [m/s].
Definition: kinematicSingleLayer.H:129
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::nCorr_
label nCorr_
Number of PISO-like correctors.
Definition: kinematicSingleLayer.H:91
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::mass
tmp< volScalarField > mass() const
Return the current film mass.
Definition: kinematicSingleLayerI.H:191
globalContErr
scalar globalContErr
Definition: compressibleContinuityErrs.H:35
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pu
virtual tmp< volScalarField > pu()
Explicit pressure source contribution.
Definition: kinematicSingleLayer.C:159
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::UPrimary_
volVectorField UPrimary_
Velocity / [m/s].
Definition: kinematicSingleLayer.H:189
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pSp_
volScalarField pSp_
Pressure / [Pa].
Definition: kinematicSingleLayer.H:166
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::cumulativeContErr_
scalar cumulativeContErr_
Cumulative continuity error.
Definition: kinematicSingleLayer.H:97
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::phi_
surfaceScalarField phi_
Mass flux (includes film thickness) / [kg.m/s].
Definition: kinematicSingleLayer.H:141
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
internalField
conserve internalField()+
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::regionModels::surfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicSingleLayer, 0)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::solve
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::regionModels::regionModel::time_
const Time & time_
Reference to the time database.
Definition: regionModel.H:90
fvcLaplacian.H
Calculate the laplacian of the given field.
Foam::sumOp
Definition: ops.H:162
mapDistribute.H
Foam::regionModels::surfaceFilmModels::injectionModelList::correct
virtual void correct(scalarField &availableMass, volScalarField &massToInject, volScalarField &diameterToInject)
Correct.
Definition: injectionModelList.C:100
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::alpha_
volScalarField alpha_
Film coverage indicator, 1 = covered, 0 = uncovered / [].
Definition: kinematicSingleLayer.H:126
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pSpPrimary_
volScalarField pSpPrimary_
Pressure / [Pa].
Definition: kinematicSingleLayer.H:179
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::mu_
volScalarField mu_
Dynamic viscosity / [Pa.s].
Definition: kinematicSingleLayer.H:114
Foam::Vector< scalar >
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::nNonOrthCorr_
label nNonOrthCorr_
Number of non-orthogonal correctors.
Definition: kinematicSingleLayer.H:94
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Uw
virtual const volVectorField & Uw() const
Return the film wall velocity [m/s].
Definition: kinematicSingleLayer.C:982
kinematicSingleLayer
Kinematic form of single-cell layer surface film model.
Foam::regionModels::singleLayerRegion::nHat
virtual const volVectorField & nHat() const
Return the patch normal vectors.
Definition: singleLayerRegion.C:207
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rhoPrimary_
volScalarField rhoPrimary_
Density / [kg/m3].
Definition: kinematicSingleLayer.H:195
fvcVolumeIntegrate.H
Volume integrate volField creating a volField.
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::updateSurfaceVelocities
virtual void updateSurfaceVelocities()
Update film surface velocities.
Definition: kinematicSingleLayer.C:275
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::regionModels::surfaceFilmModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(surfaceFilmModel, kinematicSingleLayer, mesh)
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::netMass
tmp< volScalarField > netMass() const
Return the net film mass available over the next integration.
Definition: kinematicSingleLayerI.H:197
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::readLabel
label readLabel(Istream &is)
Definition: label.H:64
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Cp
virtual const volScalarField & Cp() const
Return the film specific heat capacity [J/kg/K].
Definition: kinematicSingleLayer.C:1033
Foam::regionModels::regionModel::intCoupledPatchIDs_
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
Definition: regionModel.H:117
fvcReconstruct.H
Reconstruct volField from a face flux field.
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::kinematicSingleLayer
kinematicSingleLayer(const kinematicSingleLayer &)
Disallow default bitwise copy construct.
Foam::regionModels::surfaceFilmModels::surfaceFilmModel
Base class for surface film models.
Definition: surfaceFilmModel.H:61
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::GeometricField::null
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
Definition: GeometricFieldI.H:30
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Tw
virtual const volScalarField & Tw() const
Return the film wall temperature [K].
Definition: kinematicSingleLayer.C:1024
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::postEvolveRegion
virtual void postEvolveRegion()
Post-evolve film hook.
Definition: kinematicSingleLayer.C:927
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rho
virtual const volScalarField & rho() const
Return the film density [kg/m3].
Definition: kinematicSingleLayer.C:1000
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::deltaRho
virtual const volScalarField & deltaRho() const
Return the film thickness*density (helper field) [kg/m3].
Definition: kinematicSingleLayer.C:988
rhof
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::updateSubmodels
virtual void updateSubmodels()
Update the film sub-models.
Definition: kinematicSingleLayer.C:207
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Foam::IOobject::READ_IF_PRESENT
@ READ_IF_PRESENT
Definition: IOobject.H:110
Foam::fvMatrix::flux
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:873
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rhoSp_
volScalarField rhoSp_
Mass / [kg/m2/s].
Definition: kinematicSingleLayer.H:169
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::cloudMassTrans_
volScalarField cloudMassTrans_
Film mass available for transfer to cloud.
Definition: kinematicSingleLayer.H:150
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Sh
virtual tmp< DimensionedField< scalar, volMesh > > Sh() const
Return enthalpy source - Eulerian phase only.
Definition: kinematicSingleLayer.C:1156
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::nOuterCorr_
label nOuterCorr_
Number of outer correctors.
Definition: kinematicSingleLayer.H:88
Foam::dictionary::add
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:729
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::CourantNumber
virtual scalar CourantNumber() const
Courant number evaluation.
Definition: kinematicSingleLayer.C:939
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:61
Foam::regionModels::surfaceFilmModels::filmThermoModel::New
static autoPtr< filmThermoModel > New(surfaceFilmModel &owner, const dictionary &dict)
Return a reference to the selected phase change model.
Definition: filmThermoModelNew.C:40
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190