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-2020 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 
43 #include "mappedWallPolyPatch.H"
44 #include "mapDistribute.H"
45 #include "filmThermoModel.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 namespace regionModels
52 {
53 namespace surfaceFilmModels
54 {
55 
56 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
57 
58 defineTypeNameAndDebug(kinematicSingleLayer, 0);
59 
60 addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh);
61 
62 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
63 
65 {
67 }
68 
69 
71 {
72  rho_ == filmThermo_->rho();
73  mu_ == filmThermo_->mu();
74  sigma_ == filmThermo_->sigma();
75 }
76 
77 
79 {
81 
83  USpPrimary_ == Zero;
84  pSpPrimary_ == Zero;
85 }
86 
87 
89 {
91 
92  // Update fields from primary region via direct mapped
93  // (coupled) boundary conditions
98 }
99 
100 
102 {
104 
105  volScalarField::Boundary& rhoSpPrimaryBf =
107 
108  volVectorField::Boundary& USpPrimaryBf =
110 
111  volScalarField::Boundary& pSpPrimaryBf =
113 
114  // Convert accumulated source terms into per unit area per unit time
115  const scalar deltaT = time_.deltaTValue();
117  {
118  scalarField rpriMagSfdeltaT
119  (
120  (1/deltaT)
121  /primaryMesh().magSf().boundaryField()[patchi]
122  );
123 
124  rhoSpPrimaryBf[patchi] *= rpriMagSfdeltaT;
125  USpPrimaryBf[patchi] *= rpriMagSfdeltaT;
126  pSpPrimaryBf[patchi] *= rpriMagSfdeltaT;
127  }
128 
129  // Retrieve the source fields from the primary region
130  toRegion(rhoSp_, rhoSpPrimaryBf);
131  rhoSp_.field() /= VbyA();
132  toRegion(USp_, USpPrimaryBf);
133  USp_.field() /= VbyA();
134  toRegion(pSp_, pSpPrimaryBf);
135 
136  // update addedMassTotal counter
137  if (time().writeTime())
138  {
139  scalar addedMassTotal = 0;
140  outputProperties().readIfPresent("addedMassTotal", addedMassTotal);
141  addedMassTotal += returnReduce(addedMassTotal_, sumOp<scalar>());
142  outputProperties().add("addedMassTotal", addedMassTotal, true);
143  addedMassTotal_ = 0;
144  }
145 }
146 
147 
149 {
150  return -fvc::laplacian(sigma_, delta_);
151 }
152 
153 
155 {
157  (
159  (
160  "pSp",
161  regionMesh(),
163  zeroGradientFvPatchScalarField::typeName
164  )
165  );
166 
167  tpSp.ref().primitiveFieldRef() = pSp_;
168  tpSp.ref().correctBoundaryConditions();
169 
170  return volScalarField::New
171  (
172  IOobject::modelName("pe", typeName),
173  pPrimary_ // Pressure (mapped from primary region)
174  - tpSp // Accumulated particle impingement
175  );
176 }
177 
178 
180 {
181  return
183  (
184  max(nHat() & -g_, dimensionedScalar(g_.dimensions(), 0))*VbyA()
186 }
187 
188 
190 {
191  return
193  (
194  max(nHat() & -g_, dimensionedScalar(g_.dimensions(), 0))*VbyA()
195  )*fvc::snGrad(rho_);
196 }
197 
198 
200 {
202 }
203 
204 
206 {
208 
209  // Update injection model - mass returned is mass available for injection
211 
212  // Update transfer model - mass returned is mass available for transfer
214 
215  // Update mass source field
217 
218  turbulence_->correct();
219 }
220 
221 
223 {
225 
227 
228  // Bound film volume fraction
229  alpha_.max(0);
230 
231  delta_ == alpha_*VbyA();
232 
233  // Update continuity error caused by the delta_ bounding
235 }
236 
237 
239 {
241 }
242 
243 
245 {
247 
248  if (totalMass.value() > small)
249  {
250  const volScalarField::Internal massErr
251  (
253  );
254 
255  const scalar sumLocalContErr =
256  (
258  ).value();
259 
260  const scalar globalContErr =
261  (
263  ).value();
264 
266 
267  Info<< "time step continuity errors: sum local = "
268  << sumLocalContErr
269  << ", global = " << globalContErr
270  << ", cumulative = " << cumulativeContErr_
271  << endl;
272  }
273 }
274 
275 
277 {
278  // Push boundary film velocity values into internal field
279  for (label i=0; i<intCoupledPatchIDs_.size(); i++)
280  {
281  const label patchi = intCoupledPatchIDs_[i];
282  const polyPatch& pp = regionMesh().boundaryMesh()[patchi];
285  }
286 
287  Uw_ -= nHat()*(Uw_ & nHat());
288 
289  Us_ = turbulence_->Us();
290 }
291 
292 
294 (
295  const volScalarField& pc,
296  const volScalarField& pe
297 )
298 {
300 
301  const volScalarField::Internal rVDt
302  (
303  1/(time().deltaT()*regionMesh().V())
304  );
305 
306  // Momentum
308  (
311  ==
312  - USp_
313 
314  // Temporary treatment for the loss of momentum due to mass loss
315  // These transfers are not currently included in USp_
317 
318  + forces_.correct(U_)
319  + turbulence_->Su(U_)
320  );
321 
322  fvVectorMatrix& UEqn = tUEqn.ref();
323 
324  UEqn.relax();
325 
327  {
329 
330  solve
331  (
332  UEqn
333  ==
335  (
337  (
338  alphaf
339  *(
340  (
341  fvc::snGrad(pe + pc, "snGrad(p)")
342  + gGradRho()*alphaf
343  + rhog()*fvc::snGrad(alpha_)
344  )*regionMesh().magSf()
345  - fvc::interpolate(rho_)*(g_ & regionMesh().Sf())
346  ), 0
347  )
348  )
349  );
350 
351  // Remove any patch-normal components of velocity
352  U_ -= nHat()*(nHat() & U_);
353 
355  }
356 
357  return tUEqn;
358 }
359 
360 
362 (
363  const fvVectorMatrix& UEqn,
364  const volScalarField& pc,
365  const volScalarField& pe
366 )
367 {
369 
370  const volScalarField rAU(1/UEqn.A());
371  const volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U_, alpha_));
372 
375  const surfaceScalarField alpharAUf(fvc::interpolate(alpha_*rAU));
376  const surfaceScalarField rhogf(rhog());
377 
378  const surfaceScalarField phiu
379  (
380  "phiu",
381  (
383  (
384  (
385  fvc::snGrad(pe + pc, "snGrad(p)")
386  + gGradRho()*alphaf
387  )*regionMesh().magSf()
388  - rhof*(g_ & regionMesh().Sf()),
389  0
390  )
391  )
392  );
393 
394  surfaceScalarField phid
395  (
396  "phid",
397  constrainFilmField(rhof*(fvc::flux(HbyA) - alpharAUf*phiu), 0)
398  );
399 
400  const surfaceScalarField ddrhorAUrhogf
401  (
402  "alphaCoeff",
403  alphaf*rhof*alpharAUf*rhogf
404  );
405 
407 
408  while (pimple_.correctNonOrthogonal())
409  {
410  // Film thickness equation
411  fvScalarMatrix alphaEqn
412  (
414  + fvm::div(phid, alpha_)
415  - fvm::laplacian(ddrhorAUrhogf, alpha_)
416  ==
417  -rhoSp_
418  );
419 
420  alphaEqn.solve();
421 
423  {
424  phi_ == alphaEqn.flux();
425 
426  const surfaceScalarField phiGradAlpha
427  (
429  (
430  rhogf*fvc::snGrad(alpha_)*regionMesh().magSf(),
431  0
432  )
433  );
434 
435  phiU_ = constrainFilmField(phid/rhof - alpharAUf*phiGradAlpha, 0);
436 
437  // Update U field
438  U_ = HbyA - rAU*fvc::reconstruct(alphaf*(phiu + phiGradAlpha));
439 
440  // Remove any patch-normal components of velocity
441  U_ -= nHat()*(nHat() & U_);
442 
444  }
445  }
446 
447  // Bound film volume fraction
448  alpha_.max(0);
449 
450  delta_ == alpha_*VbyA();
451 
453 
454  // Continuity check
455  continuityCheck();
456 }
457 
458 
459 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
460 
462 (
463  const word& modelType,
464  const fvMesh& mesh,
465  const dimensionedVector& g,
466  const word& regionType,
467  const bool readFields
468 )
469 :
470  surfaceFilmRegionModel(modelType, mesh, g, regionType),
471  pimple_(regionMesh()),
472 
474 
475  deltaSmall_("deltaSmall", dimLength, small),
476  deltaCoLimit_(solution().lookupOrDefault("deltaCoLimit", 1e-4)),
477 
478  rho_
479  (
480  IOobject
481  (
482  "rho",
483  time().timeName(),
484  regionMesh(),
487  ),
488  regionMesh(),
490  zeroGradientFvPatchScalarField::typeName
491  ),
492 
493  mu_
494  (
495  IOobject
496  (
497  "mu",
498  time().timeName(),
499  regionMesh(),
502  ),
503  regionMesh(),
504  dimensionedScalar(dimPressure*dimTime, 0),
505  zeroGradientFvPatchScalarField::typeName
506  ),
507 
508  sigma_
509  (
510  IOobject
511  (
512  "sigma",
513  time().timeName(),
514  regionMesh(),
517  ),
518  regionMesh(),
519  dimensionedScalar(dimMass/sqr(dimTime), 0),
520  zeroGradientFvPatchScalarField::typeName
521  ),
522 
523  delta_
524  (
525  IOobject
526  (
527  "delta",
528  time().timeName(),
529  regionMesh(),
532  ),
533  regionMesh()
534  ),
535 
536  alpha_
537  (
538  IOobject
539  (
540  "alpha",
541  time().timeName(),
542  regionMesh(),
545  ),
546  delta_/VbyA(),
547  delta_.boundaryField().types()
548  ),
549 
550  U_
551  (
552  IOobject
553  (
554  "U",
555  time().timeName(),
556  regionMesh(),
559  ),
560  regionMesh()
561  ),
562 
563  Us_
564  (
565  IOobject
566  (
567  "Us",
568  time().timeName(),
569  regionMesh(),
572  ),
573  U_
574  ),
575 
576  Uw_
577  (
578  IOobject
579  (
580  "Uw",
581  time().timeName(),
582  regionMesh(),
585  ),
586  U_
587  ),
588 
589  phi_
590  (
591  IOobject
592  (
593  "phi",
594  time().timeName(),
595  regionMesh(),
598  ),
599  regionMesh(),
600  dimensionedScalar(dimMass/dimTime, 0)
601  ),
602 
603  phiU_
604  (
605  IOobject
606  (
607  "phiU",
608  time().timeName(),
609  regionMesh()
610  ),
611  regionMesh(),
612  dimensionedScalar(dimVolume/dimTime, 0)
613  ),
614 
616  (
617  IOobject
618  (
619  "continuityErr",
620  time().timeName(),
621  regionMesh()
622  ),
623  regionMesh(),
625  ),
626 
627  coverage_
628  (
629  IOobject
630  (
631  "coverage",
632  time().timeName(),
633  regionMesh(),
636  ),
637  regionMesh(),
639  zeroGradientFvPatchScalarField::typeName
640  ),
641 
643  (
644  IOobject
645  (
646  "primaryMassTrans",
647  time().timeName(),
648  regionMesh(),
651  ),
652  regionMesh(),
654  zeroGradientFvPatchScalarField::typeName
655  ),
656 
658  (
659  IOobject
660  (
661  "cloudMassTrans",
662  time().timeName(),
663  regionMesh(),
666  ),
667  regionMesh(),
669  zeroGradientFvPatchScalarField::typeName
670  ),
671 
673  (
674  IOobject
675  (
676  "cloudDiameterTrans",
677  time().timeName(),
678  regionMesh(),
681  ),
682  regionMesh(),
684  zeroGradientFvPatchScalarField::typeName
685  ),
686 
687  rhoSp_
688  (
689  IOobject
690  (
691  "rhoSp",
692  time_.timeName(),
693  regionMesh(),
696  ),
697  regionMesh(),
698  dimensionedScalar(dimDensity/dimTime, 0)
699  ),
700 
701  USp_
702  (
703  IOobject
704  (
705  "USp",
706  time().timeName(),
707  regionMesh(),
710  ),
711  regionMesh(),
713  ),
714 
715  pSp_
716  (
717  IOobject
718  (
719  "pSp",
720  time_.timeName(),
721  regionMesh(),
724  ),
725  regionMesh(),
727  ),
728 
730  (
731  IOobject
732  (
733  rhoSp_.name(),
734  time().timeName(),
735  primaryMesh(),
738  ),
739  primaryMesh(),
741  ),
742 
744  (
745  IOobject
746  (
747  USp_.name(),
748  time().timeName(),
749  primaryMesh(),
752  ),
753  primaryMesh(),
755  ),
756 
758  (
759  IOobject
760  (
761  pSp_.name(),
762  time().timeName(),
763  primaryMesh(),
766  ),
767  primaryMesh(),
769  ),
770 
771  UPrimary_
772  (
773  IOobject
774  (
775  "U", // must have same name as U to enable mapping
776  time().timeName(),
777  regionMesh(),
780  ),
781  regionMesh(),
783  this->mappedFieldAndInternalPatchTypes<vector>()
784  ),
785 
786  pPrimary_
787  (
788  IOobject
789  (
790  "p", // must have same name as p to enable mapping
791  time().timeName(),
792  regionMesh(),
795  ),
796  regionMesh(),
798  this->mappedFieldAndInternalPatchTypes<scalar>()
799  ),
800 
802  (
803  IOobject
804  (
805  "rho", // must have same name as rho to enable mapping
806  time().timeName(),
807  regionMesh(),
810  ),
811  regionMesh(),
813  this->mappedFieldAndInternalPatchTypes<scalar>()
814  ),
815 
816  muPrimary_
817  (
818  IOobject
819  (
820  "thermo:mu", // must have same name as mu to enable mapping
821  time().timeName(),
822  regionMesh(),
825  ),
826  regionMesh(),
827  dimensionedScalar(dimPressure*dimTime, 0),
828  this->mappedFieldAndInternalPatchTypes<scalar>()
829  ),
830 
832 
833  availableMass_(regionMesh().nCells(), 0),
834 
835  injection_(*this, coeffs_),
836 
837  transfer_(*this, coeffs_),
838 
840 
841  forces_(*this, coeffs_),
842 
843  addedMassTotal_(0)
844 {
845  alpha_ == delta_/VbyA();
846 
847  if (readFields)
848  {
850 
851  correctCoverage();
852 
854 
856  (
857  IOobject
858  (
859  "phi",
860  time().timeName(),
861  regionMesh(),
864  false
865  ),
867  );
868 
869  phi_ == phi;
870  phiU_ = fvc::flux(U_);
871  }
872 }
873 
874 
875 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
876 
878 {}
879 
880 
881 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
882 
884 (
885  const label patchi,
886  const label facei,
887  const scalar massSource,
888  const vector& momentumSource,
889  const scalar pressureSource,
890  const scalar energySource
891 )
892 {
894  << "\nSurface film: " << type() << ": adding to film source:" << nl
895  << " mass = " << massSource << nl
896  << " momentum = " << momentumSource << nl
897  << " pressure = " << pressureSource << endl;
898 
899  rhoSpPrimary_.boundaryFieldRef()[patchi][facei] -= massSource;
900  USpPrimary_.boundaryFieldRef()[patchi][facei] -= momentumSource;
901  pSpPrimary_.boundaryFieldRef()[patchi][facei] -= pressureSource;
902 
903  addedMassTotal_ += massSource;
904 }
905 
906 
908 {
910 
912 
914 
916 
918 
919  // Reset transfer fields
920  availableMass_ = mass();
924 }
925 
926 
928 {
930 
931  // Update film coverage indicator
932  correctCoverage();
933 
934  // Update film wall and surface velocities
936 
937  // Predict delta_ from continuity
938  predictDelta();
939 
940  // Update sub-models to provide updated source contributions
941  updateSubmodels();
942 
943  // Predict delta_ from continuity with updated source
944  predictDelta();
945 
946  // Capillary pressure
947  const volScalarField pc(this->pc());
948 
949  while (pimple_.loop())
950  {
951  // External pressure
952  const volScalarField pe(this->pe());
953 
954  // Solve for momentum
955  const fvVectorMatrix UEqn(solveMomentum(pc, pe));
956 
957  // Film thickness correction loop
958  while (pimple_.correct())
959  {
960  solveAlpha(UEqn, pc, pe);
961  }
962  }
963 
964  // Reset source terms for next time integration
966 }
967 
968 
970 {
971  const scalarField sumPhi(fvc::surfaceSum(mag(phiU_))().primitiveField());
972 
973  const scalar CoNum =
974  0.5*gMax(sumPhi/regionMesh().V().field())*time_.deltaTValue();
975 
976  Info<< "Film max Courant number: " << CoNum << endl;
977 
978  return CoNum;
979 }
980 
981 
983 {
984  return primaryMassTrans_;
985 }
986 
987 
989 {
990  return cloudMassTrans_;
991 }
992 
993 
995 {
996  return cloudDiameterTrans_;
997 }
998 
999 
1001 {
1002  Info<< "\nSurface film: " << type() << endl;
1003 
1004  const scalarField& deltaInternal = delta_;
1005  const vectorField& Uinternal = U_;
1006  scalar addedMassTotal = 0;
1007  outputProperties().readIfPresent("addedMassTotal", addedMassTotal);
1008  addedMassTotal += returnReduce(addedMassTotal_, sumOp<scalar>());
1009 
1010  Info<< indent << "added mass = " << addedMassTotal << nl
1011  << indent << "current mass = "
1012  << gSum((delta_*rho_*magSf())()) << nl
1013  << indent << "min/max(mag(U)) = " << gMin(mag(Uinternal)) << ", "
1014  << gMax(mag(Uinternal)) << nl
1015  << indent << "min/max(delta) = " << gMin(deltaInternal) << ", "
1016  << gMax(deltaInternal) << nl
1017  << indent << "coverage = "
1018  << gSum(coverage_.primitiveField()*magSf())/gSum(magSf()) << nl;
1019 
1020  injection_.info(Info);
1021  transfer_.info(Info);
1022 }
1023 
1024 
1026 {
1028  (
1029  IOobject::modelName("Srho", typeName),
1030  primaryMesh(),
1031  dimensionedScalar(dimMass/dimVolume/dimTime, 0)
1032  );
1033 }
1034 
1035 
1038  const label i
1039 ) const
1040 {
1042  (
1043  IOobject::modelName("Srho(" + Foam::name(i) + ")", typeName),
1044  primaryMesh(),
1045  dimensionedScalar(dimMass/dimVolume/dimTime, 0)
1046  );
1047 }
1048 
1049 
1051 {
1053  (
1054  IOobject::modelName("Sh", typeName),
1055  primaryMesh(),
1057  );
1058 }
1059 
1060 
1061 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1062 
1063 } // End namespace surfaceFilmModels
1064 } // End namespace regionModels
1065 } // End namespace Foam
1066 
1067 // ************************************************************************* //
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.
tmp< volScalarField::Internal > mass() const
Return the current film mass.
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].
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:303
virtual void correct(scalarField &availableMass, volScalarField &massToInject, volScalarField &diameterToInject)
Correct.
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
bool correctNonOrthogonal()
Non-orthogonal corrector loop.
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.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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:174
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:164
Reconstruct volField from a face flux field.
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:760
Calculate the matrix for the laplacian of the field.
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.
volVectorField::Internal USp_
Momentum [kg/m/s^2].
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
volVectorField::Internal Us_
Velocity - surface [m/s].
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:622
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 void info(Ostream &os)
Provide some info.
virtual tmp< volScalarField > primaryMassTrans() const
Return mass transfer source - Eulerian phase only.
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
tmp< surfaceScalarField > rhog() const
Hydrostatic pressure coefficient.
virtual void evolveRegion()
Evolve the film equations.
Macros for easy insertion into run-time selection tables.
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1056
tmp< fvVectorMatrix > correct(volVectorField &U)
Return (net) force system.
Definition: forceList.C:78
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.
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 const volScalarField & cloudMassTrans() const
Return the film mass available for transfer to cloud.
virtual void preEvolveRegion()
Pre-evolve region.
Definition: regionModel.C:514
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:321
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.
#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: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
Calculate the matrix for the first temporal derivative.
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
bool loop()
Pimple loop.
Definition: pimpleControl.C:81
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.
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
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].
const dimensionSet dimPressure
Calculate the divergence of the given field.
scalar CoNum
dimensionedScalar totalMass
Definition: continuityErrs.H:4
bool correct()
Piso loop.
Definition: pisoControl.C:80
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
Type gMax(const FieldField< Field, Type > &f)
virtual void resetPrimaryRegionSourceTerms()
Reset source term fields.
const Field< Type > & field() const
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevTau(U)==fvOptions(rho, U))
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:97
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:738
label patchi
const volScalarField::Internal & magSf() const
Return the face area magnitudes [m^2].
virtual void predictDelta()
Predict delta_ from the continuity equation.
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,.
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:82
void max(const dimensioned< Type > &)
volVectorField & HbyA
Definition: pEqn.H:13
fvVectorMatrix & UEqn
Definition: UEqn.H:11
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
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.
virtual void updateSurfaceVelocities()
Update film surface velocities.
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...
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const dimensionSet & dimensions() const
Return const reference to dimensions.
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:856
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
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.
const dimensionedVector & g_
Acceleration due to gravity [m/s^2].
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
rDeltaTY field()
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
Definition: regionModel.H:109
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
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
tmp< Type > constrainFilmField(const tmp< Type > &field, const typename Type::cmptType &value)
Constrain a film region master/slave boundaries of a field to a.
static autoPtr< filmMomentumTransportModel > New(surfaceFilmRegionModel &film, const dictionary &dict)
Return a reference to the selected injection model.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
defineTypeNameAndDebug(kinematicSingleLayer, 0)
virtual tmp< volScalarField::Internal > Sh() const
Return enthalpy source - Eulerian phase only.
virtual scalar CourantNumber() const
Courant number evaluation.
autoPtr< filmMomentumTransportModel > turbulence_
Turbulence model.
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
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.
virtual void correctThermoFields()
Correct the thermo fields.
const dimensionSet dimVelocity