multiphaseMixture.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 "multiphaseMixture.H"
27 #include "alphaContactAngleFvPatchScalarField.H"
28 #include "unitConversion.H"
29 #include "Time.H"
30 #include "subCycle.H"
31 #include "MULES.H"
32 #include "surfaceInterpolate.H"
33 #include "fvcGrad.H"
34 #include "fvcSnGrad.H"
35 #include "fvcDiv.H"
36 #include "fvcFlux.H"
37 
38 
39 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40 
41 void Foam::multiphaseMixture::calcAlphas()
42 {
43  scalar level = 0.0;
44  alphas_ == 0.0;
45 
46  forAllIter(PtrDictionary<phase>, phases_, iter)
47  {
48  alphas_ += level*iter();
49  level += 1.0;
50  }
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
57 (
58  const volVectorField& U,
59  const surfaceScalarField& phi
60 )
61 :
63  (
64  IOobject
65  (
66  "transportProperties",
67  U.time().constant(),
68  U.db(),
71  )
72  ),
73 
74  phases_(lookup("phases"), phase::iNew(U, phi)),
75 
76  mesh_(U.mesh()),
77  U_(U),
78  phi_(phi),
79 
80  rhoPhi_
81  (
82  IOobject
83  (
84  "rhoPhi",
85  mesh_.time().timeName(),
86  mesh_,
89  ),
90  mesh_,
92  ),
93 
94  alphas_
95  (
96  IOobject
97  (
98  "alphas",
99  mesh_.time().timeName(),
100  mesh_,
101  IOobject::NO_READ,
103  ),
104  mesh_,
106  ),
107 
108  nu_
109  (
110  IOobject
111  (
112  "nu",
113  mesh_.time().timeName(),
114  mesh_
115  ),
116  mu()/rho()
117  ),
118 
119  sigmas_(lookup("sigmas")),
120  dimSigma_(1, 0, -2, 0, 0),
121  deltaN_
122  (
123  "deltaN",
124  1e-8/pow(average(mesh_.V()), 1.0/3.0)
125  )
126 {
127  calcAlphas();
128  alphas_.write();
129 }
130 
131 
132 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
133 
136 {
137  PtrDictionary<phase>::const_iterator iter = phases_.begin();
138 
139  tmp<volScalarField> trho = iter()*iter().rho();
140  volScalarField& rho = trho.ref();
141 
142  for (++iter; iter != phases_.end(); ++iter)
143  {
144  rho += iter()*iter().rho();
145  }
146 
147  return trho;
148 }
149 
150 
152 Foam::multiphaseMixture::rho(const label patchi) const
153 {
154  PtrDictionary<phase>::const_iterator iter = phases_.begin();
155 
156  tmp<scalarField> trho = iter().boundaryField()[patchi]*iter().rho().value();
157  scalarField& rho = trho.ref();
158 
159  for (++iter; iter != phases_.end(); ++iter)
160  {
161  rho += iter().boundaryField()[patchi]*iter().rho().value();
162  }
163 
164  return trho;
165 }
166 
167 
170 {
171  PtrDictionary<phase>::const_iterator iter = phases_.begin();
172 
173  tmp<volScalarField> tmu = iter()*iter().rho()*iter().nu();
174  volScalarField& mu = tmu.ref();
175 
176  for (++iter; iter != phases_.end(); ++iter)
177  {
178  mu += iter()*iter().rho()*iter().nu();
179  }
180 
181  return tmu;
182 }
183 
184 
186 Foam::multiphaseMixture::mu(const label patchi) const
187 {
188  PtrDictionary<phase>::const_iterator iter = phases_.begin();
189 
190  tmp<scalarField> tmu =
191  iter().boundaryField()[patchi]
192  *iter().rho().value()
193  *iter().nu(patchi);
194  scalarField& mu = tmu.ref();
195 
196  for (++iter; iter != phases_.end(); ++iter)
197  {
198  mu +=
199  iter().boundaryField()[patchi]
200  *iter().rho().value()
201  *iter().nu(patchi);
202  }
203 
204  return tmu;
205 }
206 
207 
210 {
211  PtrDictionary<phase>::const_iterator iter = phases_.begin();
212 
213  tmp<surfaceScalarField> tmuf =
214  fvc::interpolate(iter())*iter().rho()*fvc::interpolate(iter().nu());
215  surfaceScalarField& muf = tmuf.ref();
216 
217  for (++iter; iter != phases_.end(); ++iter)
218  {
219  muf +=
220  fvc::interpolate(iter())*iter().rho()*fvc::interpolate(iter().nu());
221  }
222 
223  return tmuf;
224 }
225 
226 
229 {
230  return nu_;
231 }
232 
233 
235 Foam::multiphaseMixture::nu(const label patchi) const
236 {
237  return nu_.boundaryField()[patchi];
238 }
239 
240 
243 {
244  return muf()/fvc::interpolate(rho());
245 }
246 
247 
250 {
251  tmp<surfaceScalarField> tstf
252  (
254  (
255  "surfaceTensionForce",
256  mesh_,
257  dimensionedScalar(dimensionSet(1, -2, -2, 0, 0), 0)
258  )
259  );
260 
261  surfaceScalarField& stf = tstf.ref();
262 
263  forAllConstIter(PtrDictionary<phase>, phases_, iter1)
264  {
265  const phase& alpha1 = iter1();
266 
268  ++iter2;
269 
270  for (; iter2 != phases_.end(); ++iter2)
271  {
272  const phase& alpha2 = iter2();
273 
274  sigmaTable::const_iterator sigma =
275  sigmas_.find(interfacePair(alpha1, alpha2));
276 
277  if (sigma == sigmas_.end())
278  {
280  << "Cannot find interface " << interfacePair(alpha1, alpha2)
281  << " in list of sigma values"
282  << exit(FatalError);
283  }
284 
285  stf += dimensionedScalar(dimSigma_, sigma())
286  *fvc::interpolate(K(alpha1, alpha2))*
287  (
288  fvc::interpolate(alpha2)*fvc::snGrad(alpha1)
289  - fvc::interpolate(alpha1)*fvc::snGrad(alpha2)
290  );
291  }
292  }
293 
294  return tstf;
295 }
296 
297 
299 {
300  correct();
301 
302  const Time& runTime = mesh_.time();
303 
304  volScalarField& alpha = phases_.first();
305 
306  const dictionary& alphaControls = mesh_.solverDict("alpha");
307  label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles"));
308  scalar cAlpha(alphaControls.lookup<scalar>("cAlpha"));
309 
310  if (nAlphaSubCycles > 1)
311  {
312  surfaceScalarField rhoPhiSum
313  (
314  IOobject
315  (
316  "rhoPhiSum",
317  runTime.timeName(),
318  mesh_
319  ),
320  mesh_,
321  dimensionedScalar(rhoPhi_.dimensions(), 0)
322  );
323 
324  dimensionedScalar totalDeltaT = runTime.deltaT();
325 
326  for
327  (
328  subCycle<volScalarField> alphaSubCycle(alpha, nAlphaSubCycles);
329  !(++alphaSubCycle).end();
330  )
331  {
332  solveAlphas(cAlpha);
333  rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi_;
334  }
335 
336  rhoPhi_ = rhoPhiSum;
337  }
338  else
339  {
340  solveAlphas(cAlpha);
341  }
342 
343  // Update the mixture kinematic viscosity
344  nu_ = mu()/rho();
345 }
346 
347 
349 {
350  forAllIter(PtrDictionary<phase>, phases_, iter)
351  {
352  iter().correct();
353  }
354 }
355 
356 
357 Foam::tmp<Foam::surfaceVectorField> Foam::multiphaseMixture::nHatfv
358 (
359  const volScalarField& alpha1,
360  const volScalarField& alpha2
361 ) const
362 {
363  /*
364  // Cell gradient of alpha
365  volVectorField gradAlpha =
366  alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2);
367 
368  // Interpolated face-gradient of alpha
369  surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
370  */
371 
372  surfaceVectorField gradAlphaf
373  (
375  - fvc::interpolate(alpha1)*fvc::interpolate(fvc::grad(alpha2))
376  );
377 
378  // Face unit interface normal
379  return gradAlphaf/(mag(gradAlphaf) + deltaN_);
380 }
381 
382 
383 Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseMixture::nHatf
384 (
385  const volScalarField& alpha1,
386  const volScalarField& alpha2
387 ) const
388 {
389  // Face unit interface normal flux
390  return nHatfv(alpha1, alpha2) & mesh_.Sf();
391 }
392 
393 
394 // Correction for the boundary condition on the unit normal nHat on
395 // walls to produce the correct contact angle.
396 
397 // The dynamic contact angle is calculated from the component of the
398 // velocity on the direction of the interface, parallel to the wall.
399 
400 void Foam::multiphaseMixture::correctContactAngle
401 (
402  const phase& alpha1,
403  const phase& alpha2,
404  surfaceVectorField::Boundary& nHatb
405 ) const
406 {
407  const volScalarField::Boundary& gbf
408  = alpha1.boundaryField();
409 
410  const fvBoundaryMesh& boundary = mesh_.boundary();
411 
412  forAll(boundary, patchi)
413  {
414  if (isA<alphaContactAngleFvPatchScalarField>(gbf[patchi]))
415  {
416  const alphaContactAngleFvPatchScalarField& acap =
417  refCast<const alphaContactAngleFvPatchScalarField>(gbf[patchi]);
418 
419  vectorField& nHatPatch = nHatb[patchi];
420 
421  vectorField AfHatPatch
422  (
423  mesh_.Sf().boundaryField()[patchi]
424  /mesh_.magSf().boundaryField()[patchi]
425  );
426 
427  alphaContactAngleFvPatchScalarField::thetaPropsTable::
428  const_iterator tp =
429  acap.thetaProps().find(interfacePair(alpha1, alpha2));
430 
431  if (tp == acap.thetaProps().end())
432  {
434  << "Cannot find interface " << interfacePair(alpha1, alpha2)
435  << "\n in table of theta properties for patch "
436  << acap.patch().name()
437  << exit(FatalError);
438  }
439 
440  bool matched = (tp.key().first() == alpha1.name());
441 
442  scalar theta0 = degToRad(tp().theta0(matched));
443  scalarField theta(boundary[patchi].size(), theta0);
444 
445  scalar uTheta = tp().uTheta();
446 
447  // Calculate the dynamic contact angle if required
448  if (uTheta > small)
449  {
450  scalar thetaA = degToRad(tp().thetaA(matched));
451  scalar thetaR = degToRad(tp().thetaR(matched));
452 
453  // Calculated the component of the velocity parallel to the wall
454  vectorField Uwall
455  (
456  U_.boundaryField()[patchi].patchInternalField()
457  - U_.boundaryField()[patchi]
458  );
459  Uwall -= (AfHatPatch & Uwall)*AfHatPatch;
460 
461  // Find the direction of the interface parallel to the wall
462  vectorField nWall
463  (
464  nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch
465  );
466 
467  // Normalise nWall
468  nWall /= (mag(nWall) + small);
469 
470  // Calculate Uwall resolved normal to the interface parallel to
471  // the interface
472  scalarField uwall(nWall & Uwall);
473 
474  theta += (thetaA - thetaR)*tanh(uwall/uTheta);
475  }
476 
477 
478  // Reset nHatPatch to correspond to the contact angle
479 
480  scalarField a12(nHatPatch & AfHatPatch);
481 
482  scalarField b1(cos(theta));
483 
484  scalarField b2(nHatPatch.size());
485 
486  forAll(b2, facei)
487  {
488  b2[facei] = cos(acos(a12[facei]) - theta[facei]);
489  }
490 
491  scalarField det(1.0 - a12*a12);
492 
493  scalarField a((b1 - a12*b2)/det);
494  scalarField b((b2 - a12*b1)/det);
495 
496  nHatPatch = a*AfHatPatch + b*nHatPatch;
497 
498  nHatPatch /= (mag(nHatPatch) + deltaN_.value());
499  }
500  }
501 }
502 
503 
504 Foam::tmp<Foam::volScalarField> Foam::multiphaseMixture::K
505 (
506  const phase& alpha1,
507  const phase& alpha2
508 ) const
509 {
510  tmp<surfaceVectorField> tnHatfv = nHatfv(alpha1, alpha2);
511 
512  correctContactAngle(alpha1, alpha2, tnHatfv.ref().boundaryFieldRef());
513 
514  // Simple expression for curvature
515  return -fvc::div(tnHatfv & mesh_.Sf());
516 }
517 
518 
521 {
522  tmp<volScalarField> tnearInt
523  (
525  (
526  "nearInterface",
527  mesh_,
529  )
530  );
531 
532  forAllConstIter(PtrDictionary<phase>, phases_, iter)
533  {
534  tnearInt.ref() =
535  max(tnearInt(), pos0(iter() - 0.01)*pos0(0.99 - iter()));
536  }
537 
538  return tnearInt;
539 }
540 
541 
542 void Foam::multiphaseMixture::solveAlphas
543 (
544  const scalar cAlpha
545 )
546 {
547  static label nSolves=-1;
548  nSolves++;
549 
550  word alphaScheme("div(phi,alpha)");
551  word alpharScheme("div(phirb,alpha)");
552 
553  surfaceScalarField phic(mag(phi_/mesh_.magSf()));
554  phic = min(cAlpha*phic, max(phic));
555 
556  PtrList<surfaceScalarField> alphaPhiCorrs(phases_.size());
557  int phasei = 0;
558 
559  forAllIter(PtrDictionary<phase>, phases_, iter)
560  {
561  phase& alpha = iter();
562 
563  alphaPhiCorrs.set
564  (
565  phasei,
567  (
568  "phi" + alpha.name() + "Corr",
569  fvc::flux
570  (
571  phi_,
572  alpha,
574  )
575  )
576  );
577 
578  surfaceScalarField& alphaPhiCorr = alphaPhiCorrs[phasei];
579 
580  forAllIter(PtrDictionary<phase>, phases_, iter2)
581  {
582  phase& alpha2 = iter2();
583 
584  if (&alpha2 == &alpha) continue;
585 
586  surfaceScalarField phir(phic*nHatf(alpha, alpha2));
587 
588  alphaPhiCorr += fvc::flux
589  (
590  -fvc::flux(-phir, alpha2, alpharScheme),
591  alpha,
593  );
594  }
595 
597  (
598  1.0/mesh_.time().deltaT().value(),
599  geometricOneField(),
600  alpha,
601  phi_,
602  alphaPhiCorr,
603  zeroField(),
604  zeroField(),
605  oneField(),
606  zeroField(),
607  true
608  );
609 
610  phasei++;
611  }
612 
613  MULES::limitSum(alphaPhiCorrs);
614 
615  rhoPhi_ = dimensionedScalar(dimensionSet(1, 0, -1, 0, 0), 0);
616 
617  volScalarField sumAlpha
618  (
619  IOobject
620  (
621  "sumAlpha",
622  mesh_.time().timeName(),
623  mesh_
624  ),
625  mesh_,
627  );
628 
629  phasei = 0;
630 
631  forAllIter(PtrDictionary<phase>, phases_, iter)
632  {
633  phase& alpha = iter();
634 
635  surfaceScalarField& alphaPhi = alphaPhiCorrs[phasei];
636  alphaPhi += upwind<scalar>(mesh_, phi_).flux(alpha);
637 
639  (
640  geometricOneField(),
641  alpha,
642  alphaPhi
643  );
644 
645  rhoPhi_ += alphaPhi*alpha.rho();
646 
647  Info<< alpha.name() << " volume fraction, min, max = "
648  << alpha.weightedAverage(mesh_.V()).value()
649  << ' ' << min(alpha).value()
650  << ' ' << max(alpha).value()
651  << endl;
652 
653  sumAlpha += alpha;
654 
655  phasei++;
656  }
657 
658  Info<< "Phase-sum volume fraction, min, max = "
659  << sumAlpha.weightedAverage(mesh_.V()).value()
660  << ' ' << min(sumAlpha).value()
661  << ' ' << max(sumAlpha).value()
662  << endl;
663 
664  // Correct the sum of the phase-fractions to avoid 'drift'
665  volScalarField sumCorr(1.0 - sumAlpha);
666  forAllIter(PtrDictionary<phase>, phases_, iter)
667  {
668  phase& alpha = iter();
669  alpha += alpha*sumCorr;
670  }
671 
672  calcAlphas();
673 }
674 
675 
677 {
678  if (transportModel::read())
679  {
680  bool readOK = true;
681 
682  PtrList<entry> phaseData(lookup("phases"));
683  label phasei = 0;
684 
685  forAllIter(PtrDictionary<phase>, phases_, iter)
686  {
687  readOK &= iter().read(phaseData[phasei++].dict());
688  }
689 
690  lookup("sigmas") >> sigmas_;
691 
692  return readOK;
693  }
694  else
695  {
696  return false;
697  }
698 }
699 
700 
701 // ************************************************************************* //
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin, const bool returnCorr)
dimensionedScalar tanh(const dimensionedScalar &ds)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
dimensionedScalar acos(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const surfaceVectorField & Sf() const
Return cell face area vectors.
tmp< surfaceScalarField > muf() const
Return the face-interpolated dynamic laminar viscosity.
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
label phasei
Definition: pEqn.H:27
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const Boundary & boundaryField() const
Return const-reference to the boundary field.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Unit conversion functions.
const dimensionedScalar & sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:177
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, const Internal &, const PtrList< fvsPatchField< scalar >> &)
Return a temporary field constructed from name,.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< volScalarField > trho
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
const dictionary & dict() const
Return populationBalanceCoeffs dictionary.
Calculate the snGrad of the given volField.
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:622
dimensionedScalar det(const dimensionedSphericalTensor &dt)
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
Definition: solution.C:342
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
T & first()
Return the first element of the list.
Definition: UListI.H:114
CGAL::Exact_predicates_exact_constructions_kernel K
multiphaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
tmp< volScalarField > rho() const
Return the mixture density.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
bool read()
Read base transportProperties dictionary.
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:30
const dimensionedScalar & b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
Calculate the gradient of the given field.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void correct()
Correct the mixture properties.
Calculate the face-flux of the given field.
const Type & value() const
Return const reference to value.
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
word timeName
Definition: getTimeIndex.H:3
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
faceListList boundary(nPatches)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
tmp< volScalarField > nu() const
Return the kinematic laminar viscosity.
Calculate the divergence of the given field.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
dimensionedScalar pos0(const dimensionedScalar &ds)
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
const label nAlphaSubCycles(alphaControls.lookup< label >("nAlphaSubCycles"))
void solve()
Solve for the mixture phase-fractions.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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.
Internal & ref()
Return a reference to the dimensioned internal field.
const word alphaScheme(mesh.divScheme(divAlphaName)[1].wordToken())
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
tmp< surfaceScalarField > nuf() const
Return the face-interpolated dynamic laminar viscosity.
IOdictionary(const IOobject &)
Construct given an IOobject.
Definition: IOdictionary.C:30
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionedScalar & mu
Atomic mass unit.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Time & time() const
Return time.
Definition: IOobject.C:328
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
word alpharScheme("div(phirb,alpha)")
messageStream Info
virtual bool read()=0
Read transportProperties dictionary.
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
MULES: Multidimensional universal limiter for explicit solution.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:322
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:53
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
friend class const_iterator
Definition: LPtrList.H:87
tmp< surfaceScalarField > surfaceTensionForce() const
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:540