phaseSystemSolve.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) 2013-2023 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 "phaseSystem.H"
27 
28 #include "MULES.H"
29 #include "subCycle.H"
30 
31 #include "fvcDdt.H"
32 #include "fvcDiv.H"
33 #include "fvcSnGrad.H"
34 #include "fvcFlux.H"
35 #include "fvcMeshPhi.H"
36 #include "fvcSup.H"
37 
38 #include "fvmDdt.H"
39 #include "fvmLaplacian.H"
40 #include "fvmSup.H"
41 
42 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
43 
45 {
46  const dictionary& alphaControls = mesh_.solution().solverDict("alpha");
47 
48  const label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles"));
49  const label nAlphaCorr(alphaControls.lookup<label>("nAlphaCorr"));
50 
51  const bool LTS = fv::localEulerDdt::enabled(mesh_);
52 
53  // Temporary switch for testing and comparing the standard split
54  // and the new un-split phase flux discretisation
55  const bool splitPhaseFlux
56  (
57  alphaControls.lookupOrDefault<Switch>("splitPhaseFlux", false)
58  );
59 
60  // Temporary switch for testing and comparing the standard mean flux
61  // and the new phase flux reference for the phase flux correction
62  const bool meanFluxReference
63  (
64  alphaControls.lookupOrDefault<Switch>("meanFluxReference", false)
65  );
66 
67  const scalar vDotResidualAlpha
68  (
69  alphaControls.lookupOrDefault("vDotResidualAlpha", 1e-4)
70  );
71 
72  // Optional reference phase which is not solved for
73  // but obtained from the sum of the other phases
74  phaseModel* referencePhasePtr = nullptr;
75 
76  // The phases which are solved
77  // i.e. the moving phases less the optional reference phase
78  phaseModelPartialList solvePhases;
79 
81  {
82  referencePhasePtr = &phases()[referencePhaseName_];
83 
84  solvePhases.setSize(movingPhases().size() - 1);
85  label solvePhasesi = 0;
86  forAll(movingPhases(), movingPhasei)
87  {
88  if (&movingPhases()[movingPhasei] != referencePhasePtr)
89  {
90  solvePhases.set(solvePhasesi++, &movingPhases()[movingPhasei]);
91  }
92  }
93  }
94  else
95  {
96  solvePhases = movingPhases();
97  }
98 
99  forAll(phases(), phasei)
100  {
101  phases()[phasei].correctBoundaryConditions();
102  }
103 
104  // Calculate the void fraction
105  volScalarField alphaVoid
106  (
107  IOobject
108  (
109  "alphaVoid",
110  mesh_.time().name(),
111  mesh_
112  ),
113  mesh_,
115  );
116  forAll(stationaryPhases(), stationaryPhasei)
117  {
118  alphaVoid -= stationaryPhases()[stationaryPhasei];
119  }
120 
121  // Calculate the effective flux of the moving phases
122  tmp<surfaceScalarField> tphiMoving(phi_);
123  if (stationaryPhases().size())
124  {
125  tphiMoving = phi_/upwind<scalar>(mesh_, phi_).interpolate(alphaVoid);
126  }
127  const surfaceScalarField& phiMoving = tphiMoving();
128 
129  bool dilatation = false;
130  forAll(movingPhases(), movingPhasei)
131  {
132  if (movingPhases()[movingPhasei].divU().valid())
133  {
134  dilatation = true;
135  break;
136  }
137  }
138 
139  for (int acorr=0; acorr<nAlphaCorr; acorr++)
140  {
143 
144  forAll(movingPhases(), movingPhasei)
145  {
146  const phaseModel& phase = movingPhases()[movingPhasei];
147  const volScalarField& alpha = phase;
148  const label phasei = phase.index();
149 
150  Sps.set
151  (
152  phasei,
154  (
155  IOobject
156  (
157  "Sp",
158  mesh_.time().name(),
159  mesh_
160  ),
161  mesh_,
163  )
164  );
165 
166  Sus.set
167  (
168  phasei,
170  (
171  "Su",
172  min(alpha.v(), scalar(1))
173  *fvc::div(fvc::absolute(phi_, phase.U()))->v()
174  )
175  );
176 
177  if (dilatation)
178  {
179  // Construct the dilatation rate source term
181  (
183  (
184  "vDot",
185  mesh_,
187  )
188  );
189 
190  forAll(phases(), phasej)
191  {
192  const phaseModel& phase2 = phases()[phasej];
193  const volScalarField& alpha2 = phase2;
194 
195  if (&phase2 != &phase)
196  {
197  if (!phase.stationary() && phase.divU().valid())
198  {
199  vDot += alpha2()*phase.divU()()();
200  }
201 
202  if (!phase2.stationary() && phase2.divU().valid())
203  {
204  vDot -= alpha()*phase2.divU()()();
205  }
206  }
207  }
208 
209  volScalarField::Internal& Sp = Sps[phasei];
210  volScalarField::Internal& Su = Sus[phasei];
211 
212  forAll(vDot, celli)
213  {
214  if (vDot[celli] > 0)
215  {
216  Sp[celli] -=
217  vDot[celli]
218  /max(1 - alpha[celli], vDotResidualAlpha);
219  Su[celli] +=
220  vDot[celli]
221  /max(1 - alpha[celli], vDotResidualAlpha);
222  }
223  else if (vDot[celli] < 0)
224  {
225  Sp[celli] +=
226  vDot[celli]
227  /max(alpha[celli], vDotResidualAlpha);
228  }
229  }
230  }
231  }
232 
233  tmp<volScalarField> trSubDeltaT;
234 
235  if (LTS && nAlphaSubCycles > 1)
236  {
237  trSubDeltaT =
239  }
240 
241  List<volScalarField*> alphaPtrs(phases().size());
242  forAll(phases(), phasei)
243  {
244  alphaPtrs[phasei] = &phases()[phasei];
245  }
246 
247  for
248  (
250  (
251  alphaPtrs,
253  );
254  !(++alphaSubCycle).end();
255  )
256  {
257  // Create correction fluxes
258  PtrList<surfaceScalarField> alphaPhis(phases().size());
259 
261  if (implicitPhasePressure() && (rAs.size()))
262  {
263  alphaDByAf = this->alphaDByAf(rAs);
264  }
265 
266  forAll(movingPhases(), movingPhasei)
267  {
268  const phaseModel& phase = movingPhases()[movingPhasei];
269  const volScalarField& alpha = phase;
270 
271  alphaPhis.set
272  (
273  phase.index(),
275  (
276  IOobject::groupName("alphaPhiCorr", phase.name()),
277  fvc::flux
278  (
279  splitPhaseFlux ? phi_ : phase.phi()(),
280  alpha,
281  "div(phi," + alpha.name() + ')'
282  )
283  )
284  );
285 
286  surfaceScalarField& alphaPhi = alphaPhis[phase.index()];
287 
288  if (splitPhaseFlux)
289  {
290  forAll(phases(), phasei)
291  {
292  const phaseModel& phase2 = phases()[phasei];
293  const volScalarField& alpha2 = phase2;
294 
295  if (&phase2 == &phase) continue;
296 
297  surfaceScalarField phir(phase.phi() - phase2.phi());
298 
300  (
301  cAlphas_.find(phaseInterface(phase, phase2))
302  );
303 
304  if (cAlpha != cAlphas_.end())
305  {
306  surfaceScalarField phic
307  (
308  (mag(phi_) + mag(phir))/mesh_.magSf()
309  );
310 
311  phir +=
312  min(cAlpha()*phic, max(phic))
313  *nHatf(alpha, alpha2);
314  }
315 
316  const word phirScheme
317  (
318  "div(phir,"
319  + alpha2.name() + ',' + alpha.name()
320  + ')'
321  );
322 
323  alphaPhi += fvc::flux
324  (
325  -fvc::flux(-phir, alpha2, phirScheme),
326  alpha,
327  phirScheme
328  );
329  }
330  }
331  else if (!cAlphas_.empty())
332  {
333  forAll(phases(), phasei)
334  {
335  const phaseModel& phase2 = phases()[phasei];
336  const volScalarField& alpha2 = phase2;
337 
338  if (&phase2 == &phase) continue;
339 
341  (
342  cAlphas_.find(phaseInterface(phase, phase2))
343  );
344 
345  if (cAlpha != cAlphas_.end())
346  {
347  const surfaceScalarField phir
348  (
349  phase.phi() - phase2.phi()
350  );
351 
352  const surfaceScalarField phic
353  (
354  (mag(phi_) + mag(phir))/mesh_.magSf()
355  );
356 
357  const surfaceScalarField phirc
358  (
359  min(cAlpha()*phic, max(phic))
360  *nHatf(alpha, alpha2)
361  );
362 
363  const word phirScheme
364  (
365  "div(phir,"
366  + alpha2.name() + ',' + alpha.name()
367  + ')'
368  );
369 
370  alphaPhi += fvc::flux
371  (
372  -fvc::flux(-phirc, alpha2, phirScheme),
373  alpha,
374  phirScheme
375  );
376  }
377  }
378  }
379 
380  if (alphaDByAf.valid())
381  {
382  alphaPhi +=
383  alphaDByAf()
384  *fvc::snGrad(alpha, "bounded")*mesh_.magSf();
385  }
386 
387  phase.correctInflowOutflow(alphaPhi);
388 
390  (
392  alpha,
393  meanFluxReference
394  ? phiMoving // Guarantees boundedness but less accurate
395  : phase.phi()(), // Less robust but more accurate
396  alphaPhi,
397  Sps[phase.index()],
398  Sus[phase.index()],
399  min(alphaVoid.primitiveField(), phase.alphaMax())(),
400  zeroField(),
401  false
402  );
403  }
404 
405  // Limit the flux corrections to ensure the phase fractions sum to 1
406  {
407  // Generate alphas for the moving phases
408  UPtrList<const volScalarField> alphasMoving
409  (
410  movingPhases().size()
411  );
412 
413  UPtrList<surfaceScalarField> alphaPhisMoving
414  (
415  movingPhases().size()
416  );
417 
418  forAll(movingPhases(), movingPhasei)
419  {
420  const phaseModel& phase = movingPhases()[movingPhasei];
421 
422  alphasMoving.set(movingPhasei, &phase);
423 
424  alphaPhisMoving.set
425  (
426  movingPhasei,
427  &alphaPhis[phase.index()]
428  );
429  }
430 
431  MULES::limitSum(alphasMoving, alphaPhisMoving, phiMoving);
432  }
433 
434  forAll(solvePhases, solvePhasei)
435  {
436  phaseModel& phase = solvePhases[solvePhasei];
437  volScalarField& alpha = phase;
438 
439  surfaceScalarField& alphaPhi = alphaPhis[phase.index()];
440  phase.correctInflowOutflow(alphaPhi);
441 
443  (
445  alpha,
446  alphaPhi,
447  Sps[phase.index()],
448  Sus[phase.index()]
449  );
450 
451  if (alphaSubCycle.index() == 1)
452  {
453  phase.alphaPhiRef() = alphaPhi;
454  }
455  else
456  {
457  phase.alphaPhiRef() += alphaPhi;
458  }
459  }
460 
461  if (alphaDByAf.valid())
462  {
463  // Update alphaDByAf due to changes in alpha
464  alphaDByAf = this->alphaDByAf(rAs);
465 
466  forAll(solvePhases, solvePhasei)
467  {
468  phaseModel& phase = solvePhases[solvePhasei];
469  volScalarField& alpha = phase;
470 
471  fvScalarMatrix alphaEqn
472  (
474  - fvm::laplacian(alphaDByAf(), alpha, "bounded")
475  );
476 
477  alphaEqn.solve();
478 
479  phase.alphaPhiRef() += alphaEqn.flux();
480  }
481  }
482 
483  // Report the phase fractions and the phase fraction sum
484  forAll(solvePhases, solvePhasei)
485  {
486  phaseModel& phase = solvePhases[solvePhasei];
487 
488  Info<< phase.name() << " fraction, min, max = "
489  << phase.weightedAverage(mesh_.V()).value()
490  << ' ' << min(phase).value()
491  << ' ' << max(phase).value()
492  << endl;
493  }
494 
495  if (referencePhasePtr)
496  {
497  volScalarField& referenceAlpha = *referencePhasePtr;
498  referenceAlpha = alphaVoid;
499 
500  forAll(solvePhases, solvePhasei)
501  {
502  referenceAlpha -= solvePhases[solvePhasei];
503  }
504  }
505  else
506  {
508  (
509  IOobject
510  (
511  "sumAlphaMoving",
512  mesh_.time().name(),
513  mesh_
514  ),
515  mesh_,
517  );
518  forAll(movingPhases(), movingPhasei)
519  {
520  sumAlphaMoving += movingPhases()[movingPhasei];
521  }
522 
523  Info<< "Phase-sum volume fraction, min, max = "
524  << (sumAlphaMoving + 1 - alphaVoid)()
525  .weightedAverage(mesh_.V()).value()
526  << ' ' << min(sumAlphaMoving + 1 - alphaVoid).value()
527  << ' ' << max(sumAlphaMoving + 1 - alphaVoid).value()
528  << endl;
529 
530  // Correct the sum of the phase fractions to avoid drift
531  forAll(movingPhases(), movingPhasei)
532  {
533  movingPhases()[movingPhasei] *= alphaVoid/sumAlphaMoving;
534  }
535  }
536  }
537 
538  if (nAlphaSubCycles > 1)
539  {
540  forAll(solvePhases, solvePhasei)
541  {
542  phaseModel& phase = solvePhases[solvePhasei];
543 
544  phase.alphaPhiRef() /= nAlphaSubCycles;
545  }
546  }
547 
548  forAll(solvePhases, solvePhasei)
549  {
550  phaseModel& phase = solvePhases[solvePhasei];
551 
552  phase.alphaRhoPhiRef() =
553  fvc::interpolate(phase.rho())*phase.alphaPhi();
554 
555  phase.maxMin(0, 1);
556  }
557 
558  if (referencePhasePtr)
559  {
560  phaseModel& referencePhase = *referencePhasePtr;
561 
562  referencePhase.alphaPhiRef() = phi_;
563 
564  forAll(solvePhases, solvePhasei)
565  {
566  phaseModel& phase = solvePhases[solvePhasei];
567  referencePhase.alphaPhiRef() -= phase.alphaPhi();
568  }
569 
570  referencePhase.alphaRhoPhiRef() =
571  fvc::interpolate(referencePhase.rho())
572  *referencePhase.alphaPhi();
573 
574  volScalarField& referenceAlpha = referencePhase;
575  referenceAlpha = alphaVoid;
576 
577  forAll(solvePhases, solvePhasei)
578  {
579  referenceAlpha -= solvePhases[solvePhasei];
580  }
581  }
582  }
583 }
584 
585 
586 // ************************************************************************* //
MULES: Multidimensional universal limiter for explicit solution.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const label nAlphaSubCycles(alphaControls.lookup< label >("nAlphaSubCycles"))
const dictionary & alphaControls
Definition: alphaControls.H:1
const label nAlphaCorr(alphaControls.lookup< label >("nAlphaCorr"))
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &, const Field< Type > &)
Return a temporary field constructed from name, mesh,.
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &) const
Calculate and return weighted average.
Generic GeometricField class.
void maxMin(const dimensioned< Type > &minDt, const dimensioned< Type > &maxDt)
const Internal::FieldType & primitiveField() const
Return a const-reference to the primitive field.
bool empty() const
Return true if the hash table is empty.
Definition: HashTableI.H:72
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:167
friend class const_iterator
Declare friendship with the const_iterator.
Definition: HashTable.H:197
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
static word groupName(Name name, const word &group)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
const iterator & end()
Definition: UILList.H:223
bool set(const label) const
Is element set.
Definition: UPtrListI.H:87
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void setSize(const label)
Reset size of UPtrList. This can only be used to set the size.
Definition: UPtrList.C:54
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const word & name() const
Return const reference to name.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
tmp< SurfaceField< Type > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:964
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const fvSolution & solution() const
Return the fvSolution.
Definition: fvMesh.C:1762
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
static tmp< volScalarField > localRSubDeltaT(const fvMesh &mesh, const label nAlphaSubCycles)
Calculate and return the reciprocal of the local sub-cycling.
Definition: localEulerDdt.C:70
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:37
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Class to represent an interface between phases. Derivations can further specify the configuration of ...
scalar alphaMax() const
Return the maximum phase-fraction (e.g. packing limit)
Definition: phaseModel.C:175
virtual bool stationary() const =0
Return whether the phase is stationary.
virtual const autoPtr< volScalarField > & divU() const =0
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual tmp< surfaceScalarField > phi() const =0
Return the volumetric flux.
label index() const
Return the index of the phase.
Definition: phaseModel.C:157
virtual surfaceScalarField & alphaRhoPhiRef()=0
Access the mass flux of the phase.
virtual tmp< volVectorField > U() const =0
Return the velocity.
virtual surfaceScalarField & alphaPhiRef()=0
Access the volumetric flux of the phase.
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Definition: phaseModel.C:245
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
virtual const volScalarField & rho() const =0
Return the density field.
virtual tmp< surfaceScalarField > alphaPhi() const =0
Return the volumetric flux of the phase.
virtual tmp< surfaceScalarField > alphaDByAf(const PtrList< volScalarField > &rAs) const =0
Return the phase diffusivity.
const fvMesh & mesh_
Reference to the mesh.
Definition: phaseSystem.H:128
const phaseModelPartialList & movingPhases() const
Return the models for phases that are moving.
Definition: phaseSystemI.H:116
virtual void solve(const PtrList< volScalarField > &rAs)
Solve for the phase fractions.
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:102
cAlphaTable cAlphas_
Interface compression coefficients.
Definition: phaseSystem.H:162
const phaseModelPartialList & stationaryPhases() const
Return the models for phases that are stationary.
Definition: phaseSystemI.H:130
surfaceScalarField phi_
Total volumetric flux.
Definition: phaseSystem.H:156
virtual bool implicitPhasePressure() const
Returns true if the phase pressure is treated implicitly.
Definition: phaseSystem.C:506
tmp< volScalarField > sumAlphaMoving() const
Return the sum of the phase fractions of the moving phases.
Definition: phaseSystem.C:54
tmp< surfaceScalarField > nHatf(const volScalarField &alpha1, const volScalarField &alpha2) const
Normal to interface between two phases dotted with face areas.
Definition: phaseSystem.C:151
word referencePhaseName_
Name of optional reference phase which is not solved for.
Definition: phaseSystem.H:138
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
Definition: solution.C:348
Perform a subCycleTime on a field or list of fields.
Definition: subCycle.H:235
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
Return the face-interpolate of the given cell field.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition: zeroField.H:53
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Calculate the snGrad of the given volField.
Calculate the field for explicit evaluation of implicit and explicit sources.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the laplacian of the field.
Calculate the matrix for implicit and explicit sources.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
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)
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:30
bool valid(const PtrList< ModelType > &l)
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > flux(const VolField< Type > &vf)
Return the face-flux field obtained from the given volVectorField.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
const doubleScalar e
Definition: doubleScalar.H:106
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
const dimensionSet dimless
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112