thermophysicalPredictor.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) 2022-2024 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 "XiFluid.H"
27 #include "fvcSnGrad.H"
28 #include "fvcLaplacian.H"
29 #include "fvcDdt.H"
30 #include "fvcMeshPhi.H"
31 #include "fvmDiv.H"
32 #include "fvmSup.H"
33 
34 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
35 
37 (
39 )
40 {
41  volScalarField& ft = thermo_.Y("ft");
42 
43  fvScalarMatrix ftEqn
44  (
45  fvm::ddt(rho, ft)
46  + mvConvection.fvmDiv(phi, ft)
48  ==
49  fvModels().source(rho, ft)
50  );
51 
52  fvConstraints().constrain(ftEqn);
53 
54  ftEqn.solve();
55 
57 }
58 
59 
61 (
62  const volScalarField& c,
63  const surfaceScalarField& nf,
64  const dimensionedScalar& dMgb
65 ) const
66 {
67  dimensionedScalar StCorr("StCorr", dimless, 1.0);
68 
69  if (ign.igniting())
70  {
71  // Calculate volume of ignition kernel
72  const dimensionedScalar Vk
73  (
74  "Vk",
75  dimVolume,
76  gSum(c*mesh.V().primitiveField())
77  );
78  dimensionedScalar Ak("Ak", dimArea, 0.0);
79 
80  if (Vk.value() > small)
81  {
82  // Calculate kernel area from its volume
83  // and the dimensionality of the case
84 
85  switch(mesh.nGeometricD())
86  {
87  case 3:
88  {
89  // Assume it is part-spherical
90  const scalar sphereFraction
91  (
92  combustionProperties.lookup<scalar>
93  (
94  "ignitionSphereFraction"
95  )
96  );
97 
98  Ak = sphereFraction*4.0*constant::mathematical::pi
99  *pow
100  (
101  3.0*Vk
102  /(sphereFraction*4.0*constant::mathematical::pi),
103  2.0/3.0
104  );
105  }
106  break;
107 
108  case 2:
109  {
110  // Assume it is part-circular
111  const dimensionedScalar thickness
112  (
113  combustionProperties.lookup("ignitionThickness")
114  );
115 
116  const scalar circleFraction
117  (
118  combustionProperties.lookup<scalar>
119  (
120  "ignitionCircleFraction"
121  )
122  );
123 
124  Ak = circleFraction*constant::mathematical::pi*thickness
125  *sqrt
126  (
127  4.0*Vk
128  /(
129  circleFraction
130  *thickness
132  )
133  );
134  }
135  break;
136 
137  case 1:
138  // Assume it is plane or two planes
139  Ak = dimensionedScalar
140  (
141  combustionProperties.lookup("ignitionKernelArea")
142  );
143  break;
144  }
145 
146  // Calculate kernel area from b field consistent with the
147  // discretisation of the b equation.
148  const volScalarField mgb
149  (
150  fvc::div(nf, b, "div(phiSt,b)") - b*fvc::div(nf) + dMgb
151  );
152  const dimensionedScalar AkEst = gSum(mgb*mesh.V().primitiveField());
153 
154  StCorr.value() = max(min((Ak/AkEst).value(), 10.0), 1.0);
155 
156  Info<< "StCorr = " << StCorr.value() << endl;
157  }
158  }
159 
160  return StCorr;
161 }
162 
163 
165 (
167 )
168 {
169  volScalarField& b(b_);
170  volScalarField& Xi(Xi_);
171 
172  // progress variable
173  // ~~~~~~~~~~~~~~~~~
174  const volScalarField c("c", scalar(1) - b);
175 
176  // Unburnt gas density
177  // ~~~~~~~~~~~~~~~~~~~
178  const volScalarField rhou(thermo.rhou());
179 
180  // Calculate flame normal etc.
181  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
182 
183  volVectorField n("n", fvc::grad(b));
184 
185  volScalarField mgb(mag(n));
186 
187  const dimensionedScalar dMgb = 1.0e-3*
188  (b*c*mgb)().weightedAverage(mesh.V())
189  /((b*c)().weightedAverage(mesh.V()) + small)
190  + dimensionedScalar(mgb.dimensions(), small);
191 
192  mgb += dMgb;
193 
194  const surfaceVectorField SfHat(mesh.Sf()/mesh.magSf());
196  nfVec += SfHat*(fvc::snGrad(b) - (SfHat & nfVec));
197  nfVec /= (mag(nfVec) + dMgb);
198  surfaceScalarField nf((mesh.Sf() & nfVec));
199  n /= mgb;
200 
201  // Calculate turbulent flame speed flux
202  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
203  const surfaceScalarField phiSt
204  (
205  "phiSt",
206  fvc::interpolate(rhou*StCorr(c, nf, dMgb)*Su*Xi)*nf
207  );
208 
209  const scalar StCoNum = max
210  (
211  mesh.surfaceInterpolation::deltaCoeffs()
212  *mag(phiSt)/(fvc::interpolate(rho)*mesh.magSf())
213  ).value()*runTime.deltaTValue();
214 
215  Info<< "Max St-Courant Number = " << StCoNum << endl;
216 
217  // Create b equation
218  // ~~~~~~~~~~~~~~~~~
219  fvScalarMatrix bEqn
220  (
221  fvm::ddt(rho, b)
222  + mvConvection.fvmDiv(phi, b)
223  + fvm::div(phiSt, b)
224  - fvm::Sp(fvc::div(phiSt), b)
226  ==
227  fvModels().source(rho, b)
228  );
229 
230 
231  // Add ignition cell contribution to b-equation
232  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
233  forAll(ign.sites(), i)
234  {
235  const ignitionSite& ignSite = ign.sites()[i];
236 
237  if (ignSite.igniting())
238  {
239  forAll(ignSite.cells(), icelli)
240  {
241  label ignCell = ignSite.cells()[icelli];
242  Info<< "Igniting cell " << ignCell;
243 
244  Info<< " state :"
245  << ' ' << b[ignCell]
246  << ' ' << Xi[ignCell]
247  << ' ' << Su[ignCell]
248  << ' ' << mgb[ignCell]
249  << endl;
250 
251  bEqn.diag()[ignSite.cells()[icelli]] +=
252  (
253  ignSite.strength()*ignSite.cellVolumes()[icelli]
254  *rhou[ignSite.cells()[icelli]]/ignSite.duration()
255  )/(b[ignSite.cells()[icelli]] + 0.001);
256  }
257  }
258  }
259 
260 
261  // Solve for b
262  // ~~~~~~~~~~~
263  bEqn.relax();
264 
265  fvConstraints().constrain(bEqn);
266 
267  bEqn.solve();
268 
270 
271  Info<< "min(b) = " << min(b).value() << endl;
272 
273 
274  // Calculate coefficients for Gulder's flame speed correlation
275  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
276 
277  const volScalarField up(uPrimeCoef*sqrt((2.0/3.0)*momentumTransport->k()));
278  // volScalarField up(sqrt(mag(diag(n * n) & diag(momentumTransport->r()))));
279 
280  const volScalarField epsilon
281  (
282  pow(uPrimeCoef, 3)*momentumTransport->epsilon()
283  );
284 
285  const volScalarField tauEta(sqrt(thermo.muu()/(rhou*epsilon)));
286 
287  const volScalarField Reta
288  (
289  up
290  / (
291  sqrt(epsilon*tauEta)
292  + dimensionedScalar(up.dimensions(), 1e-8)
293  )
294  );
295 
296  // volScalarField l = 0.337*k*sqrt(k)/epsilon;
297  // Reta *= max((l - dimensionedScalar(dimLength, 1.5e-3))/l, 0);
298 
299  // Calculate Xi flux
300  // ~~~~~~~~~~~~~~~~~
301  const surfaceScalarField phiXi
302  (
303  phiSt
305  (
307  )*nf
308  + fvc::interpolate(rho)*fvc::interpolate(Su*(1.0/Xi - Xi))*nf
309  );
310 
311 
312  // Calculate mean and turbulent strain rates
313  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
314 
315  const volVectorField Ut(U + Su*Xi*n);
316  const volScalarField sigmat((n & n)*fvc::div(Ut) - (n & fvc::grad(Ut) & n));
317 
318  const volScalarField sigmas
319  (
320  ((n & n)*fvc::div(U) - (n & fvc::grad(U) & n))/Xi
321  + (
322  (n & n)*fvc::div(Su*n)
323  - (n & fvc::grad(Su*n) & n)
324  )*(Xi + scalar(1))/(2*Xi)
325  );
326 
327 
328  // Calculate the unstrained laminar flame speed
329  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
331 
332 
333  // Calculate the laminar flame speed in equilibrium with the applied strain
334  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
335  const volScalarField SuInf
336  (
337  Su0*max(scalar(1) - sigmas/sigmaExt, scalar(0.01))
338  );
339 
340  if (SuModel == "unstrained")
341  {
342  Su == Su0;
343  }
344  else if (SuModel == "equilibrium")
345  {
346  Su == SuInf;
347  }
348  else if (SuModel == "transport")
349  {
350  // Solve for the strained laminar flame speed
351  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
352 
353  const volScalarField Rc
354  (
355  (sigmas*SuInf*(Su0 - SuInf) + sqr(SuMin)*sigmaExt)
356  /(sqr(Su0 - SuInf) + sqr(SuMin))
357  );
358 
359  fvScalarMatrix SuEqn
360  (
361  fvm::ddt(rho, Su)
362  + fvm::div(phi + phiXi, Su, "div(phiXi,Su)")
363  - fvm::Sp(fvc::div(phiXi), Su)
364  ==
365  - fvm::SuSp(-rho*Rc*Su0/Su, Su)
366  - fvm::SuSp(rho*(sigmas + Rc), Su)
367  + fvModels().source(rho, Su)
368  );
369 
370  SuEqn.relax();
371 
372  fvConstraints().constrain(SuEqn);
373 
374  SuEqn.solve();
375 
377 
378  // Limit the maximum Su
379  // ~~~~~~~~~~~~~~~~~~~~
380  Su.min(SuMax);
381  Su.max(SuMin);
382  }
383  else
384  {
386  << "Unknown Su model " << SuModel
387  << abort(FatalError);
388  }
389 
390 
391  // Calculate Xi according to the selected flame wrinkling model
392  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
393 
394  if (XiModel == "fixed")
395  {
396  // Do nothing, Xi is fixed!
397  }
398  else if (XiModel == "algebraic")
399  {
400  // Simple algebraic model for Xi based on Gulders correlation
401  // with a linear correction function to give a plausible profile for Xi
402  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
403  Xi == scalar(1) +
404  (scalar(1) + (2*XiShapeCoef)*(scalar(0.5) - b))
405  *XiCoef*sqrt(up/(Su + SuMin))*Reta;
406  }
407  else if (XiModel == "transport")
408  {
409  // Calculate Xi transport coefficients based on Gulders correlation
410  // and DNS data for the rate of generation
411  // with a linear correction function to give a plausible profile for Xi
412  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
413 
414  const volScalarField XiEqStar
415  (
416  scalar(1.001) + XiCoef*sqrt(up/(Su + SuMin))*Reta
417  );
418 
419  const volScalarField XiEq
420  (
421  scalar(1.001)
422  + (
423  scalar(1)
424  + (2*XiShapeCoef)
425  *(scalar(0.5) - min(max(b, scalar(0)), scalar(1)))
426  )*(XiEqStar - scalar(1.001))
427  );
428 
429  const volScalarField Gstar(0.28/tauEta);
430  const volScalarField R(Gstar*XiEqStar/(XiEqStar - scalar(1)));
431  const volScalarField G(R*(XiEq - scalar(1.001))/XiEq);
432 
433  // R *= (Gstar + 2*mag(dev(symm(fvc::grad(U)))))/Gstar;
434 
435  // Solve for the flame wrinkling
436  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
437  fvScalarMatrix XiEqn
438  (
439  fvm::ddt(rho, Xi)
440  + fvm::div(phi + phiXi, Xi, "div(phiXi,Xi)")
441  - fvm::Sp(fvc::div(phiXi), Xi)
442  ==
443  rho*R
444  - fvm::Sp(rho*(R - G), Xi)
445  - fvm::Sp
446  (
447  rho*max
448  (
449  sigmat - sigmas,
450  dimensionedScalar(sigmat.dimensions(), 0)
451  ),
452  Xi
453  )
454  + fvModels().source(rho, Xi)
455  );
456 
457  XiEqn.relax();
458 
459  fvConstraints().constrain(XiEqn);
460 
461  XiEqn.solve();
462 
463  fvConstraints().constrain(Xi);
464 
465  // Correct boundedness of Xi
466  // ~~~~~~~~~~~~~~~~~~~~~~~~~
467  Xi.max(1.0);
468  Info<< "max(Xi) = " << max(Xi).value() << endl;
469  Info<< "max(XiEq) = " << max(XiEq).value() << endl;
470  }
471  else
472  {
474  << "Unknown Xi model " << XiModel
475  << abort(FatalError);
476  }
477 
478  Info<< "Combustion progress = "
479  << 100*(scalar(1) - b)().weightedAverage(mesh.V()).value() << "%"
480  << endl;
481 
482  St = Xi*Su;
483 }
484 
485 
487 (
489 )
490 {
491  volScalarField& heau = thermo_.heu();
492 
493  const volScalarField::Internal rhoByRhou(rho()/thermo.rhou()());
494 
495  fvScalarMatrix heauEqn
496  (
497  fvm::ddt(rho, heau) + mvConvection.fvmDiv(phi, heau)
498  + rhoByRhou
499  *(
500  (fvc::ddt(rho, K) + fvc::div(phi, K))()
501  + pressureWork
502  (
503  heau.name() == "eau"
504  ? mvConvection.fvcDiv(phi, p/rho)()
505  : -dpdt
506  )
507  )
508  + thermophysicalTransport.divq(heau)
509 
510  // These terms cannot be used in partially-premixed combustion due to
511  // the resultant inconsistency between ft and heau transport.
512  // A possible solution would be to solve for ftu as well as ft.
513  //- fvm::div(muEff*fvc::grad(b)/(b + 0.001), heau)
514  //+ fvm::Sp(fvc::div(muEff*fvc::grad(b)/(b + 0.001)), heau)
515 
516  ==
517  fvModels().source(rho, heau)
518  );
519 
520  fvConstraints().constrain(heauEqn);
521 
522  heauEqn.solve();
523 
524  fvConstraints().constrain(heau);
525 }
526 
527 
529 (
531 )
532 {
533  volScalarField& hea = thermo_.he();
534 
536  (
537  fvm::ddt(rho, hea) + mvConvection.fvmDiv(phi, hea)
538  + fvc::ddt(rho, K) + fvc::div(phi, K)
539  + pressureWork
540  (
541  hea.name() == "ea"
542  ? mvConvection.fvcDiv(phi, p/rho)()
543  : -dpdt
544  )
545  + thermophysicalTransport.divq(hea)
546  ==
547  (
548  buoyancy.valid()
549  ? fvModels().source(rho, hea) + rho*(U & buoyancy->g)
550  : fvModels().source(rho, hea)
551  )
552  );
553 
554  EaEqn.relax();
555 
557 
558  EaEqn.solve();
559 
560  fvConstraints().constrain(hea);
561 }
562 
563 
565 {
567  (
569  (
570  mesh,
571  fields,
572  phi,
573  mesh.schemes().div("div(phi,ft_b_ha_hau)")
574  )
575  );
576 
577  if (thermo_.containsSpecie("ft"))
578  {
579  ftSolve(mvConvection());
580  }
581 
582  if (ign.ignited())
583  {
584  bSolve(mvConvection());
585  EauSolve(mvConvection());
586  }
587 
588  EaSolve(mvConvection());
589 
590  if (!ign.ignited())
591  {
592  thermo_.heu() == thermo.he();
593  }
594 
595  thermo_.correct();
596 }
597 
598 
599 // ************************************************************************* //
fvScalarMatrix EaEqn(betav *fvm::ddt(rho, hea)+mvConvection->fvmDiv(phi, hea)+betav *fvc::ddt(rho, K)+fvc::div(phi, K)+(hea.name()=="ea" ? fvc::div(fvc::absolute(phi, rho, U), p/rho) :-betav *dpdt) - fvm::laplacian(Db, hea)+betav *fvModels.source(rho, hea))
dimensionedScalar StCorr("StCorr", dimless, 1.0)
StCoNum
Definition: StCourantNo.H:36
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
tmp< fv::convectionScheme< scalar > > mvConvection(fv::convectionScheme< scalar >::New(mesh, fields, phi, mesh.schemes().div("div(phi,ft_b_ha_hau)")))
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
void max(const dimensioned< Type > &)
const Internal::FieldType & primitiveField() const
Return a const-reference to the primitive field.
const word & name() const
Return name.
Definition: IOobject.H:310
Base-class for all Xi models used by the b-Xi combustion model. See Technical Report SH/RE/01R for de...
Definition: XiModel.H:109
const Type & value() const
Return const reference to value.
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:603
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
Finite volume models.
Definition: fvModels.H:65
Abstract base class for convection schemes.
Foam::ignitionSite.
Definition: ignitionSite.H:58
const labelList & cells() const
Return the ignition cells updated if the mesh moved.
Definition: ignitionSite.C:86
scalar duration() const
Definition: ignitionSite.H:145
const scalarList & cellVolumes() const
Definition: ignitionSite.H:158
scalar strength() const
Definition: ignitionSite.H:150
bool igniting() const
Definition: ignitionSite.C:98
scalarField & diag()
Definition: lduMatrix.C:186
virtual PtrList< volScalarField > & Y()=0
Access the mass-fraction fields.
Foam::fvConstraints & fvConstraints() const
Return the fvConstraints that are created on demand.
Definition: solver.C:107
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
void bSolve(const fv::convectionScheme< scalar > &mvConvection)
Solve the Xi and regress variable equations.
void ftSolve(const fv::convectionScheme< scalar > &mvConvection)
Solve the ft equation for partially-premixed mixtures.
turbulenceThermophysicalTransportModels::unityLewisEddyDiffusivity< RASThermophysicalTransportModel< ThermophysicalTransportModel< compressibleMomentumTransportModel, fluidThermo > > > thermophysicalTransport
Definition: XiFluid.H:172
psiuMulticomponentThermo & thermo_
Definition: XiFluid.H:103
void EauSolve(const fv::convectionScheme< scalar > &mvConvection)
Solve the unburnt energy equation.
void EaSolve(const fv::convectionScheme< scalar > &mvConvection)
Solve the energy equation.
dimensionedScalar StCorr(const volScalarField &c, const surfaceScalarField &nf, const dimensionedScalar &dMgb) const
Calculate and return the turbulent flame-speed kernel correction.
Buoyancy related data for the Foam::solvers::isothermalFluid solver module when solving buoyant cases...
Definition: buoyancy.H:70
uniformDimensionedVectorField g
Gravitational acceleration.
Definition: buoyancy.H:83
const surfaceScalarField & phi
Mass-flux field.
const volScalarField & rho
Reference to the continuity density field.
A class for managing temporary objects.
Definition: tmp.H:55
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const scalar epsilon
Calculate the first temporal derivative.
Calculate the laplacian 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 matrix for the divergence of the given field and flux.
Calculate the matrix for implicit and explicit sources.
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.name(), mesh), mesh, dimensionedScalar(p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));Info<< "Creating the unstrained laminar flame speed\n"<< endl;autoPtr< laminarFlameSpeed > unstrainedLaminarFlameSpeed(laminarFlameSpeed::New(thermo))
volScalarField & b
Definition: createFields.H:25
Info<< "Creating thermophysical transport model\n"<< endl;turbulenceThermophysicalTransportModels::unityLewisEddyDiffusivity< RASThermophysicalTransportModel< ThermophysicalTransportModel< compressibleMomentumTransportModel, fluidThermo > >> thermophysicalTransport(turbulence(), thermo, true)
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:228
K
Definition: pEqn.H:75
U
Definition: pEqn.H:72
const dimensionedScalar G
Newtonian constant of gravitation.
const dimensionedScalar c
Speed of light in a vacuum.
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 > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
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 > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Type gSum(const FieldField< Field, Type > &f)
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimVolume
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
error FatalError
const dimensionSet dimArea
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31