CoEulerDdtScheme.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2018 OpenFOAM Foundation
9  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "CoEulerDdtScheme.H"
30 #include "surfaceInterpolate.H"
31 #include "fvcDiv.H"
32 #include "fvMatrices.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace fv
42 {
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 template<class Type>
47 tmp<volScalarField> CoEulerDdtScheme<Type>::CorDeltaT() const
48 {
49  const surfaceScalarField cofrDeltaT(CofrDeltaT());
50 
51  tmp<volScalarField> tcorDeltaT
52  (
53  new volScalarField
54  (
55  IOobject
56  (
57  "CorDeltaT",
58  cofrDeltaT.instance(),
59  mesh()
60  ),
61  mesh(),
62  dimensionedScalar(cofrDeltaT.dimensions(), Zero),
64  )
65  );
66 
67  volScalarField& corDeltaT = tcorDeltaT.ref();
68 
69  const labelUList& owner = mesh().owner();
70  const labelUList& neighbour = mesh().neighbour();
71 
72  forAll(owner, facei)
73  {
74  corDeltaT[owner[facei]] =
75  max(corDeltaT[owner[facei]], cofrDeltaT[facei]);
76 
77  corDeltaT[neighbour[facei]] =
78  max(corDeltaT[neighbour[facei]], cofrDeltaT[facei]);
79  }
80 
81  const surfaceScalarField::Boundary& cofrDeltaTbf =
82  cofrDeltaT.boundaryField();
83 
84  forAll(cofrDeltaTbf, patchi)
85  {
86  const fvsPatchScalarField& pcofrDeltaT = cofrDeltaTbf[patchi];
87  const fvPatch& p = pcofrDeltaT.patch();
88  const labelUList& faceCells = p.patch().faceCells();
89 
90  forAll(pcofrDeltaT, patchFacei)
91  {
92  corDeltaT[faceCells[patchFacei]] = max
93  (
94  corDeltaT[faceCells[patchFacei]],
95  pcofrDeltaT[patchFacei]
96  );
97  }
98  }
99 
100  corDeltaT.correctBoundaryConditions();
101 
102  return tcorDeltaT;
103 }
104 
105 
106 template<class Type>
107 tmp<surfaceScalarField> CoEulerDdtScheme<Type>::CofrDeltaT() const
108 {
109  const dimensionedScalar& deltaT = mesh().time().deltaT();
110 
111  const surfaceScalarField& phi =
112  static_cast<const objectRegistry&>(mesh())
113  .lookupObject<surfaceScalarField>(phiName_);
114 
115  if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
116  {
118  (
120  *(mag(phi)/mesh().magSf())
121  *deltaT
122  );
123 
124  return max(Co/maxCo_, scalar(1))/deltaT;
125  }
126  else if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
127  {
128  const volScalarField& rho =
129  static_cast<const objectRegistry&>(mesh())
130  .lookupObject<volScalarField>(rhoName_).oldTime();
131 
133  (
135  *(mag(phi)/(fvc::interpolate(rho)*mesh().magSf()))
136  *deltaT
137  );
138 
139  return max(Co/maxCo_, scalar(1))/deltaT;
140  }
141 
143  << "Incorrect dimensions of phi: " << phi.dimensions()
144  << abort(FatalError);
145 
146  return nullptr;
147 }
148 
149 
150 template<class Type>
153 (
154  const dimensioned<Type>& dt
155 )
156 {
157  const volScalarField rDeltaT(CorDeltaT());
158 
159  IOobject ddtIOobject
160  (
161  "ddt("+dt.name()+')',
162  mesh().time().timeName(),
163  mesh().thisDb()
164  );
165 
166  if (mesh().moving())
167  {
169  (
171  (
172  ddtIOobject,
173  mesh(),
175  )
176  );
177 
178  tdtdt.ref().primitiveFieldRef() =
179  rDeltaT.primitiveField()*dt.value()
180  *(1.0 - mesh().Vsc0()/mesh().Vsc());
181 
182  // Different operation on boundary v.s. internal so re-evaluate
183  // coupled boundaries
184  tdtdt.ref().boundaryFieldRef().
185  template evaluateCoupled<coupledFvPatch>();
186 
187  return tdtdt;
188  }
189  else
190  {
192  (
194  (
195  ddtIOobject,
196  mesh(),
199  )
200  );
201  }
202 }
203 
204 
205 template<class Type>
208 (
210 )
211 {
212  const volScalarField rDeltaT(CorDeltaT());
213 
214  IOobject ddtIOobject
215  (
216  "ddt("+vf.name()+')',
217  mesh().time().timeName(),
218  mesh().thisDb()
219  );
220 
221  if (mesh().moving())
222  {
224  (
226  (
227  ddtIOobject,
228  mesh(),
229  rDeltaT.dimensions()*vf.dimensions(),
230  rDeltaT.primitiveField()*
231  (
232  vf.primitiveField()
233  - vf.oldTime().primitiveField()*mesh().Vsc0()/mesh().Vsc()
234  ),
235  rDeltaT.boundaryField()*
236  (
237  vf.boundaryField() - vf.oldTime().boundaryField()
238  )
239  )
240  );
241 
242  // Different operation on boundary v.s. internal so re-evaluate
243  // coupled boundaries
244  tdtdt.ref().boundaryFieldRef().
245  template evaluateCoupled<coupledFvPatch>();
246 
247  return tdtdt;
248  }
249  else
250  {
251  return tmp<GeometricField<Type, fvPatchField, volMesh>>
252  (
253  new GeometricField<Type, fvPatchField, volMesh>
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(CorDeltaT());
272 
273  IOobject ddtIOobject
274  (
275  "ddt("+rho.name()+','+vf.name()+')',
276  mesh().time().timeName(),
277  mesh().thisDb()
278  );
279 
280  if (mesh().moving())
281  {
283  (
285  (
286  ddtIOobject,
287  mesh(),
288  rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
289  rDeltaT.primitiveField()*rho.value()*
290  (
291  vf.primitiveField()
292  - vf.oldTime().primitiveField()*mesh().Vsc0()/mesh().Vsc()
293  ),
294  rDeltaT.boundaryField()*rho.value()*
295  (
296  vf.boundaryField() - vf.oldTime().boundaryField()
297  )
298  )
299  );
300 
301  // Different operation on boundary v.s. internal so re-evaluate
302  // coupled boundaries
303  tdtdt.ref().boundaryFieldRef().
304  template evaluateCoupled<coupledFvPatch>();
305 
306  return tdtdt;
307  }
308  else
309  {
310  return tmp<GeometricField<Type, fvPatchField, volMesh>>
311  (
312  new GeometricField<Type, fvPatchField, volMesh>
313  (
314  ddtIOobject,
315  rDeltaT*rho*(vf - vf.oldTime())
316  )
317  );
318  }
319 }
320 
321 
322 template<class Type>
325 (
326  const volScalarField& rho,
328 )
329 {
330  const volScalarField rDeltaT(CorDeltaT());
331 
332  IOobject ddtIOobject
333  (
334  "ddt("+rho.name()+','+vf.name()+')',
335  mesh().time().timeName(),
336  mesh().thisDb()
337  );
338 
339  if (mesh().moving())
340  {
342  (
344  (
345  ddtIOobject,
346  mesh(),
347  rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
348  rDeltaT.primitiveField()*
349  (
350  rho.primitiveField()*vf.primitiveField()
351  - rho.oldTime().primitiveField()
352  *vf.oldTime().primitiveField()*mesh().Vsc0()/mesh().Vsc()
353  ),
354  rDeltaT.boundaryField()*
355  (
356  rho.boundaryField()*vf.boundaryField()
357  - rho.oldTime().boundaryField()
358  *vf.oldTime().boundaryField()
359  )
360  )
361  );
362 
363  // Different operation on boundary v.s. internal so re-evaluate
364  // coupled boundaries
365  tdtdt.ref().boundaryFieldRef().
366  template evaluateCoupled<coupledFvPatch>();
367 
368  return tdtdt;
369  }
370  else
371  {
372  return tmp<GeometricField<Type, fvPatchField, volMesh>>
373  (
374  new GeometricField<Type, fvPatchField, volMesh>
375  (
376  ddtIOobject,
377  rDeltaT*(rho*vf - rho.oldTime()*vf.oldTime())
378  )
379  );
380  }
381 }
382 
383 
384 template<class Type>
387 (
388  const volScalarField& alpha,
389  const volScalarField& rho,
391 )
392 {
393  const volScalarField rDeltaT(CorDeltaT());
394 
395  IOobject ddtIOobject
396  (
397  "ddt("+alpha.name()+','+rho.name()+','+vf.name()+')',
398  mesh().time().timeName(),
399  mesh().thisDb()
400  );
401 
402  if (mesh().moving())
403  {
405  (
407  (
408  ddtIOobject,
409  mesh(),
410  rDeltaT.dimensions()
411  *alpha.dimensions()*rho.dimensions()*vf.dimensions(),
412  rDeltaT.primitiveField()*
413  (
414  alpha.primitiveField()
415  *rho.primitiveField()
416  *vf.primitiveField()
417 
418  - alpha.oldTime().primitiveField()
419  *rho.oldTime().primitiveField()
420  *vf.oldTime().primitiveField()*mesh().Vsc0()/mesh().Vsc()
421  ),
422  rDeltaT.boundaryField()*
423  (
424  alpha.boundaryField()
425  *rho.boundaryField()
426  *vf.boundaryField()
427 
428  - alpha.oldTime().boundaryField()
429  *rho.oldTime().boundaryField()
430  *vf.oldTime().boundaryField()
431  )
432  )
433  );
434 
435  // Different operation on boundary v.s. internal so re-evaluate
436  // coupled boundaries
437  tdtdt.ref().boundaryFieldRef().
438  template evaluateCoupled<coupledFvPatch>();
439 
440  return tdtdt;
441  }
442  else
443  {
444  return tmp<GeometricField<Type, fvPatchField, volMesh>>
445  (
446  new GeometricField<Type, fvPatchField, volMesh>
447  (
448  ddtIOobject,
449  rDeltaT
450  *(
451  alpha*rho*vf
452  - alpha.oldTime()*rho.oldTime()*vf.oldTime()
453  )
454  )
455  );
456  }
457 }
458 
459 
460 template<class Type>
463 (
465 )
466 {
467  tmp<fvMatrix<Type>> tfvm
468  (
469  new fvMatrix<Type>
470  (
471  vf,
473  )
474  );
475 
476  fvMatrix<Type>& fvm = tfvm.ref();
477 
478  scalarField rDeltaT(CorDeltaT()().primitiveField());
479 
480  fvm.diag() = rDeltaT*mesh().Vsc();
481 
482  if (mesh().moving())
483  {
484  fvm.source() = rDeltaT*vf.oldTime().primitiveField()*mesh().Vsc0();
485  }
486  else
487  {
488  fvm.source() = rDeltaT*vf.oldTime().primitiveField()*mesh().Vsc();
489  }
490 
491  return tfvm;
492 }
493 
494 
495 template<class Type>
498 (
499  const dimensionedScalar& rho,
501 )
502 {
503  tmp<fvMatrix<Type>> tfvm
504  (
505  new fvMatrix<Type>
506  (
507  vf,
508  rho.dimensions()*vf.dimensions()*dimVol/dimTime
509  )
510  );
511  fvMatrix<Type>& fvm = tfvm.ref();
512 
513  scalarField rDeltaT(CorDeltaT()().primitiveField());
514 
515  fvm.diag() = rDeltaT*rho.value()*mesh().Vsc();
516 
517  if (mesh().moving())
518  {
519  fvm.source() = rDeltaT
520  *rho.value()*vf.oldTime().primitiveField()*mesh().Vsc0();
521  }
522  else
523  {
524  fvm.source() = rDeltaT
525  *rho.value()*vf.oldTime().primitiveField()*mesh().Vsc();
526  }
527 
528  return tfvm;
529 }
530 
531 
532 template<class Type>
535 (
536  const volScalarField& rho,
538 )
539 {
540  tmp<fvMatrix<Type>> tfvm
541  (
542  new fvMatrix<Type>
543  (
544  vf,
545  rho.dimensions()*vf.dimensions()*dimVol/dimTime
546  )
547  );
548  fvMatrix<Type>& fvm = tfvm.ref();
549 
550  scalarField rDeltaT(CorDeltaT()().primitiveField());
551 
552  fvm.diag() = rDeltaT*rho.primitiveField()*mesh().Vsc();
553 
554  if (mesh().moving())
555  {
556  fvm.source() = rDeltaT
557  *rho.oldTime().primitiveField()
558  *vf.oldTime().primitiveField()*mesh().Vsc0();
559  }
560  else
561  {
562  fvm.source() = rDeltaT
563  *rho.oldTime().primitiveField()
564  *vf.oldTime().primitiveField()*mesh().Vsc();
565  }
566 
567  return tfvm;
568 }
569 
570 
571 template<class Type>
574 (
575  const volScalarField& alpha,
576  const volScalarField& rho,
578 )
579 {
580  tmp<fvMatrix<Type>> tfvm
581  (
582  new fvMatrix<Type>
583  (
584  vf,
585  alpha.dimensions()*rho.dimensions()*vf.dimensions()*dimVol/dimTime
586  )
587  );
588  fvMatrix<Type>& fvm = tfvm.ref();
589 
590  scalarField rDeltaT(CorDeltaT()().primitiveField());
591 
592  fvm.diag() =
593  rDeltaT*alpha.primitiveField()*rho.primitiveField()*mesh().Vsc();
594 
595  if (mesh().moving())
596  {
597  fvm.source() = rDeltaT
598  *alpha.oldTime().primitiveField()
599  *rho.oldTime().primitiveField()
600  *vf.oldTime().primitiveField()*mesh().Vsc0();
601  }
602  else
603  {
604  fvm.source() = rDeltaT
605  *alpha.oldTime().primitiveField()
606  *rho.oldTime().primitiveField()
607  *vf.oldTime().primitiveField()*mesh().Vsc();
608  }
609 
610  return tfvm;
611 }
612 
613 
614 template<class Type>
617 (
620 )
621 {
622  const surfaceScalarField rDeltaT(fvc::interpolate(CorDeltaT()));
623 
624  fluxFieldType phiUf0(mesh().Sf() & Uf.oldTime());
625  fluxFieldType phiCorr
626  (
627  phiUf0 - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
628  );
629 
630  return tmp<fluxFieldType>
631  (
632  new fluxFieldType
633  (
634  IOobject
635  (
636  "ddtCorr(" + U.name() + ',' + Uf.name() + ')',
637  mesh().time().timeName(),
638  mesh().thisDb()
639  ),
640  this->fvcDdtPhiCoeff(U.oldTime(), phiUf0, phiCorr)
641  *rDeltaT*phiCorr
642  )
643  );
644 }
645 
646 
647 template<class Type>
650 (
652  const fluxFieldType& phi
653 )
654 {
655  const surfaceScalarField rDeltaT(fvc::interpolate(CorDeltaT()));
656 
657  fluxFieldType phiCorr
658  (
659  phi.oldTime() - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
660  );
661 
662  return tmp<fluxFieldType>
663  (
664  new fluxFieldType
665  (
666  IOobject
667  (
668  "ddtCorr(" + U.name() + ',' + phi.name() + ')',
669  mesh().time().timeName(),
670  mesh().thisDb()
671  ),
672  this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime(), phiCorr)
673  *rDeltaT*phiCorr
674  )
675  );
676 }
677 
678 
679 template<class Type>
682 (
683  const volScalarField& rho,
686 )
687 {
688  const surfaceScalarField rDeltaT(fvc::interpolate(CorDeltaT()));
689 
690  if
691  (
692  U.dimensions() == dimVelocity
693  && Uf.dimensions() == dimDensity*dimVelocity
694  )
695  {
697  (
698  rho.oldTime()*U.oldTime()
699  );
700 
701  fluxFieldType phiUf0(mesh().Sf() & Uf.oldTime());
702  fluxFieldType phiCorr(phiUf0 - fvc::dotInterpolate(mesh().Sf(), rhoU0));
703 
704  return tmp<fluxFieldType>
705  (
706  new fluxFieldType
707  (
708  IOobject
709  (
710  "ddtCorr("
711  + rho.name() + ',' + U.name() + ',' + Uf.name() + ')',
712  mesh().time().timeName(),
713  mesh().thisDb()
714  ),
715  this->fvcDdtPhiCoeff(rhoU0, phiUf0, phiCorr, rho.oldTime())
716  *rDeltaT*phiCorr
717  )
718  );
719  }
720  else if
721  (
722  U.dimensions() == dimDensity*dimVelocity
723  && Uf.dimensions() == dimDensity*dimVelocity
724  )
725  {
726  fluxFieldType phiUf0(mesh().Sf() & Uf.oldTime());
727  fluxFieldType phiCorr
728  (
729  phiUf0 - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
730  );
731 
732  return tmp<fluxFieldType>
733  (
734  new fluxFieldType
735  (
736  IOobject
737  (
738  "ddtCorr("
739  + rho.name() + ',' + U.name() + ',' + Uf.name() + ')',
740  mesh().time().timeName(),
741  mesh().thisDb()
742  ),
743  this->fvcDdtPhiCoeff
744  (
745  U.oldTime(),
746  phiUf0,
747  phiCorr,
748  rho.oldTime()
749  )*rDeltaT*phiCorr
750  )
751  );
752  }
753  else
754  {
756  << "dimensions of Uf are not correct"
757  << abort(FatalError);
758 
759  return fluxFieldType::null();
760  }
761 }
762 
763 
764 template<class Type>
767 (
768  const volScalarField& rho,
770  const fluxFieldType& phi
771 )
772 {
773  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
774 
775  if
776  (
777  U.dimensions() == dimVelocity
778  && phi.dimensions() == rho.dimensions()*dimVelocity*dimArea
779  )
780  {
782  (
783  rho.oldTime()*U.oldTime()
784  );
785 
786  fluxFieldType phiCorr
787  (
788  phi.oldTime() - fvc::dotInterpolate(mesh().Sf(), rhoU0)
789  );
790 
791  return tmp<fluxFieldType>
792  (
793  new fluxFieldType
794  (
795  IOobject
796  (
797  "ddtCorr("
798  + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
799  mesh().time().timeName(),
800  mesh().thisDb()
801  ),
802  this->fvcDdtPhiCoeff
803  (
804  rhoU0,
805  phi.oldTime(),
806  phiCorr,
807  rho.oldTime()
808  )*rDeltaT*phiCorr
809  )
810  );
811  }
812  else if
813  (
814  U.dimensions() == rho.dimensions()*dimVelocity
815  && phi.dimensions() == rho.dimensions()*dimVelocity*dimArea
816  )
817  {
818  fluxFieldType phiCorr
819  (
820  phi.oldTime() - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
821  );
822 
823  return tmp<fluxFieldType>
824  (
825  new fluxFieldType
826  (
827  IOobject
828  (
829  "ddtCorr("
830  + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
831  mesh().time().timeName(),
832  mesh().thisDb()
833  ),
834  this->fvcDdtPhiCoeff
835  (
836  U.oldTime(),
837  phi.oldTime(),
838  phiCorr,
839  rho.oldTime()
840  )*rDeltaT*phiCorr
841  )
842  );
843  }
844  else
845  {
847  << "dimensions of phi are not correct"
848  << abort(FatalError);
849 
850  return fluxFieldType::null();
851  }
852 }
853 
854 
855 template<class Type>
857 (
859 )
860 {
861  tmp<surfaceScalarField> tmeshPhi
862  (
864  (
865  IOobject
866  (
867  "meshPhi",
868  mesh().time().timeName(),
869  mesh().thisDb(),
873  ),
874  mesh(),
876  )
877  );
878 
879  tmeshPhi.ref().setOriented();
880 
881  return tmeshPhi;
882 }
883 
884 
885 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
886 
887 } // End namespace fv
888 
889 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
890 
891 } // End namespace Foam
892 
893 // ************************************************************************* //
fvsPatchField< scalar > fvsPatchScalarField
const Type & value() const noexcept
Return const reference to value.
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
static tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &tvf)
Interpolate field onto faces.
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:65
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
Type of boundary fields.
Generic dimensioned Type class.
Ignore writing from objectRegistry::writeObject()
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
word timeName
Definition: getTimeIndex.H:3
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
dynamicFvMesh & mesh
labelList fv(nPoints)
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
autoPtr< surfaceVectorField > Uf
const dimensionSet & dimensions() const noexcept
Return const reference to dimensions.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Calculate the divergence of the given field.
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.
const dimensionSet dimDensity
const word & name() const noexcept
Return const reference to name.
Field< Type > & source() noexcept
Definition: fvMatrix.H:533
U
Definition: pEqn.H:72
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Nothing to be read.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:204
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
scalarField & diag()
Definition: lduMatrix.C:197
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Definition: fvPatchField.H:213
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Do not request registration (bool: false)
virtual const surfaceScalarField & deltaCoeffs() const
Return reference to cell-centre difference coefficients.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
const dimensionSet & dimensions() const noexcept
Return dimensions.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity