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 PtrList<volScalarField>& rAUs,
47  const PtrList<surfaceScalarField>& rAUfs
48 )
49 {
50  const dictionary& alphaControls = mesh_.solution().solverDict("alpha");
51 
52  const label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles"));
53  const label nAlphaCorr(alphaControls.lookup<label>("nAlphaCorr"));
54 
55  const bool LTS = fv::localEulerDdt::enabled(mesh_);
56 
57  // Temporary switch for testing and comparing the standard split
58  // and the new un-split phase flux discretisation
59  const bool splitPhaseFlux
60  (
61  alphaControls.lookupOrDefault<Switch>("splitPhaseFlux", false)
62  );
63 
64  // Temporary switch for testing and comparing the standard mean flux
65  // and the new phase flux reference for the phase flux correction
66  const bool meanFluxReference
67  (
68  alphaControls.lookupOrDefault<Switch>("meanFluxReference", false)
69  );
70 
71  // Optional reference phase which is not solved for
72  // but obtained from the sum of the other phases
73  phaseModel* referencePhasePtr = nullptr;
74 
75  // The phases which are solved
76  // i.e. the moving phases less the optional reference phase
77  phaseModelPartialList solvePhases;
78 
80  {
81  referencePhasePtr = &phases()[referencePhaseName_];
82 
83  solvePhases.setSize(movingPhases().size() - 1);
84  label solvePhasesi = 0;
85  forAll(movingPhases(), movingPhasei)
86  {
87  if (&movingPhases()[movingPhasei] != referencePhasePtr)
88  {
89  solvePhases.set(solvePhasesi++, &movingPhases()[movingPhasei]);
90  }
91  }
92  }
93  else
94  {
95  solvePhases = movingPhases();
96  }
97 
98  forAll(phases(), phasei)
99  {
100  phases()[phasei].correctBoundaryConditions();
101  }
102 
103  // Calculate the void fraction
104  volScalarField alphaVoid
105  (
106  IOobject
107  (
108  "alphaVoid",
109  mesh_.time().name(),
110  mesh_
111  ),
112  mesh_,
114  );
115  forAll(stationaryPhases(), stationaryPhasei)
116  {
117  alphaVoid -= stationaryPhases()[stationaryPhasei];
118  }
119 
120  // Calculate the effective flux of the moving phases
121  tmp<surfaceScalarField> tphiMoving(phi_);
122  if (stationaryPhases().size())
123  {
124  tphiMoving = phi_/upwind<scalar>(mesh_, phi_).interpolate(alphaVoid);
125  }
126  const surfaceScalarField& phiMoving = tphiMoving();
127 
128  bool dilatation = false;
129  forAll(movingPhases(), movingPhasei)
130  {
131  if (movingPhases()[movingPhasei].divU().valid())
132  {
133  dilatation = true;
134  break;
135  }
136  }
137 
138  for (int acorr=0; acorr<nAlphaCorr; acorr++)
139  {
142 
143  forAll(movingPhases(), movingPhasei)
144  {
145  const phaseModel& phase = movingPhases()[movingPhasei];
146  const volScalarField& alpha = phase;
147  const label phasei = phase.index();
148 
149  Sps.set
150  (
151  phasei,
153  (
154  IOobject
155  (
156  "Sp",
157  mesh_.time().name(),
158  mesh_
159  ),
160  mesh_,
162  )
163  );
164 
165  Sus.set
166  (
167  phasei,
169  (
170  "Su",
171  min(alpha.v(), scalar(1))
172  *fvc::div(fvc::absolute(phi_, phase.U()))->v()
173  )
174  );
175 
176  if (dilatation)
177  {
178  // Construct the dilatation rate source term
180  (
182  (
183  "dgdt",
184  mesh_,
186  )
187  );
188 
189  forAll(phases(), phasej)
190  {
191  const phaseModel& phase2 = phases()[phasej];
192  const volScalarField& alpha2 = phase2;
193 
194  if (&phase2 != &phase)
195  {
196  if (phase.divU().valid())
197  {
198  dgdt += alpha2()*phase.divU()()();
199  }
200 
201  if (phase2.divU().valid())
202  {
203  dgdt -= alpha()*phase2.divU()()();
204  }
205  }
206  }
207 
208  volScalarField::Internal& Sp = Sps[phasei];
209  volScalarField::Internal& Su = Sus[phasei];
210 
211  forAll(dgdt, celli)
212  {
213  if (dgdt[celli] > 0)
214  {
215  Sp[celli] -= dgdt[celli]/max(1 - alpha[celli], 1e-4);
216  Su[celli] += dgdt[celli]/max(1 - alpha[celli], 1e-4);
217  }
218  else if (dgdt[celli] < 0)
219  {
220  Sp[celli] += dgdt[celli]/max(alpha[celli], 1e-4);
221  }
222  }
223  }
224  }
225 
226  tmp<volScalarField> trSubDeltaT;
227 
228  if (LTS && nAlphaSubCycles > 1)
229  {
230  trSubDeltaT =
232  }
233 
234  List<volScalarField*> alphaPtrs(phases().size());
235  forAll(phases(), phasei)
236  {
237  alphaPtrs[phasei] = &phases()[phasei];
238  }
239 
240  for
241  (
243  (
244  alphaPtrs,
246  );
247  !(++alphaSubCycle).end();
248  )
249  {
250  // Create correction fluxes
251  PtrList<surfaceScalarField> alphaPhis(phases().size());
252 
254  if (implicitPhasePressure() && (rAUs.size() || rAUfs.size()))
255  {
256  alphaDByAf = this->alphaDByAf(rAUs, rAUfs);
257  }
258 
259  forAll(movingPhases(), movingPhasei)
260  {
261  const phaseModel& phase = movingPhases()[movingPhasei];
262  const volScalarField& alpha = phase;
263 
264  alphaPhis.set
265  (
266  phase.index(),
268  (
269  IOobject::groupName("alphaPhiCorr", phase.name()),
270  fvc::flux
271  (
272  splitPhaseFlux ? phi_ : phase.phi()(),
273  alpha,
274  "div(phi," + alpha.name() + ')'
275  )
276  )
277  );
278 
279  surfaceScalarField& alphaPhi = alphaPhis[phase.index()];
280 
281  if (splitPhaseFlux)
282  {
283  forAll(phases(), phasei)
284  {
285  const phaseModel& phase2 = phases()[phasei];
286  const volScalarField& alpha2 = phase2;
287 
288  if (&phase2 == &phase) continue;
289 
290  surfaceScalarField phir(phase.phi() - phase2.phi());
291 
293  (
294  cAlphas_.find(phaseInterface(phase, phase2))
295  );
296 
297  if (cAlpha != cAlphas_.end())
298  {
299  surfaceScalarField phic
300  (
301  (mag(phi_) + mag(phir))/mesh_.magSf()
302  );
303 
304  phir +=
305  min(cAlpha()*phic, max(phic))
306  *nHatf(alpha, alpha2);
307  }
308 
309  const word phirScheme
310  (
311  "div(phir,"
312  + alpha2.name() + ',' + alpha.name()
313  + ')'
314  );
315 
316  alphaPhi += fvc::flux
317  (
318  -fvc::flux(-phir, alpha2, phirScheme),
319  alpha,
320  phirScheme
321  );
322  }
323  }
324  else if (!cAlphas_.empty())
325  {
326  forAll(phases(), phasei)
327  {
328  const phaseModel& phase2 = phases()[phasei];
329  const volScalarField& alpha2 = phase2;
330 
331  if (&phase2 == &phase) continue;
332 
334  (
335  cAlphas_.find(phaseInterface(phase, phase2))
336  );
337 
338  if (cAlpha != cAlphas_.end())
339  {
340  const surfaceScalarField phir
341  (
342  phase.phi() - phase2.phi()
343  );
344 
345  const surfaceScalarField phic
346  (
347  (mag(phi_) + mag(phir))/mesh_.magSf()
348  );
349 
350  const surfaceScalarField phirc
351  (
352  min(cAlpha()*phic, max(phic))
353  *nHatf(alpha, alpha2)
354  );
355 
356  const word phirScheme
357  (
358  "div(phir,"
359  + alpha2.name() + ',' + alpha.name()
360  + ')'
361  );
362 
363  alphaPhi += fvc::flux
364  (
365  -fvc::flux(-phirc, alpha2, phirScheme),
366  alpha,
367  phirScheme
368  );
369  }
370  }
371  }
372 
373  if (alphaDByAf.valid())
374  {
375  alphaPhi +=
376  alphaDByAf()
377  *fvc::snGrad(alpha, "bounded")*mesh_.magSf();
378  }
379 
380  phase.correctInflowOutflow(alphaPhi);
381 
383  (
385  alpha,
386  meanFluxReference
387  ? phiMoving // Guarantees boundedness but less accurate
388  : phase.phi()(), // Less robust but more accurate
389  alphaPhi,
390  Sps[phase.index()],
391  Sus[phase.index()],
392  min(alphaVoid.primitiveField(), phase.alphaMax())(),
393  zeroField(),
394  false
395  );
396  }
397 
398  // Limit the flux corrections to ensure the phase fractions sum to 1
399  {
400  // Generate alphas for the moving phases
401  UPtrList<const volScalarField> alphasMoving
402  (
403  movingPhases().size()
404  );
405 
406  UPtrList<surfaceScalarField> alphaPhisMoving
407  (
408  movingPhases().size()
409  );
410 
411  forAll(movingPhases(), movingPhasei)
412  {
413  const phaseModel& phase = movingPhases()[movingPhasei];
414 
415  alphasMoving.set(movingPhasei, &phase);
416 
417  alphaPhisMoving.set
418  (
419  movingPhasei,
420  &alphaPhis[phase.index()]
421  );
422  }
423 
424  MULES::limitSum(alphasMoving, alphaPhisMoving, phiMoving);
425  }
426 
427  forAll(solvePhases, solvePhasei)
428  {
429  phaseModel& phase = solvePhases[solvePhasei];
430  volScalarField& alpha = phase;
431 
432  surfaceScalarField& alphaPhi = alphaPhis[phase.index()];
433  phase.correctInflowOutflow(alphaPhi);
434 
436  (
438  alpha,
439  alphaPhi,
440  Sps[phase.index()],
441  Sus[phase.index()]
442  );
443 
444  if (alphaSubCycle.index() == 1)
445  {
446  phase.alphaPhiRef() = alphaPhi;
447  }
448  else
449  {
450  phase.alphaPhiRef() += alphaPhi;
451  }
452  }
453 
454  if (alphaDByAf.valid())
455  {
456  // Update alphaDByAf due to changes in alpha
457  alphaDByAf = this->alphaDByAf(rAUs, rAUfs);
458 
459  forAll(solvePhases, solvePhasei)
460  {
461  phaseModel& phase = solvePhases[solvePhasei];
462  volScalarField& alpha = phase;
463 
464  fvScalarMatrix alphaEqn
465  (
467  - fvm::laplacian(alphaDByAf(), alpha, "bounded")
468  );
469 
470  alphaEqn.solve();
471 
472  phase.alphaPhiRef() += alphaEqn.flux();
473  }
474  }
475 
476  // Report the phase fractions and the phase fraction sum
477  forAll(solvePhases, solvePhasei)
478  {
479  phaseModel& phase = solvePhases[solvePhasei];
480 
481  Info<< phase.name() << " fraction, min, max = "
482  << phase.weightedAverage(mesh_.V()).value()
483  << ' ' << min(phase).value()
484  << ' ' << max(phase).value()
485  << endl;
486  }
487 
488  if (referencePhasePtr)
489  {
490  volScalarField& referenceAlpha = *referencePhasePtr;
491  referenceAlpha = alphaVoid;
492 
493  forAll(solvePhases, solvePhasei)
494  {
495  referenceAlpha -= solvePhases[solvePhasei];
496  }
497  }
498  else
499  {
501  (
502  IOobject
503  (
504  "sumAlphaMoving",
505  mesh_.time().name(),
506  mesh_
507  ),
508  mesh_,
510  );
511  forAll(movingPhases(), movingPhasei)
512  {
513  sumAlphaMoving += movingPhases()[movingPhasei];
514  }
515 
516  Info<< "Phase-sum volume fraction, min, max = "
517  << (sumAlphaMoving + 1 - alphaVoid)()
518  .weightedAverage(mesh_.V()).value()
519  << ' ' << min(sumAlphaMoving + 1 - alphaVoid).value()
520  << ' ' << max(sumAlphaMoving + 1 - alphaVoid).value()
521  << endl;
522 
523  // Correct the sum of the phase fractions to avoid drift
524  forAll(movingPhases(), movingPhasei)
525  {
526  movingPhases()[movingPhasei] *= alphaVoid/sumAlphaMoving;
527  }
528  }
529  }
530 
531  if (nAlphaSubCycles > 1)
532  {
533  forAll(solvePhases, solvePhasei)
534  {
535  phaseModel& phase = solvePhases[solvePhasei];
536 
537  phase.alphaPhiRef() /= nAlphaSubCycles;
538  }
539  }
540 
541  forAll(solvePhases, solvePhasei)
542  {
543  phaseModel& phase = solvePhases[solvePhasei];
544 
545  phase.alphaRhoPhiRef() =
546  fvc::interpolate(phase.rho())*phase.alphaPhi();
547 
548  phase.maxMin(0, 1);
549  }
550 
551  if (referencePhasePtr)
552  {
553  phaseModel& referencePhase = *referencePhasePtr;
554 
555  referencePhase.alphaPhiRef() = phi_;
556 
557  forAll(solvePhases, solvePhasei)
558  {
559  phaseModel& phase = solvePhases[solvePhasei];
560  referencePhase.alphaPhiRef() -= phase.alphaPhi();
561  }
562 
563  referencePhase.alphaRhoPhiRef() =
564  fvc::interpolate(referencePhase.rho())
565  *referencePhase.alphaPhi();
566 
567  volScalarField& referenceAlpha = referencePhase;
568  referenceAlpha = alphaVoid;
569 
570  forAll(solvePhases, solvePhasei)
571  {
572  referenceAlpha -= solvePhases[solvePhasei];
573  }
574  }
575  }
576 }
577 
578 
579 // ************************************************************************* //
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 internal 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:142
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:65
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:78
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:160
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:963
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:406
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const fvSolution & solution() const
Return the fvSchemes.
Definition: fvMesh.C:1759
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:139
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:121
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:209
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:109
virtual const volScalarField & rho() const =0
Return the density field.
virtual tmp< surfaceScalarField > alphaPhi() const =0
Return the volumetric flux of the phase.
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:55
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:41
virtual void solve(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs)
Solve for the phase fractions.
cAlphaTable cAlphas_
Interface compression coefficients.
Definition: phaseSystem.H:162
const phaseModelPartialList & stationaryPhases() const
Return the models for phases that are stationary.
Definition: phaseSystemI.H:69
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:513
tmp< volScalarField > sumAlphaMoving() const
Return the sum of the phase fractions of the moving phases.
Definition: phaseSystem.C:78
tmp< surfaceScalarField > nHatf(const volScalarField &alpha1, const volScalarField &alpha2) const
Normal to interface between two phases dotted with face areas.
Definition: phaseSystem.C:175
virtual tmp< surfaceScalarField > alphaDByAf(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs) const =0
Return the phase diffusivity.
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
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const doubleScalar e
Definition: doubleScalar.H:105
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:251
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)
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112