Loading...
Searching...
No Matches
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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2018 OpenFOAM Foundation
9 Copyright (C) 2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "SLTSDdtScheme.H"
30#include "surfaceInterpolate.H"
31#include "fvMatrices.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37
38// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39
40namespace fv
41{
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45template<class Type>
46void SLTSDdtScheme<Type>::relaxedDiag
47(
48 scalarField& rD,
50) const
51{
52 const labelUList& owner = mesh().owner();
53 const labelUList& neighbour = mesh().neighbour();
54 scalarField diag(rD.size(), Zero);
55
56 forAll(owner, facei)
57 {
58 if (phi[facei] > 0.0)
59 {
60 diag[owner[facei]] += phi[facei];
61 rD[neighbour[facei]] += phi[facei];
62 }
63 else
64 {
65 diag[neighbour[facei]] -= phi[facei];
66 rD[owner[facei]] -= phi[facei];
67 }
68 }
69
70 forAll(phi.boundaryField(), patchi)
71 {
72 const fvsPatchScalarField& pphi = phi.boundaryField()[patchi];
73 const labelUList& faceCells = pphi.patch().patch().faceCells();
74
75 forAll(pphi, patchFacei)
76 {
77 if (pphi[patchFacei] > 0.0)
78 {
79 diag[faceCells[patchFacei]] += pphi[patchFacei];
80 }
81 else
82 {
83 rD[faceCells[patchFacei]] -= pphi[patchFacei];
84 }
85 }
86 }
87
88 rD += (1.0/alpha_ - 2.0)*diag;
89}
90
91
92template<class Type>
93tmp<volScalarField> SLTSDdtScheme<Type>::SLrDeltaT() const
94{
95 const surfaceScalarField& phi =
96 mesh().objectRegistry::template
97 lookupObject<surfaceScalarField>(phiName_);
98
99 const dimensionedScalar& deltaT = mesh().time().deltaT();
100
101 tmp<volScalarField> trDeltaT
102 (
104 (
105 IOobject
106 (
107 "rDeltaT",
108 phi.instance(),
109 mesh()
110 ),
111 mesh(),
114 )
115 );
116
117 volScalarField& rDeltaT = trDeltaT.ref();
118
119 relaxedDiag(rDeltaT, phi);
120
121 if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
122 {
123 rDeltaT.primitiveFieldRef() = max
124 (
125 rDeltaT.primitiveField()/mesh().V(),
126 scalar(1)/deltaT.value()
127 );
128 }
129 else if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
130 {
131 const volScalarField& rho =
132 mesh().objectRegistry::template lookupObject<volScalarField>
133 (
134 rhoName_
135 ).oldTime();
136
137 rDeltaT.primitiveFieldRef() = max
138 (
139 rDeltaT.primitiveField()/(rho.primitiveField()*mesh().V()),
140 scalar(1)/deltaT.value()
141 );
142 }
143 else
144 {
146 << "Incorrect dimensions of phi: " << phi.dimensions()
147 << abort(FatalError);
148 }
149
150 rDeltaT.correctBoundaryConditions();
152 return trDeltaT;
153}
154
155
156template<class Type>
159(
160 const dimensioned<Type>& dt
161)
162{
163 const volScalarField rDeltaT(SLrDeltaT());
164
165 IOobject ddtIOobject
166 (
167 "ddt("+dt.name()+')',
168 mesh().time().timeName(),
169 mesh().thisDb()
170 );
171
172 if (mesh().moving())
173 {
175 (
177 (
178 ddtIOobject,
179 mesh(),
181 )
182 );
183
184 tdtdt.ref().primitiveFieldRef() =
185 rDeltaT.primitiveField()*dt.value()*(1.0 - mesh().V0()/mesh().V());
186
187 // Different operation on boundary v.s. internal so re-evaluate
188 // coupled boundaries
189 tdtdt.ref().boundaryFieldRef().
190 template evaluateCoupled<coupledFvPatch>();
191
192 return tdtdt;
193 }
194 else
195 {
197 (
199 (
200 ddtIOobject,
201 mesh(),
204 )
205 );
206 }
207}
208
209
210template<class Type>
213(
215)
216{
217 const volScalarField rDeltaT(SLrDeltaT());
218
219 IOobject ddtIOobject
220 (
221 "ddt("+vf.name()+')',
222 mesh().time().timeName(),
223 mesh().thisDb()
224 );
225
226 if (mesh().moving())
227 {
229 (
231 (
232 ddtIOobject,
233 mesh(),
234 rDeltaT.dimensions()*vf.dimensions(),
235 rDeltaT.primitiveField()*
236 (
237 vf.primitiveField()
238 - vf.oldTime().primitiveField()*mesh().V0()/mesh().V()
239 ),
240 rDeltaT.boundaryField()*
241 (
242 vf.boundaryField() - vf.oldTime().boundaryField()
243 )
244 )
245 );
246
247 // Different operation on boundary v.s. internal so re-evaluate
248 // coupled boundaries
249 tdtdt.ref().boundaryFieldRef().
250 template evaluateCoupled<coupledFvPatch>();
251
252 return tdtdt;
253 }
254 else
255 {
256 return tmp<GeometricField<Type, fvPatchField, volMesh>>
257 (
258 new GeometricField<Type, fvPatchField, volMesh>
259 (
260 ddtIOobject,
261 rDeltaT*(vf - vf.oldTime())
262 )
263 );
264 }
265}
266
267
268template<class Type>
271(
272 const dimensionedScalar& rho,
274)
275{
276 const volScalarField rDeltaT(SLrDeltaT());
277
278 IOobject ddtIOobject
279 (
280 "ddt("+rho.name()+','+vf.name()+')',
281 mesh().time().timeName(),
282 mesh().thisDb()
283 );
284
285 if (mesh().moving())
286 {
288 (
290 (
291 ddtIOobject,
292 mesh(),
293 rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
294 rDeltaT.primitiveField()*rho.value()*
295 (
296 vf.primitiveField()
297 - vf.oldTime().primitiveField()*mesh().V0()/mesh().V()
298 ),
299 rDeltaT.boundaryField()*rho.value()*
300 (
301 vf.boundaryField() - vf.oldTime().boundaryField()
302 )
303 )
304 );
305
306 // Different operation on boundary v.s. internal so re-evaluate
307 // coupled boundaries
308 tdtdt.ref().boundaryFieldRef().
309 template evaluateCoupled<coupledFvPatch>();
310
311 return tdtdt;
312 }
313 else
314 {
315 return tmp<GeometricField<Type, fvPatchField, volMesh>>
316 (
317 new GeometricField<Type, fvPatchField, volMesh>
318 (
319 ddtIOobject,
320 rDeltaT*rho*(vf - vf.oldTime())
321 )
322 );
323 }
324}
325
326
327template<class Type>
330(
331 const volScalarField& rho,
333)
334{
335 const volScalarField rDeltaT(SLrDeltaT());
336
337 IOobject ddtIOobject
338 (
339 "ddt("+rho.name()+','+vf.name()+')',
340 mesh().time().timeName(),
341 mesh().thisDb()
342 );
343
344 if (mesh().moving())
345 {
347 (
349 (
350 ddtIOobject,
351 mesh(),
352 rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
353 rDeltaT.primitiveField()*
354 (
355 rho.primitiveField()*vf.primitiveField()
356 - rho.oldTime().primitiveField()
357 *vf.oldTime().primitiveField()*mesh().V0()/mesh().V()
358 ),
359 rDeltaT.boundaryField()*
360 (
361 rho.boundaryField()*vf.boundaryField()
362 - rho.oldTime().boundaryField()
363 *vf.oldTime().boundaryField()
364 )
365 )
366 );
367
368 // Different operation on boundary v.s. internal so re-evaluate
369 // coupled boundaries
370 tdtdt.ref().boundaryFieldRef().
371 template evaluateCoupled<coupledFvPatch>();
372
373 return tdtdt;
374 }
375 else
376 {
377 return tmp<GeometricField<Type, fvPatchField, volMesh>>
378 (
379 new GeometricField<Type, fvPatchField, volMesh>
380 (
381 ddtIOobject,
382 rDeltaT*(rho*vf - rho.oldTime()*vf.oldTime())
383 )
384 );
385 }
386}
387
388
389template<class Type>
392(
393 const volScalarField& alpha,
394 const volScalarField& rho,
396)
397{
398 const volScalarField rDeltaT(SLrDeltaT());
399
400 IOobject ddtIOobject
401 (
402 "ddt("+alpha.name()+','+rho.name()+','+vf.name()+')',
403 mesh().time().timeName(),
404 mesh().thisDb()
405 );
406
407 if (mesh().moving())
408 {
410 (
412 (
413 ddtIOobject,
414 mesh(),
415 rDeltaT.dimensions()
416 *alpha.dimensions()*rho.dimensions()*vf.dimensions(),
417 rDeltaT.primitiveField()*
418 (
419 alpha.primitiveField()
420 *rho.primitiveField()
421 *vf.primitiveField()
422
423 - alpha.oldTime().primitiveField()
424 *rho.oldTime().primitiveField()
425 *vf.oldTime().primitiveField()*mesh().Vsc0()/mesh().Vsc()
426 ),
427 rDeltaT.boundaryField()*
428 (
429 alpha.boundaryField()
430 *rho.boundaryField()
431 *vf.boundaryField()
432
433 - alpha.oldTime().boundaryField()
434 *rho.oldTime().boundaryField()
435 *vf.oldTime().boundaryField()
436 )
437 )
438 );
439
440 // Different operation on boundary v.s. internal so re-evaluate
441 // coupled boundaries
442 tdtdt.ref().boundaryFieldRef().
443 template evaluateCoupled<coupledFvPatch>();
444
445 return tdtdt;
446 }
447 else
448 {
449 return tmp<GeometricField<Type, fvPatchField, volMesh>>
450 (
451 new GeometricField<Type, fvPatchField, volMesh>
452 (
453 ddtIOobject,
454 rDeltaT
455 *(
456 alpha*rho*vf
457 - alpha.oldTime()*rho.oldTime()*vf.oldTime()
458 )
459 )
460 );
461 }
462}
463
464
465template<class Type>
468(
470)
471{
473 (
475 (
476 vf,
478 )
479 );
480
481 fvMatrix<Type>& fvm = tfvm.ref();
482
483 scalarField rDeltaT(SLrDeltaT()().primitiveField());
484
485 fvm.diag() = rDeltaT*mesh().V();
486
487 if (mesh().moving())
488 {
489 fvm.source() = rDeltaT*vf.oldTime().primitiveField()*mesh().V0();
490 }
491 else
492 {
493 fvm.source() = rDeltaT*vf.oldTime().primitiveField()*mesh().V();
494 }
496 return tfvm;
497}
498
499
500template<class Type>
503(
504 const dimensionedScalar& rho,
506)
507{
509 (
511 (
512 vf,
513 rho.dimensions()*vf.dimensions()*dimVol/dimTime
514 )
515 );
516 fvMatrix<Type>& fvm = tfvm.ref();
517
518 scalarField rDeltaT(SLrDeltaT()().primitiveField());
519
520 fvm.diag() = rDeltaT*rho.value()*mesh().V();
521
522 if (mesh().moving())
523 {
524 fvm.source() = rDeltaT
525 *rho.value()*vf.oldTime().primitiveField()*mesh().V0();
526 }
527 else
528 {
529 fvm.source() = rDeltaT
530 *rho.value()*vf.oldTime().primitiveField()*mesh().V();
531 }
533 return tfvm;
534}
535
536
537template<class Type>
540(
541 const volScalarField& rho,
543)
544{
546 (
548 (
549 vf,
550 rho.dimensions()*vf.dimensions()*dimVol/dimTime
551 )
552 );
553 fvMatrix<Type>& fvm = tfvm.ref();
554
555 scalarField rDeltaT(SLrDeltaT()().primitiveField());
556
557 fvm.diag() = rDeltaT*rho.primitiveField()*mesh().V();
558
559 if (mesh().moving())
560 {
561 fvm.source() = rDeltaT
562 *rho.oldTime().primitiveField()
563 *vf.oldTime().primitiveField()*mesh().V0();
564 }
565 else
566 {
567 fvm.source() = rDeltaT
568 *rho.oldTime().primitiveField()
569 *vf.oldTime().primitiveField()*mesh().V();
570 }
572 return tfvm;
573}
574
575
576template<class Type>
579(
580 const volScalarField& alpha,
581 const volScalarField& rho,
583)
584{
586 (
588 (
589 vf,
590 alpha.dimensions()*rho.dimensions()*vf.dimensions()*dimVol/dimTime
591 )
592 );
593 fvMatrix<Type>& fvm = tfvm.ref();
594
595 scalarField rDeltaT(SLrDeltaT()().primitiveField());
596
597 fvm.diag() =
598 rDeltaT*alpha.primitiveField()*rho.primitiveField()*mesh().Vsc();
599
600 if (mesh().moving())
601 {
602 fvm.source() = rDeltaT
603 *alpha.oldTime().primitiveField()
604 *rho.oldTime().primitiveField()
605 *vf.oldTime().primitiveField()*mesh().Vsc0();
606 }
607 else
608 {
609 fvm.source() = rDeltaT
610 *alpha.oldTime().primitiveField()
611 *rho.oldTime().primitiveField()
612 *vf.oldTime().primitiveField()*mesh().Vsc();
613 }
615 return tfvm;
616}
617
618
619template<class Type>
622(
625)
626{
627 const surfaceScalarField rDeltaT(fvc::interpolate(SLrDeltaT()));
628
629 fluxFieldType phiUf0(mesh().Sf() & Uf.oldTime());
630 fluxFieldType phiCorr
631 (
632 phiUf0 - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
633 );
634
635 return tmp<fluxFieldType>
636 (
637 new fluxFieldType
638 (
640 (
641 "ddtCorr(" + U.name() + ',' + Uf.name() + ')',
642 mesh().time().timeName(),
643 mesh().thisDb()
644 ),
645 this->fvcDdtPhiCoeff(U.oldTime(), phiUf0, phiCorr)
646 *rDeltaT*phiCorr
648 );
649}
650
651
652template<class Type>
655(
657 const fluxFieldType& phi
658)
659{
660 const surfaceScalarField rDeltaT(fvc::interpolate(SLrDeltaT()));
661
662 fluxFieldType phiCorr
663 (
664 phi.oldTime() - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
665 );
666
667 return tmp<fluxFieldType>
668 (
669 new fluxFieldType
670 (
672 (
673 "ddtCorr(" + U.name() + ',' + phi.name() + ')',
674 mesh().time().timeName(),
675 mesh().thisDb()
676 ),
677 this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime(), phiCorr)
678 *rDeltaT*phiCorr
680 );
681}
682
683
684template<class Type>
687(
688 const volScalarField& rho,
691)
692{
693 const surfaceScalarField rDeltaT(fvc::interpolate(SLrDeltaT()));
694
695 if
696 (
697 U.dimensions() == dimVelocity
698 && Uf.dimensions() == dimDensity*dimVelocity
699 )
700 {
702 (
703 rho.oldTime()*U.oldTime()
704 );
705
706 fluxFieldType phiUf0(mesh().Sf() & Uf.oldTime());
707 fluxFieldType phiCorr(phiUf0 - fvc::dotInterpolate(mesh().Sf(), rhoU0));
708
709 return tmp<fluxFieldType>
710 (
711 new fluxFieldType
712 (
714 (
715 "ddtCorr("
716 + rho.name() + ',' + U.name() + ',' + Uf.name() + ')',
717 mesh().time().timeName(),
718 mesh().thisDb()
719 ),
720 this->fvcDdtPhiCoeff(rhoU0, phiUf0, phiCorr, rho.oldTime())
721 *rDeltaT*phiCorr
722 )
723 );
724 }
725 else if
726 (
727 U.dimensions() == dimDensity*dimVelocity
728 && Uf.dimensions() == dimDensity*dimVelocity
729 )
730 {
731 fluxFieldType phiUf0(mesh().Sf() & Uf.oldTime());
732 fluxFieldType phiCorr
733 (
734 phiUf0 - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
735 );
736
737 return tmp<fluxFieldType>
738 (
739 new fluxFieldType
740 (
742 (
743 "ddtCorr("
744 + rho.name() + ',' + U.name() + ',' + Uf.name() + ')',
745 mesh().time().timeName(),
746 mesh().thisDb()
747 ),
748 this->fvcDdtPhiCoeff
749 (
750 U.oldTime(),
751 phiUf0,
752 phiCorr,
753 rho.oldTime()
754 )*rDeltaT*phiCorr
755 )
756 );
757 }
758 else
759 {
761 << "dimensions of Uf are not correct"
762 << abort(FatalError);
763
765 }
766}
767
768
769template<class Type>
772(
773 const volScalarField& rho,
775 const fluxFieldType& phi
776)
777{
778 dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
779
780 if
781 (
782 U.dimensions() == dimVelocity
783 && phi.dimensions() == rho.dimensions()*dimVelocity*dimArea
784 )
785 {
787 (
788 rho.oldTime()*U.oldTime()
789 );
790
791 fluxFieldType phiCorr
792 (
793 phi.oldTime() - fvc::dotInterpolate(mesh().Sf(), rhoU0)
794 );
795
796 return tmp<fluxFieldType>
797 (
798 new fluxFieldType
799 (
801 (
802 "ddtCorr("
803 + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
804 mesh().time().timeName(),
805 mesh().thisDb()
806 ),
807 this->fvcDdtPhiCoeff
808 (
809 rhoU0,
810 phi.oldTime(),
811 phiCorr,
812 rho.oldTime()
813 )*rDeltaT*phiCorr
814 )
815 );
816 }
817 else if
818 (
819 U.dimensions() == rho.dimensions()*dimVelocity
820 && phi.dimensions() == rho.dimensions()*dimVelocity*dimArea
821 )
822 {
823 fluxFieldType phiCorr
824 (
825 phi.oldTime() - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
826 );
827
828 return tmp<fluxFieldType>
829 (
830 new fluxFieldType
831 (
832 IOobject
833 (
834 "ddtCorr("
835 + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
836 mesh().time().timeName(),
837 mesh().thisDb()
838 ),
839 this->fvcDdtPhiCoeff
840 (
841 U.oldTime(),
842 phi.oldTime(),
843 phiCorr,
844 rho.oldTime()
845 )*rDeltaT*phiCorr
846 )
847 );
848 }
849 else
850 {
852 << "dimensions of phi are not correct"
853 << abort(FatalError);
855 return fluxFieldType::null();
856 }
857}
858
859
860template<class Type>
862(
864)
865{
867 (
869 (
871 (
872 "meshPhi",
873 mesh().time().timeName(),
874 mesh().thisDb(),
878 ),
879 mesh(),
881 )
882 );
883
884 tmeshPhi.ref().setOriented();
885
886 return tmeshPhi;
887}
888
889
890// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
891
892} // End namespace fv
893
894// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
895
896} // End namespace Foam
897
898// ************************************************************************* //
const dimensionSet & dimensions() const noexcept
Return dimensions.
Generic GeometricField class.
static const this_type & null() noexcept
Return a null GeometricField (reference to a nullObject).
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
Generic dimensioned Type class.
const dimensionSet & dimensions() const noexcept
Return const reference to dimensions.
const word & name() const noexcept
Return const reference to name.
const Type & value() const noexcept
Return const reference to value.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
Field< Type > & source() noexcept
Definition fvMatrix.H:535
static const word & calculatedType() noexcept
The type name for calculated patch fields.
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
ddtScheme< Type >::fluxFieldType fluxFieldType
const fvMesh & mesh() const
Return mesh reference.
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
tmp< surfaceScalarField > fvcDdtPhiCoeff(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
Definition ddtScheme.C:122
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
U
Definition pEqn.H:72
dynamicFvMesh & mesh
tmp< volScalarField > trDeltaT
autoPtr< surfaceVectorField > Uf
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
A special matrix type and solver, designed for finite volume solutions of scalar equations.
word timeName
Definition getTimeIndex.H:3
const expr V(m.psi().mesh().V())
Namespace for finite-volume.
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.
static tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &tvf)
Interpolate field onto faces.
Namespace of functions to calculate implicit derivatives returning a matrix.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
const dimensionSet dimVelocity
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
errorManip< error > abort(error &err)
Definition errorManip.H:139
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
const dimensionSet dimDensity
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvsPatchField< scalar > fvsPatchScalarField
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
volScalarField & alpha
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299