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-2018 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"
33 #include "fvcFlux.H"
35 #include "mappedWallPolyPatch.H"
36 #include "mapDistribute.H"
37 #include "filmThermoModel.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 namespace regionModels
44 {
45 namespace surfaceFilmModels
46 {
47 
48 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
49 
50 defineTypeNameAndDebug(kinematicSingleLayer, 0);
51 
52 addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh);
53 
54 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
55 
57 {
59  {
60  const dictionary& solution = this->solution().subDict("PISO");
61  solution.lookup("momentumPredictor") >> momentumPredictor_;
62  solution.readIfPresent("nOuterCorr", nOuterCorr_);
63  solution.lookup("nCorr") >> nCorr_;
64  solution.lookup("nNonOrthCorr") >> nNonOrthCorr_;
65 
66  return true;
67  }
68  else
69  {
70  return false;
71  }
72 }
73 
74 
76 {
77  rho_ == filmThermo_->rho();
78  mu_ == filmThermo_->mu();
79  sigma_ == filmThermo_->sigma();
80 }
81 
82 
84 {
85  if (debug)
86  {
88  }
89 
93 }
94 
95 
97 {
98  if (debug)
99  {
100  InfoInFunction << 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  InfoInFunction << endl;
117  }
118 
119  volScalarField::Boundary& rhoSpPrimaryBf =
121 
122  volVectorField::Boundary& USpPrimaryBf =
124 
125  volScalarField::Boundary& pSpPrimaryBf =
127 
128  // Convert accummulated source terms into per unit area per unit time
129  const scalar deltaT = time_.deltaTValue();
131  {
132  scalarField rpriMagSfdeltaT
133  (
134  (1.0/deltaT)
135  /primaryMesh().magSf().boundaryField()[patchi]
136  );
137 
138  rhoSpPrimaryBf[patchi] *= rpriMagSfdeltaT;
139  USpPrimaryBf[patchi] *= rpriMagSfdeltaT;
140  pSpPrimaryBf[patchi] *= rpriMagSfdeltaT;
141  }
142 
143  // Retrieve the source fields from the primary region via direct mapped
144  // (coupled) boundary conditions
145  // - fields require transfer of values for both patch AND to push the
146  // values into the first layer of internal cells
150 
151  // update addedMassTotal counter
152  if (time().writeTime())
153  {
154  scalar addedMassTotal = 0.0;
155  outputProperties().readIfPresent("addedMassTotal", addedMassTotal);
156  addedMassTotal += returnReduce(addedMassTotal_, sumOp<scalar>());
157  outputProperties().add("addedMassTotal", addedMassTotal, true);
158  addedMassTotal_ = 0.0;
159  }
160 }
161 
162 
164 {
165  return volScalarField::New
166  (
167  typeName + ":pu",
168  pPrimary_ // pressure (mapped from primary region)
169  - pSp_ // accumulated particle impingement
170  - fvc::laplacian(sigma_, delta_) // surface tension
171  );
172 }
173 
174 
176 {
177  return volScalarField::New
178  (
179  typeName + ":pp",
180  -rho_*gNormClipped() // hydrostatic effect only
181  );
182 }
183 
184 
186 {
188 }
189 
190 
192 {
193  if (debug)
194  {
195  InfoInFunction << endl;
196  }
197 
198  // Update injection model - mass returned is mass available for injection
200 
201  // Update transfer model - mass returned is mass available for transfer
203 
204  // Update mass source field
206 
207  turbulence_->correct();
208 }
209 
210 
212 {
213  const volScalarField deltaRho0(deltaRho_);
214 
215  solveContinuity();
216 
217  if (debug)
218  {
222  + dimensionedScalar(dimMass*dimVolume, rootVSmall);
223 
224  const scalar sumLocalContErr =
225  (
226  fvc::domainIntegrate(mag(mass - magSf()*deltaRho0))/totalMass
227  ).value();
228 
229  const scalar globalContErr =
230  (
231  fvc::domainIntegrate(mass - magSf()*deltaRho0)/totalMass
232  ).value();
233 
235 
237  << "Surface film: " << type() << nl
238  << " time step continuity errors: sum local = "
239  << sumLocalContErr << ", global = " << globalContErr
240  << ", cumulative = " << cumulativeContErr_ << endl;
241  }
242 }
243 
244 
246 {
247  if (debug)
248  {
249  InfoInFunction << endl;
250  }
251 
252  solve
253  (
255  + fvc::div(phi_)
256  ==
257  - rhoSp_
258  );
259 }
260 
261 
263 {
264  // Push boundary film velocity values into internal field
265  for (label i=0; i<intCoupledPatchIDs_.size(); i++)
266  {
268  const polyPatch& pp = regionMesh().boundaryMesh()[patchi];
271  }
272  Uw_ -= nHat()*(Uw_ & nHat());
274 
275  Us_ = turbulence_->Us();
276 }
277 
278 
280 (
281  const volScalarField& pu,
282  const volScalarField& pp
283 )
284 {
285  if (debug)
286  {
287  InfoInFunction << endl;
288  }
289 
290  // Momentum
292  (
294  + fvm::div(phi_, U_)
295  ==
296  - USp_
297  // - fvm::SuSp(rhoSp_, U_)
298  - rhoSp_*U_
299  + forces_.correct(U_)
300  + turbulence_->Su(U_)
301  );
302 
303  fvVectorMatrix& UEqn = tUEqn.ref();
304 
305  UEqn.relax();
306 
307  if (momentumPredictor_)
308  {
309  solve
310  (
311  UEqn
312  ==
314  (
316  * (
317  regionMesh().magSf()
318  * (
319  fvc::snGrad(pu, "snGrad(p)")
320  + fvc::snGrad(pp, "snGrad(p)")*fvc::interpolate(delta_)
322  )
323  - fvc::flux(rho_*gTan())
324  )
325  )
326  );
327 
328  // Remove any patch-normal components of velocity
329  U_ -= nHat()*(nHat() & U_);
331  }
332 
333  return tUEqn;
334 }
335 
336 
338 (
339  const volScalarField& pu,
340  const volScalarField& pp,
341  const fvVectorMatrix& UEqn
342 )
343 {
344  if (debug)
345  {
346  InfoInFunction << endl;
347  }
348 
349  volScalarField rUA(1.0/UEqn.A());
350  U_ = rUA*UEqn.H();
351 
354 
355  surfaceScalarField phiAdd
356  (
357  "phiAdd",
358  regionMesh().magSf()
359  * (
360  fvc::snGrad(pu, "snGrad(p)")
361  + fvc::snGrad(pp, "snGrad(p)")*fvc::interpolate(delta_)
362  )
363  - fvc::flux(rho_*gTan())
364  );
365  constrainFilmField(phiAdd, 0.0);
366 
367  surfaceScalarField phid
368  (
369  "phid",
370  fvc::flux(U_*rho_) - deltarUAf*phiAdd*rhof
371  );
372  constrainFilmField(phid, 0.0);
373 
374  surfaceScalarField ddrhorUAppf
375  (
376  "deltaCoeff",
377  fvc::interpolate(delta_)*deltarUAf*rhof*fvc::interpolate(pp)
378  );
379 
381 
382  for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
383  {
384  // Film thickness equation
385  fvScalarMatrix deltaEqn
386  (
387  fvm::ddt(rho_, delta_)
388  + fvm::div(phid, delta_)
389  - fvm::laplacian(ddrhorUAppf, delta_)
390  ==
391  - rhoSp_
392  );
393 
394  deltaEqn.solve();
395 
396  if (nonOrth == nNonOrthCorr_)
397  {
398  phiAdd +=
399  fvc::interpolate(pp)
401  * regionMesh().magSf();
402 
403  phi_ == deltaEqn.flux();
404  }
405  }
406 
407  // Bound film thickness by a minimum of zero
408  delta_.max(0.0);
409 
410  // Update U field
411  U_ -= fvc::reconstruct(deltarUAf*phiAdd);
412 
413  // Remove any patch-normal components of velocity
414  U_ -= nHat()*(nHat() & U_);
415 
417 
418  // Continuity check
419  continuityCheck();
420 }
421 
422 
423 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
424 
426 (
427  const word& modelType,
428  const fvMesh& mesh,
429  const dimensionedVector& g,
430  const word& regionType,
431  const bool readFields
432 )
433 :
434  surfaceFilmRegionModel(modelType, mesh, g, regionType),
435 
436  momentumPredictor_(solution().subDict("PISO").lookup("momentumPredictor")),
437  nOuterCorr_(solution().subDict("PISO").lookupOrDefault("nOuterCorr", 1)),
438  nCorr_(readLabel(solution().subDict("PISO").lookup("nCorr"))),
440  (
441  readLabel(solution().subDict("PISO").lookup("nNonOrthCorr"))
442  ),
443 
444  cumulativeContErr_(0.0),
445 
446  deltaSmall_("deltaSmall", dimLength, small),
447  deltaCoLimit_(solution().lookupOrDefault("deltaCoLimit", 1e-4)),
448 
449  rho_
450  (
451  IOobject
452  (
453  "rhof",
454  time().timeName(),
455  regionMesh(),
458  ),
459  regionMesh(),
461  zeroGradientFvPatchScalarField::typeName
462  ),
463  mu_
464  (
465  IOobject
466  (
467  "muf",
468  time().timeName(),
469  regionMesh(),
472  ),
473  regionMesh(),
474  dimensionedScalar(dimPressure*dimTime, 0),
475  zeroGradientFvPatchScalarField::typeName
476  ),
477  sigma_
478  (
479  IOobject
480  (
481  "sigmaf",
482  time().timeName(),
483  regionMesh(),
486  ),
487  regionMesh(),
488  dimensionedScalar(dimMass/sqr(dimTime), 0),
489  zeroGradientFvPatchScalarField::typeName
490  ),
491 
492  delta_
493  (
494  IOobject
495  (
496  "deltaf",
497  time().timeName(),
498  regionMesh(),
501  ),
502  regionMesh()
503  ),
504  alpha_
505  (
506  IOobject
507  (
508  "alpha",
509  time().timeName(),
510  regionMesh(),
513  ),
514  regionMesh(),
516  zeroGradientFvPatchScalarField::typeName
517  ),
518  U_
519  (
520  IOobject
521  (
522  "Uf",
523  time().timeName(),
524  regionMesh(),
527  ),
528  regionMesh()
529  ),
530  Us_
531  (
532  IOobject
533  (
534  "Usf",
535  time().timeName(),
536  regionMesh(),
539  ),
540  U_,
541  zeroGradientFvPatchScalarField::typeName
542  ),
543  Uw_
544  (
545  IOobject
546  (
547  "Uwf",
548  time().timeName(),
549  regionMesh(),
552  ),
553  U_,
554  zeroGradientFvPatchScalarField::typeName
555  ),
556  deltaRho_
557  (
558  IOobject
559  (
560  delta_.name() + "*" + rho_.name(),
561  time().timeName(),
562  regionMesh(),
565  ),
566  regionMesh(),
568  zeroGradientFvPatchScalarField::typeName
569  ),
570 
571  phi_
572  (
573  IOobject
574  (
575  "phi",
576  time().timeName(),
577  regionMesh(),
580  ),
581  regionMesh(),
583  ),
584 
586  (
587  IOobject
588  (
589  "primaryMassTrans",
590  time().timeName(),
591  regionMesh(),
594  ),
595  regionMesh(),
597  zeroGradientFvPatchScalarField::typeName
598  ),
600  (
601  IOobject
602  (
603  "cloudMassTrans",
604  time().timeName(),
605  regionMesh(),
608  ),
609  regionMesh(),
611  zeroGradientFvPatchScalarField::typeName
612  ),
614  (
615  IOobject
616  (
617  "cloudDiameterTrans",
618  time().timeName(),
619  regionMesh(),
622  ),
623  regionMesh(),
625  zeroGradientFvPatchScalarField::typeName
626  ),
627 
628  USp_
629  (
630  IOobject
631  (
632  "USpf",
633  time().timeName(),
634  regionMesh(),
637  ),
638  regionMesh(),
640  (
641  "zero", dimMass*dimVelocity/dimArea/dimTime, Zero
642  ),
643  this->mappedPushedFieldPatchTypes<vector>()
644  ),
645  pSp_
646  (
647  IOobject
648  (
649  "pSpf",
650  time_.timeName(),
651  regionMesh(),
654  ),
655  regionMesh(),
657  this->mappedPushedFieldPatchTypes<scalar>()
658  ),
659  rhoSp_
660  (
661  IOobject
662  (
663  "rhoSpf",
664  time_.timeName(),
665  regionMesh(),
668  ),
669  regionMesh(),
670  dimensionedScalar(dimMass/dimTime/dimArea, 0),
671  this->mappedPushedFieldPatchTypes<scalar>()
672  ),
673 
675  (
676  IOobject
677  (
678  USp_.name(), // must have same name as USp_ to enable mapping
679  time().timeName(),
680  primaryMesh(),
683  ),
684  primaryMesh(),
686  ),
688  (
689  IOobject
690  (
691  pSp_.name(), // must have same name as pSp_ to enable mapping
692  time().timeName(),
693  primaryMesh(),
696  ),
697  primaryMesh(),
699  ),
701  (
702  IOobject
703  (
704  rhoSp_.name(), // must have same name as rhoSp_ to enable mapping
705  time().timeName(),
706  primaryMesh(),
709  ),
710  primaryMesh(),
712  ),
713 
714  UPrimary_
715  (
716  IOobject
717  (
718  "U", // must have same name as U to enable mapping
719  time().timeName(),
720  regionMesh(),
723  ),
724  regionMesh(),
726  this->mappedFieldAndInternalPatchTypes<vector>()
727  ),
728  pPrimary_
729  (
730  IOobject
731  (
732  "p", // must have same name as p to enable mapping
733  time().timeName(),
734  regionMesh(),
737  ),
738  regionMesh(),
740  this->mappedFieldAndInternalPatchTypes<scalar>()
741  ),
743  (
744  IOobject
745  (
746  "rho", // must have same name as rho to enable mapping
747  time().timeName(),
748  regionMesh(),
751  ),
752  regionMesh(),
754  this->mappedFieldAndInternalPatchTypes<scalar>()
755  ),
756  muPrimary_
757  (
758  IOobject
759  (
760  "thermo:mu", // must have same name as mu to enable mapping
761  time().timeName(),
762  regionMesh(),
765  ),
766  regionMesh(),
767  dimensionedScalar(dimPressure*dimTime, 0),
768  this->mappedFieldAndInternalPatchTypes<scalar>()
769  ),
770 
772 
773  availableMass_(regionMesh().nCells(), 0.0),
774 
775  injection_(*this, coeffs_),
776 
777  transfer_(*this, coeffs_),
778 
780 
781  forces_(*this, coeffs_),
782 
783  addedMassTotal_(0.0)
784 {
785  if (readFields)
786  {
788 
789  correctAlpha();
790 
792 
793  deltaRho_ == delta_*rho_;
794 
796  (
797  IOobject
798  (
799  "phi",
800  time().timeName(),
801  regionMesh(),
804  false
805  ),
807  );
808 
809  phi_ == phi0;
810  }
811 }
812 
813 
814 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
815 
817 {}
818 
819 
820 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
821 
823 (
824  const label patchi,
825  const label facei,
826  const scalar massSource,
827  const vector& momentumSource,
828  const scalar pressureSource,
829  const scalar energySource
830 )
831 {
832  if (debug)
833  {
835  << "\nSurface film: " << type() << ": adding to film source:" << nl
836  << " mass = " << massSource << nl
837  << " momentum = " << momentumSource << nl
838  << " pressure = " << pressureSource << endl;
839  }
840 
841  rhoSpPrimary_.boundaryFieldRef()[patchi][facei] -= massSource;
842  USpPrimary_.boundaryFieldRef()[patchi][facei] -= momentumSource;
843  pSpPrimary_.boundaryFieldRef()[patchi][facei] -= pressureSource;
844 
845  addedMassTotal_ += massSource;
846 }
847 
848 
850 {
851  if (debug)
852  {
853  InfoInFunction << endl;
854  }
855 
857 
859 
861 
863 
864  // Reset transfer fields
865  availableMass_ = mass();
869 }
870 
871 
873 {
874  if (debug)
875  {
876  InfoInFunction << endl;
877  }
878 
879  // Update film coverage indicator
880  correctAlpha();
881 
882  // Update film wall and surface velocities
884 
885  // Update sub-models to provide updated source contributions
886  updateSubmodels();
887 
888  // Solve continuity for deltaRho_
889  solveContinuity();
890 
891  // Implicit pressure source coefficient - constant
892  tmp<volScalarField> tpp(this->pp());
893 
894  for (int oCorr=1; oCorr<=nOuterCorr_; oCorr++)
895  {
896  // Explicit pressure source contribution - varies with delta_
897  tmp<volScalarField> tpu(this->pu());
898 
899  // Solve for momentum for U_
900  tmp<fvVectorMatrix> UEqn = solveMomentum(tpu(), tpp());
901 
902  // Film thickness correction loop
903  for (int corr=1; corr<=nCorr_; corr++)
904  {
905  // Solve thickness for delta_
906  solveThickness(tpu(), tpp(), UEqn());
907  }
908  }
909 
910  // Update deltaRho_ with new delta_
911  deltaRho_ == delta_*rho_;
912 
913  // Reset source terms for next time integration
915 }
916 
917 
919 {
920  scalar CoNum = 0.0;
921 
922  if (regionMesh().nInternalFaces() > 0)
923  {
924  const scalarField sumPhi
925  (
926  fvc::surfaceSum(mag(phi_))().primitiveField()
927  / (deltaRho_.primitiveField() + rootVSmall)
928  );
929 
930  forAll(delta_, i)
931  {
932  if (delta_[i] > deltaCoLimit_)
933  {
934  CoNum = max(CoNum, sumPhi[i]/(delta_[i]*magSf()[i]));
935  }
936  }
937 
938  CoNum *= 0.5*time_.deltaTValue();
939  }
940 
941  reduce(CoNum, maxOp<scalar>());
942 
943  Info<< "Film max Courant number: " << CoNum << endl;
944 
945  return CoNum;
946 }
947 
948 
950 {
951  return U_;
952 }
953 
954 
956 {
957  return Us_;
958 }
959 
960 
962 {
963  return Uw_;
964 }
965 
966 
968 {
969  return deltaRho_;
970 }
971 
972 
974 {
975  return phi_;
976 }
977 
978 
980 {
981  return rho_;
982 }
983 
984 
986 {
988  << "T field not available for " << type() << abort(FatalError);
989 
990  return volScalarField::null();
991 }
992 
993 
995 {
997  << "Ts field not available for " << type() << abort(FatalError);
998 
999  return volScalarField::null();
1000 }
1001 
1002 
1004 {
1006  << "Tw field not available for " << type() << abort(FatalError);
1007 
1008  return volScalarField::null();
1009 }
1010 
1011 
1013 {
1015  << "hs field not available for " << type() << abort(FatalError);
1016 
1017  return volScalarField::null();
1018 }
1019 
1020 
1022 {
1024  << "Cp field not available for " << type() << abort(FatalError);
1025 
1026  return volScalarField::null();
1027 }
1028 
1029 
1031 {
1033  << "kappa field not available for " << type() << abort(FatalError);
1034 
1035  return volScalarField::null();
1036 }
1037 
1038 
1040 {
1041  return primaryMassTrans_;
1042 }
1043 
1044 
1046 {
1047  return cloudMassTrans_;
1048 }
1049 
1050 
1052 {
1053  return cloudDiameterTrans_;
1054 }
1055 
1056 
1058 {
1059  Info<< "\nSurface film: " << type() << endl;
1060 
1061  const scalarField& deltaInternal = delta_;
1062  const vectorField& Uinternal = U_;
1063  scalar addedMassTotal = 0.0;
1064  outputProperties().readIfPresent("addedMassTotal", addedMassTotal);
1065  addedMassTotal += returnReduce(addedMassTotal_, sumOp<scalar>());
1066 
1067  Info<< indent << "added mass = " << addedMassTotal << nl
1068  << indent << "current mass = "
1069  << gSum((deltaRho_*magSf())()) << nl
1070  << indent << "min/max(mag(U)) = " << gMin(mag(Uinternal)) << ", "
1071  << gMax(mag(Uinternal)) << nl
1072  << indent << "min/max(delta) = " << gMin(deltaInternal) << ", "
1073  << gMax(deltaInternal) << nl
1074  << indent << "coverage = "
1075  << gSum(alpha_.primitiveField()*magSf())/gSum(magSf()) << nl;
1076 
1077  injection_.info(Info);
1078  transfer_.info(Info);
1079 }
1080 
1081 
1083 {
1085  (
1086  typeName + ":Srho",
1087  primaryMesh(),
1088  dimensionedScalar(dimMass/dimVolume/dimTime, 0)
1089  );
1090 }
1091 
1092 
1095  const label i
1096 ) const
1097 {
1099  (
1100  typeName + ":Srho(" + Foam::name(i) + ")",
1101  primaryMesh(),
1102  dimensionedScalar(dimMass/dimVolume/dimTime, 0)
1103  );
1104 }
1105 
1106 
1108 {
1110  (
1111  typeName + ":Sh",
1112  primaryMesh(),
1114  );
1115 }
1116 
1117 
1118 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1119 
1120 } // End namespace surfaceFilmModels
1121 } // End namespace regionModels
1122 } // End namespace Foam
1123 
1124 // ************************************************************************* //
virtual bool read()
Read control parameters from dictionary.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
autoPtr< filmThermoModel > filmThermo_
Film thermo model.
virtual bool read()
Read control parameters from dictionary.
surfaceFilmRegionModel(const word &modelType, const fvMesh &mesh, const dimensionedVector &g, const word &regionType)
Construct from type name, mesh and gravity vector.
virtual const volVectorField & U() const
Return the film velocity [m/s].
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
autoPtr< filmTurbulenceModel > turbulence_
Turbulence model.
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
const word & name() const
Return name.
Definition: IOobject.H:295
virtual void correct(scalarField &availableMass, volScalarField &massToInject, volScalarField &diameterToInject)
Correct.
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:226
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
virtual const volScalarField & deltaRho() const
Return the film thickness*density (helper field) [kg/m^3].
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
volScalarField primaryMassTrans_
Film mass available for transfer to the primary region.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Type gMin(const FieldField< Field, Type > &f)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Reconstruct volField from a face flux field.
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:760
virtual const volVectorField & Us() const
Return the film surface velocity [m/s].
const dimensionedScalar phi0
Magnetic flux quantum: default SI units: [Wb].
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:497
static autoPtr< filmThermoModel > New(surfaceFilmRegionModel &film, const dictionary &dict)
Return a reference to the selected phase change model.
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
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.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
scalar addedMassTotal_
Cumulative mass added via sources [kg].
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
virtual void info(Ostream &os)
Provide some info.
virtual tmp< volScalarField > primaryMassTrans() const
Return mass transfer source - Eulerian phase only.
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
virtual void solveThickness(const volScalarField &pu, const volScalarField &pp, const fvVectorMatrix &UEqn)
Solve coupled velocity-thickness equations.
virtual void evolveRegion()
Evolve the film equations.
Macros for easy insertion into run-time selection tables.
virtual void solveContinuity()
Solve continuity equation.
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:821
virtual tmp< volScalarField > pp()
Implicit pressure source coefficient.
tmp< fvVectorMatrix > correct(volVectorField &U)
Return (net) force system.
Definition: forceList.C:78
volScalarField deltaRho_
Film thickness*density (helper field) [kg/m^2].
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:699
virtual void transferPrimaryRegionSourceFields()
Transfer source fields from the primary region to the film region.
virtual void correct(scalarField &availableMass, volScalarField &massToTransfer)
Correct kinematic transfers.
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const dimensionSet & dimensions() const
Return dimensions.
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionModelI.H:31
Type gSum(const FieldField< Field, Type > &f)
virtual tmp< fvVectorMatrix > solveMomentum(const volScalarField &pu, const volScalarField &pp)
Solve for film velocity.
virtual const volScalarField & cloudMassTrans() const
Return the film mass available for transfer to cloud.
virtual void preEvolveRegion()
Pre-evolve region.
Definition: regionModel.C:519
static autoPtr< filmTurbulenceModel > New(surfaceFilmRegionModel &film, const dictionary &dict)
Return a reference to the selected injection model.
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:315
tmp< volVectorField > gTan() const
Return the gravity tangential component contributions.
tmp< volScalarField > mass() const
Return the current film mass.
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.
Calculate the laplacian of the given field.
tmp< volScalarField > gNormClipped() const
Return the gravity normal-to-patch component contribution.
const IOdictionary & outputProperties() const
Return const access to the output properties dictionary.
Definition: regionModelI.H:110
scalar sumLocalContErr
const dictionary & solution() const
Return the solution dictionary.
Definition: regionModelI.H:103
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
virtual const volVectorField & nHat() const
Return the patch normal vectors.
virtual const surfaceScalarField & phi() const
Return the film flux [kg m/s].
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:61
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
static const zero Zero
Definition: zero.H:97
const dimensionSet dimPressure
label nNonOrthCorr_
Number of non-orthogonal correctors.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
label readLabel(Istream &is)
Definition: label.H:64
virtual const volScalarField & kappa() const
Return the film thermal conductivity [W/m/K].
virtual tmp< volScalarField > pu()
Explicit pressure source contribution.
Calculate the divergence of the given field.
scalar CoNum
dimensionedScalar totalMass
Definition: continuityErrs.H:4
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
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:265
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Type gMax(const FieldField< Field, Type > &f)
virtual void correctAlpha()
Correct film coverage field.
virtual const volScalarField & hs() const
Return the film surface enthalpy [J/kg].
virtual void resetPrimaryRegionSourceTerms()
Reset source term fields.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimEnergy
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.
const dimensionSet dimDensity
dictionary coeffs_
Model coefficients dictionary.
Definition: regionModel.H:99
static const GeometricField< scalar, fvPatchField, volMesh > & null()
Return a null geometric field.
virtual void updateSubmodels()
Update the film sub-models.
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m^2].
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:738
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual const volVectorField & Uw() const
Return the film wall velocity [m/s].
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
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
const Time & time_
Reference to the time database.
Definition: regionModel.H:84
void max(const dimensioned< Type > &)
surfaceScalarField phi_
Mass flux (includes film thickness) [kg m/s].
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
volScalarField alpha_
Film coverage indicator, 1 = covered, 0 = uncovered [].
A List with indirect addressing.
Definition: fvMatrix.H:106
virtual void updateSurfaceVelocities()
Update film surface velocities.
virtual const volScalarField & Tw() const
Return the film wall temperature [K].
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
void correctBoundaryConditions()
Correct boundary field.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
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:863
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
virtual const volScalarField & rho() const
Return the film density [kg/m^3].
volScalarField cloudMassTrans_
Film mass available for transfer to cloud.
void constrainFilmField(Type &field, const typename Type::cmptType &value)
Constrain a film region master/slave boundaries of a field to a.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
Definition: regionModel.H:111
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
virtual const volScalarField & Cp() const
Return the film specific heat capacity [J/kg/K].
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
virtual const volScalarField & T() const
Return the film mean temperature [K].
defineTypeNameAndDebug(kinematicSingleLayer, 0)
virtual tmp< volScalarField::Internal > Sh() const
Return enthalpy source - Eulerian phase only.
virtual scalar CourantNumber() const
Courant number evaluation.
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:53
scalar deltaCoLimit_
Film thickness above which Courant number calculation in valid.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
virtual void correctThermoFields()
Correct the thermo fields.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583
#define InfoInFunction
Report an information message using Foam::Info.
const dimensionSet dimVelocity