SLTSDdtScheme.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 "SLTSDdtScheme.H"
27 #include "surfaceInterpolate.H"
28 #include "fvMatrices.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace fv
38 {
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 template<class Type>
44 (
45  scalarField& rD,
46  const surfaceScalarField& phi
47 ) const
48 {
49  const labelUList& owner = mesh().owner();
50  const labelUList& neighbour = mesh().neighbour();
51  scalarField diag(rD.size(), 0.0);
52 
53  forAll(owner, faceI)
54  {
55  if (phi[faceI] > 0.0)
56  {
57  diag[owner[faceI]] += phi[faceI];
58  rD[neighbour[faceI]] += phi[faceI];
59  }
60  else
61  {
62  diag[neighbour[faceI]] -= phi[faceI];
63  rD[owner[faceI]] -= phi[faceI];
64  }
65  }
66 
67  forAll(phi.boundaryField(), patchi)
68  {
69  const fvsPatchScalarField& pphi = phi.boundaryField()[patchi];
70  const labelUList& faceCells = pphi.patch().patch().faceCells();
71 
72  forAll(pphi, patchFacei)
73  {
74  if (pphi[patchFacei] > 0.0)
75  {
76  diag[faceCells[patchFacei]] += pphi[patchFacei];
77  }
78  else
79  {
80  rD[faceCells[patchFacei]] -= pphi[patchFacei];
81  }
82  }
83  }
84 
85  rD += (1.0/alpha_ - 2.0)*diag;
86 }
87 
88 
89 template<class Type>
91 {
92  const surfaceScalarField& phi =
93  mesh().objectRegistry::template
94  lookupObject<surfaceScalarField>(phiName_);
95 
96  const dimensionedScalar& deltaT = mesh().time().deltaT();
97 
99  (
100  new volScalarField
101  (
102  IOobject
103  (
104  "rDeltaT",
105  phi.instance(),
106  mesh()
107  ),
108  mesh(),
109  dimensionedScalar("rDeltaT", dimless/dimTime, 0.0),
110  zeroGradientFvPatchScalarField::typeName
111  )
112  );
113 
114  volScalarField& rDeltaT = trDeltaT();
115 
116  relaxedDiag(rDeltaT, phi);
117 
118  if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
119  {
120  rDeltaT.internalField() = max
121  (
122  rDeltaT.internalField()/mesh().V(),
123  scalar(1)/deltaT.value()
124  );
125  }
126  else if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
127  {
128  const volScalarField& rho =
129  mesh().objectRegistry::template lookupObject<volScalarField>
130  (
131  rhoName_
132  ).oldTime();
133 
134  rDeltaT.internalField() = max
135  (
136  rDeltaT.internalField()/(rho.internalField()*mesh().V()),
137  scalar(1)/deltaT.value()
138  );
139  }
140  else
141  {
143  << "Incorrect dimensions of phi: " << phi.dimensions()
144  << abort(FatalError);
145  }
146 
147  rDeltaT.correctBoundaryConditions();
148 
149  return trDeltaT;
150 }
151 
152 
153 template<class Type>
156 (
157  const dimensioned<Type>& dt
158 )
159 {
160  const volScalarField rDeltaT(SLrDeltaT());
161 
162  IOobject ddtIOobject
163  (
164  "ddt("+dt.name()+')',
165  mesh().time().timeName(),
166  mesh()
167  );
168 
169  if (mesh().moving())
170  {
172  (
174  (
175  ddtIOobject,
176  mesh(),
178  (
179  "0",
180  dt.dimensions()/dimTime,
182  )
183  )
184  );
185 
186  tdtdt().internalField() =
187  rDeltaT.internalField()*dt.value()*(1.0 - mesh().V0()/mesh().V());
188 
189  return tdtdt;
190  }
191  else
192  {
194  (
196  (
197  ddtIOobject,
198  mesh(),
200  (
201  "0",
202  dt.dimensions()/dimTime,
204  ),
206  )
207  );
208  }
209 }
210 
211 
212 template<class Type>
215 (
217 )
218 {
219  const volScalarField rDeltaT(SLrDeltaT());
220 
221  IOobject ddtIOobject
222  (
223  "ddt("+vf.name()+')',
224  mesh().time().timeName(),
225  mesh()
226  );
227 
228  if (mesh().moving())
229  {
231  (
233  (
234  ddtIOobject,
235  mesh(),
236  rDeltaT.dimensions()*vf.dimensions(),
237  rDeltaT.internalField()*
238  (
239  vf.internalField()
240  - vf.oldTime().internalField()*mesh().V0()/mesh().V()
241  ),
242  rDeltaT.boundaryField()*
243  (
244  vf.boundaryField() - vf.oldTime().boundaryField()
245  )
246  )
247  );
248  }
249  else
250  {
252  (
254  (
255  ddtIOobject,
256  rDeltaT*(vf - vf.oldTime())
257  )
258  );
259  }
260 }
261 
262 
263 template<class Type>
266 (
267  const dimensionedScalar& rho,
269 )
270 {
271  const volScalarField rDeltaT(SLrDeltaT());
272 
273  IOobject ddtIOobject
274  (
275  "ddt("+rho.name()+','+vf.name()+')',
276  mesh().time().timeName(),
277  mesh()
278  );
279 
280  if (mesh().moving())
281  {
283  (
285  (
286  ddtIOobject,
287  mesh(),
288  rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
289  rDeltaT.internalField()*rho.value()*
290  (
291  vf.internalField()
292  - vf.oldTime().internalField()*mesh().V0()/mesh().V()
293  ),
294  rDeltaT.boundaryField()*rho.value()*
295  (
296  vf.boundaryField() - vf.oldTime().boundaryField()
297  )
298  )
299  );
300  }
301  else
302  {
304  (
306  (
307  ddtIOobject,
308  rDeltaT*rho*(vf - vf.oldTime())
309  )
310  );
311  }
312 }
313 
314 
315 template<class Type>
318 (
319  const volScalarField& rho,
321 )
322 {
323  const volScalarField rDeltaT(SLrDeltaT());
324 
325  IOobject ddtIOobject
326  (
327  "ddt("+rho.name()+','+vf.name()+')',
328  mesh().time().timeName(),
329  mesh()
330  );
331 
332  if (mesh().moving())
333  {
335  (
337  (
338  ddtIOobject,
339  mesh(),
340  rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
341  rDeltaT.internalField()*
342  (
343  rho.internalField()*vf.internalField()
344  - rho.oldTime().internalField()
345  *vf.oldTime().internalField()*mesh().V0()/mesh().V()
346  ),
347  rDeltaT.boundaryField()*
348  (
349  rho.boundaryField()*vf.boundaryField()
350  - rho.oldTime().boundaryField()
351  *vf.oldTime().boundaryField()
352  )
353  )
354  );
355  }
356  else
357  {
359  (
361  (
362  ddtIOobject,
363  rDeltaT*(rho*vf - rho.oldTime()*vf.oldTime())
364  )
365  );
366  }
367 }
368 
369 
370 template<class Type>
373 (
374  const volScalarField& alpha,
375  const volScalarField& rho,
377 )
378 {
379  const volScalarField rDeltaT(SLrDeltaT());
380 
381  IOobject ddtIOobject
382  (
383  "ddt("+alpha.name()+','+rho.name()+','+vf.name()+')',
384  mesh().time().timeName(),
385  mesh()
386  );
387 
388  if (mesh().moving())
389  {
391  (
393  (
394  ddtIOobject,
395  mesh(),
396  rDeltaT.dimensions()
397  *alpha.dimensions()*rho.dimensions()*vf.dimensions(),
398  rDeltaT.internalField()*
399  (
400  alpha.internalField()
401  *rho.internalField()
402  *vf.internalField()
403 
404  - alpha.oldTime().internalField()
405  *rho.oldTime().internalField()
406  *vf.oldTime().internalField()*mesh().Vsc0()/mesh().Vsc()
407  ),
408  rDeltaT.boundaryField()*
409  (
410  alpha.boundaryField()
411  *rho.boundaryField()
412  *vf.boundaryField()
413 
414  - alpha.oldTime().boundaryField()
415  *rho.oldTime().boundaryField()
416  *vf.oldTime().boundaryField()
417  )
418  )
419  );
420  }
421  else
422  {
424  (
426  (
427  ddtIOobject,
428  rDeltaT
429  *(
430  alpha*rho*vf
431  - alpha.oldTime()*rho.oldTime()*vf.oldTime()
432  )
433  )
434  );
435  }
436 }
437 
438 
439 template<class Type>
442 (
444 )
445 {
446  tmp<fvMatrix<Type> > tfvm
447  (
448  new fvMatrix<Type>
449  (
450  vf,
451  vf.dimensions()*dimVol/dimTime
452  )
453  );
454 
455  fvMatrix<Type>& fvm = tfvm();
456 
457  scalarField rDeltaT(SLrDeltaT()().internalField());
458 
459  Info<< "SLTSDdtScheme<Type>::fvmDdt: max/min rDeltaT "
460  << gMax(rDeltaT) << " " << gMin(rDeltaT) << endl;
461 
462  fvm.diag() = rDeltaT*mesh().V();
463 
464  if (mesh().moving())
465  {
466  fvm.source() = rDeltaT*vf.oldTime().internalField()*mesh().V0();
467  }
468  else
469  {
470  fvm.source() = rDeltaT*vf.oldTime().internalField()*mesh().V();
471  }
472 
473  return tfvm;
474 }
475 
476 
477 template<class Type>
480 (
481  const dimensionedScalar& rho,
483 )
484 {
485  tmp<fvMatrix<Type> > tfvm
486  (
487  new fvMatrix<Type>
488  (
489  vf,
490  rho.dimensions()*vf.dimensions()*dimVol/dimTime
491  )
492  );
493  fvMatrix<Type>& fvm = tfvm();
494 
495  scalarField rDeltaT(SLrDeltaT()().internalField());
496 
497  fvm.diag() = rDeltaT*rho.value()*mesh().V();
498 
499  if (mesh().moving())
500  {
501  fvm.source() = rDeltaT
502  *rho.value()*vf.oldTime().internalField()*mesh().V0();
503  }
504  else
505  {
506  fvm.source() = rDeltaT
507  *rho.value()*vf.oldTime().internalField()*mesh().V();
508  }
509 
510  return tfvm;
511 }
512 
513 
514 template<class Type>
517 (
518  const volScalarField& rho,
520 )
521 {
522  tmp<fvMatrix<Type> > tfvm
523  (
524  new fvMatrix<Type>
525  (
526  vf,
527  rho.dimensions()*vf.dimensions()*dimVol/dimTime
528  )
529  );
530  fvMatrix<Type>& fvm = tfvm();
531 
532  scalarField rDeltaT(SLrDeltaT()().internalField());
533 
534  fvm.diag() = rDeltaT*rho.internalField()*mesh().V();
535 
536  if (mesh().moving())
537  {
538  fvm.source() = rDeltaT
539  *rho.oldTime().internalField()
540  *vf.oldTime().internalField()*mesh().V0();
541  }
542  else
543  {
544  fvm.source() = rDeltaT
545  *rho.oldTime().internalField()
546  *vf.oldTime().internalField()*mesh().V();
547  }
548 
549  return tfvm;
550 }
551 
552 
553 template<class Type>
556 (
557  const volScalarField& alpha,
558  const volScalarField& rho,
560 )
561 {
562  tmp<fvMatrix<Type> > tfvm
563  (
564  new fvMatrix<Type>
565  (
566  vf,
567  alpha.dimensions()*rho.dimensions()*vf.dimensions()*dimVol/dimTime
568  )
569  );
570  fvMatrix<Type>& fvm = tfvm();
571 
572  scalarField rDeltaT(SLrDeltaT()().internalField());
573 
574  fvm.diag() = rDeltaT*alpha.internalField()*rho.internalField()*mesh().Vsc();
575 
576  if (mesh().moving())
577  {
578  fvm.source() = rDeltaT
579  *alpha.oldTime().internalField()
580  *rho.oldTime().internalField()
581  *vf.oldTime().internalField()*mesh().Vsc0();
582  }
583  else
584  {
585  fvm.source() = rDeltaT
586  *alpha.oldTime().internalField()
587  *rho.oldTime().internalField()
588  *vf.oldTime().internalField()*mesh().Vsc();
589  }
590 
591  return tfvm;
592 }
593 
594 
595 template<class Type>
598 (
601 )
602 {
603  const surfaceScalarField rDeltaT(fvc::interpolate(SLrDeltaT()));
604 
605  fluxFieldType phiCorr
606  (
607  mesh().Sf() & (Uf.oldTime() - fvc::interpolate(U.oldTime()))
608  );
609 
610  return tmp<fluxFieldType>
611  (
612  new fluxFieldType
613  (
614  IOobject
615  (
616  "ddtCorr(" + U.name() + ',' + Uf.name() + ')',
617  mesh().time().timeName(),
618  mesh()
619  ),
620  this->fvcDdtPhiCoeff
621  (
622  U.oldTime(),
623  (mesh().Sf() & Uf.oldTime()),
624  phiCorr
625  )
626  *rDeltaT*phiCorr
627  )
628  );
629 }
630 
631 
632 template<class Type>
635 (
637  const fluxFieldType& phi
638 )
639 {
640  const surfaceScalarField rDeltaT(fvc::interpolate(SLrDeltaT()));
641 
642  fluxFieldType phiCorr
643  (
644  phi.oldTime() - (mesh().Sf() & fvc::interpolate(U.oldTime()))
645  );
646 
647  return tmp<fluxFieldType>
648  (
649  new fluxFieldType
650  (
651  IOobject
652  (
653  "ddtCorr(" + U.name() + ',' + phi.name() + ')',
654  mesh().time().timeName(),
655  mesh()
656  ),
657  this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime(), phiCorr)
658  *rDeltaT*phiCorr
659  )
660  );
661 }
662 
663 
664 template<class Type>
667 (
668  const volScalarField& rho,
671 )
672 {
673  if
674  (
675  U.dimensions() == dimVelocity
676  && Uf.dimensions() == dimDensity*dimVelocity
677  )
678  {
679  const surfaceScalarField rDeltaT(fvc::interpolate(SLrDeltaT()));
680 
682  (
683  rho.oldTime()*U.oldTime()
684  );
685 
686  fluxFieldType phiCorr
687  (
688  mesh().Sf() & (Uf.oldTime() - fvc::interpolate(rhoU0))
689  );
690 
691  return tmp<fluxFieldType>
692  (
693  new fluxFieldType
694  (
695  IOobject
696  (
697  "ddtCorr("
698  + rho.name() + ',' + U.name() + ',' + Uf.name() + ')',
699  mesh().time().timeName(),
700  mesh()
701  ),
702  this->fvcDdtPhiCoeff
703  (
704  rhoU0,
705  mesh().Sf() & Uf.oldTime(),
706  phiCorr
707  )
708  *rDeltaT*phiCorr
709  )
710  );
711  }
712  else if
713  (
714  U.dimensions() == dimDensity*dimVelocity
715  && Uf.dimensions() == dimDensity*dimVelocity
716  )
717  {
718  return fvcDdtUfCorr(U, Uf);
719  }
720  else
721  {
723  << "dimensions of Uf are not correct"
724  << abort(FatalError);
725 
726  return fluxFieldType::null();
727  }
728 }
729 
730 
731 template<class Type>
734 (
735  const volScalarField& rho,
737  const fluxFieldType& phi
738 )
739 {
740  if
741  (
742  U.dimensions() == dimVelocity
743  && phi.dimensions() == rho.dimensions()*dimVelocity*dimArea
744  )
745  {
746  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
747 
749  (
750  rho.oldTime()*U.oldTime()
751  );
752 
753  fluxFieldType phiCorr
754  (
755  phi.oldTime() - (mesh().Sf() & fvc::interpolate(rhoU0))
756  );
757 
758  return tmp<fluxFieldType>
759  (
760  new fluxFieldType
761  (
762  IOobject
763  (
764  "ddtCorr("
765  + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
766  mesh().time().timeName(),
767  mesh()
768  ),
769  this->fvcDdtPhiCoeff(rhoU0, phi.oldTime(), phiCorr)
770  *rDeltaT*phiCorr
771  )
772  );
773  }
774  else if
775  (
776  U.dimensions() == rho.dimensions()*dimVelocity
777  && phi.dimensions() == rho.dimensions()*dimVelocity*dimArea
778  )
779  {
780  return fvcDdtPhiCorr(U, phi);
781  }
782  else
783  {
785  << "dimensions of phi are not correct"
786  << abort(FatalError);
787 
788  return fluxFieldType::null();
789  }
790 }
791 
792 
793 template<class Type>
795 (
797 )
798 {
800  (
802  (
803  IOobject
804  (
805  "meshPhi",
806  mesh().time().timeName(),
807  mesh(),
810  false
811  ),
812  mesh(),
814  )
815  );
816 }
817 
818 
819 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
820 
821 } // End namespace fv
822 
823 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
824 
825 } // End namespace Foam
826 
827 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::fv::SLTSDdtScheme::fvmDdt
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
Definition: SLTSDdtScheme.C:442
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::fv::SLTSDdtScheme::SLrDeltaT
tmp< volScalarField > SLrDeltaT() const
Return the reciprocal of the stabilised local time-step.
Definition: SLTSDdtScheme.C:90
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
Foam::diag
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
Definition: pointPatchFieldFunctions.H:262
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::fvc::interpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Definition: surfaceInterpolate.C:110
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Foam::GeometricField::oldTime
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
Definition: GeometricField.C:804
Foam::dimensioned::name
const word & name() const
Return const reference to name.
Definition: dimensionedType.C:235
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
SLTSDdtScheme.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:116
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
U
U
Definition: pEqn.H:46
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::Info
messageStream Info
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fv::SLTSDdtScheme::fvcDdtUfCorr
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: SLTSDdtScheme.C:598
surfaceInterpolate.H
Surface Interpolation.
Foam::calculatedFvPatchField
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
Definition: calculatedFvPatchField.H:65
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::polyPatch::faceCells
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:340
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
rho
rho
Definition: pEqn.H:3
internalField
conserve internalField()+
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
fv
labelList fv(nPoints)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::fvsPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvsPatchField.H:278
Foam::fv::SLTSDdtScheme::meshPhi
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
Definition: SLTSDdtScheme.C:795
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::fvMatrix::source
Field< Type > & source()
Definition: fvMatrix.H:291
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
patchi
label patchi
Definition: getPatchFieldScalar.H:1
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::fvMatrix< Type >
Foam::fvPatch::patch
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
Foam::fv::SLTSDdtScheme::relaxedDiag
void relaxedDiag(scalarField &rD, const surfaceScalarField &phi) const
Calculate a relaxed diagonal from the given flux field.
Definition: SLTSDdtScheme.C:44
trDeltaT
tmp< volScalarField > trDeltaT
Definition: createRDeltaT.H:3
Foam::dimVol
const dimensionSet dimVol(dimVolume)
Definition: dimensionSets.H:59
Foam::fv::SLTSDdtScheme::fvcDdt
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
Definition: SLTSDdtScheme.C:156
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
Foam::dimensioned::dimensions
const dimensionSet & dimensions() const
Return const reference to dimensions.
Definition: dimensionedType.C:248
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Uf
Uf
Definition: pEqn.H:78
Foam::fv::SLTSDdtScheme::fvcDdtPhiCorr
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
Definition: SLTSDdtScheme.C:635
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562