Loading...
Searching...
No Matches
fusedGaussLaplacianSchemes.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-2016 OpenFOAM Foundation
9 Copyright (C) 2024 M. Janssens
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\*---------------------------------------------------------------------------*/
30#include "fvMesh.H"
31#include "fvMatrixExpression.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35makeFvLaplacianScheme(fusedGaussLaplacianScheme)
36
37#define declareFvmLaplacianScalarGamma(Type) \
38 \
39template<> \
40Foam::tmp<Foam::fvMatrix<Foam::Type>> \
41Foam::fv::fusedGaussLaplacianScheme<Foam::Type, Foam::scalar>:: \
42fvmLaplacian \
43( \
44 const GeometricField<scalar, fvPatchField, volMesh>& gamma, \
45 const GeometricField<Type, fvPatchField, volMesh>& vf \
46) \
47{ \
48 DebugPout<< "fusedGaussLaplacianScheme::fvmLaplacian on " << vf.name() \
49 << " with scalar gamma " << gamma.name() << " and ET" << endl; \
50 \
51 const fvMesh& mesh = this->mesh(); \
52 \
53 const auto weights = this->tinterpGammaScheme_().weights(gamma).expr(); \
54 const auto gammaMagSf = \
55 Expression::interpolate(gamma.expr(), weights, mesh) \
56 * mesh.magSf().expr(); \
57 /* For compatibility with linearInterpolate : avoid orientation. TBD. */ \
58 const_cast<orientedType&>(gammaMagSf.oriented()).setOriented(false); \
59 const auto deltaCoeffs = this->tsnGradScheme_().deltaCoeffs(vf).expr(); \
60 \
61 tmp<fvMatrix<Type>> tfvm \
62 ( \
63 new fvMatrix<Type> \
64 ( \
65 vf, \
66 gamma.dimensions()*mesh.magSf().dimensions()*vf.dimensions() \
67 ) \
68 ); \
69 fvMatrix<Type>& fvm = tfvm.ref(); \
70 \
71 Expression::fvmLaplacianUncorrected(fvm, gammaMagSf, deltaCoeffs); \
72 \
73 if (this->tsnGradScheme_().corrected()) \
74 { \
75 const auto corr(this->tsnGradScheme_().correction(vf).expr()); \
76 fvmCorrection(fvm, gamma.dimensions(), gammaMagSf, corr); \
77 } \
78 \
79 return tfvm; \
80} \
81 \
82template<> \
83Foam::tmp<Foam::fvMatrix<Foam::Type>> \
84Foam::fv::fusedGaussLaplacianScheme<Foam::Type, Foam::scalar>:: \
85fvmLaplacian \
86( \
87 const GeometricField<scalar, fvsPatchField, surfaceMesh>& gamma, \
88 const GeometricField<Type, fvPatchField, volMesh>& vf \
89) \
90{ \
91 DebugPout<< "fusedGaussLaplacianScheme::fvmLaplacian on " << vf.name() \
92 << " with interpolated gamma " << gamma.name() << " and ET" << endl; \
93 \
94 const fvMesh& mesh = this->mesh(); \
95 \
96 const auto gammaMagSf = gamma.expr()* mesh.magSf().expr(); \
97 const auto deltaCoeffs = this->tsnGradScheme_().deltaCoeffs(vf).expr(); \
98 \
99 tmp<fvMatrix<Type>> tfvm \
100 ( \
101 new fvMatrix<Type> \
102 ( \
103 vf, \
104 gamma.dimensions()*mesh.magSf().dimensions()*vf.dimensions() \
105 ) \
106 ); \
107 fvMatrix<Type>& fvm = tfvm.ref(); \
108 \
109 Expression::fvmLaplacianUncorrected(fvm, gammaMagSf, deltaCoeffs); \
110 \
111 if (this->tsnGradScheme_().corrected()) \
112 { \
113 const auto corr(this->tsnGradScheme_().correction(vf).expr()); \
114 fvmCorrection(fvm, gamma.dimensions(), gammaMagSf, corr); \
115 } \
116 \
117 return tfvm; \
118} \
119 \
120template<> \
121Foam::tmp<Foam::GeometricField<Foam::Type, Foam::fvPatchField, Foam::volMesh>> \
122Foam::fv::fusedGaussLaplacianScheme<Foam::Type, Foam::scalar>:: \
123fvcLaplacian \
124( \
125 const GeometricField<scalar, fvsPatchField, surfaceMesh>& gamma, \
126 const GeometricField<Type, fvPatchField, volMesh>& vf \
127) \
128{ \
129 DebugPout<< "fvcLaplacian on " << vf.name() \
130 << " with scalar gamma " << gamma.name() << endl; \
131 \
132 const fvMesh& mesh = this->mesh(); \
133 \
134 tmp<GeometricField<Type, fvPatchField, volMesh>> tLaplacian \
135 ( \
136 fvc::div(gamma*this->tsnGradScheme_().snGrad(vf)*mesh.magSf()) \
137 ); \
138 \
139 tLaplacian.ref().rename \
140 ( \
141 "laplacian(" + gamma.name() + ',' + vf.name() + ')' \
142 ); \
143 \
144 return tLaplacian; \
145}
146
147
148template<>
150<
152>
154(
155 const GeometricField<scalar, fvPatchField, volMesh>& gamma,
156 const GeometricField<scalar, fvPatchField, volMesh>& vf
157)
158{
159 typedef scalar Type;
160 typedef GeometricField<Type, fvPatchField, volMesh> FieldType;
161 typedef GeometricField<Type, fvsPatchField, surfaceMesh> SurfaceFieldType;
162 typedef typename outerProduct<vector, Type>::type GradType;
163 typedef GeometricField<GradType, fvPatchField, volMesh> GradFieldType;
164
166 << "fusedGaussLaplacianScheme<scalar, scalar>::fvcLaplacian"
167 << " on " << vf.name() << " with gamma " << gamma.name() << endl;
168
169 const fvMesh& mesh = vf.mesh();
170
171 tmp<FieldType> tresult
172 (
173 new FieldType
174 (
175 IOobject
176 (
177 "laplacian(" + vf.name() + ')',
178 vf.instance(),
179 mesh,
180 IOobject::NO_READ,
181 IOobject::NO_WRITE
182 ),
183 mesh,
184 dimensioned<Type>
185 (
186 gamma.dimensions()*vf.dimensions()/dimArea, Zero
187 ),
188 fvPatchFieldBase::extrapolatedCalculatedType()
189 )
190 );
191 FieldType& result = tresult.ref();
192
193 const auto tweights(this->tinterpGammaScheme_().weights(gamma));
194 const auto& weights = tweights();
195 const auto tdeltaCoeffs(this->tsnGradScheme_().deltaCoeffs(vf));
196 const auto& deltaCoeffs = tdeltaCoeffs();
197
198 if (this->tsnGradScheme_().corrected())
199 {
200 // Calculate sn gradient
201 tmp<SurfaceFieldType> tfaceGrad
202 (
203 new SurfaceFieldType
204 (
205 IOobject
206 (
207 "snGradCorr("+vf.name()+')',
208 vf.instance(),
209 mesh,
210 IOobject::NO_READ,
211 IOobject::NO_WRITE
212 ),
213 mesh,
214 vf.dimensions()
215 )
216 );
217
218 {
219 // Calculate gradient
220 tmp<GradFieldType> tgGrad
221 (
222 gradScheme<Type>::New
223 (
224 mesh,
225 mesh.gradScheme("grad(" + vf.name() + ')')
226 )().grad(vf, "grad(" + vf.name() + ')')
227 );
228 const auto& gGrad = tgGrad();
229
230 // Doing a dotinterpolate with nonOrthCorrectionVectors
231 const auto dotInterpolate = [&]
232 (
233 const vector& area,
234 const scalar lambda,
235
236 const GradType& ownVal,
237 const GradType& neiVal,
238
239 const vector& dotVector,
240
241 Type& result
242 )
243 {
244 result = dotVector&(lambda*(ownVal - neiVal) + neiVal);
245 };
246
247 fvc::interpolate
248 (
249 mesh.surfaceInterpolation::weights(), // linear interpolation
250 gGrad, // volume field
251 mesh.nonOrthCorrectionVectors(),// surface multiplier
252 dotInterpolate,
253 tfaceGrad.ref()
254 );
255 }
256 const auto& faceGrad = tfaceGrad();
257
258 const auto snGrad = [&]
259 (
260 const vector& Sf,
261
262 const scalar weight,
263 const scalar ownGamma,
264 const scalar neiGamma,
265
266 const scalar dc,
267 const Type& ownVal,
268 const Type& neiVal,
269
270 const Type& correction
271 ) -> Type
272 {
273 const auto snGrad(dc*(neiVal-ownVal) + correction);
274 const scalar faceGamma(weight*(ownGamma-neiGamma)+neiGamma);
275 return mag(Sf)*faceGamma*snGrad;
276 };
277
278 fvc::surfaceSnSum
279 (
280 weights, // gamma weights
281 gamma,
282
283 deltaCoeffs,
284 vf,
285
286 faceGrad, // face-based addition
287
288 snGrad,
289
290 result,
291 false // avoid boundary evaluation until volume division
292 );
293 }
294 else
295 {
296 const auto snGrad = [&]
297 (
298 const vector& Sf,
299
300 const scalar weight,
301 const scalar ownGamma,
302 const scalar neiGamma,
303
304 const scalar dc,
305 const Type& ownVal,
306 const Type& neiVal
307 ) -> Type
308 {
309 const auto snGrad(dc*(neiVal-ownVal));
310 const scalar faceGamma(weight*(ownGamma-neiGamma)+neiGamma);
311 return mag(Sf)*faceGamma*snGrad;
312 };
313
314 fvc::surfaceSnSum
315 (
316 weights,
317 gamma,
318
319 deltaCoeffs,
320 vf,
321
322 snGrad,
323
324 result,
325 false // avoid boundary evaluation until volume division
326 );
327 }
328
329 result.primitiveFieldRef() /= mesh.V();
330 result.correctBoundaryConditions();
331
332 return tresult;
333}
335
336template<>
338<
340>
342(
345)
346{
348 << "fusedGaussLaplacianScheme<vector, scalar>::fvcLaplacian"
349 << " on " << vf.name() << " with gamma " << gamma.name() << endl;
350
351 typedef vector Type;
354 typedef typename outerProduct<vector, Type>::type GradType;
356
357 const fvMesh& mesh = vf.mesh();
358
359 tmp<FieldType> tresult
360 (
361 new FieldType
362 (
364 (
365 "laplacian(" + vf.name() + ')',
366 vf.instance(),
367 mesh,
370 ),
371 mesh,
373 (
374 gamma.dimensions()*vf.dimensions()/dimArea, Zero
375 ),
377 )
378 );
379 FieldType& result = tresult.ref();
380
381 const auto tweights(this->tinterpGammaScheme_().weights(gamma));
382 const auto& weights = tweights();
383 const auto tdeltaCoeffs(this->tsnGradScheme_().deltaCoeffs(vf));
384 const auto& deltaCoeffs = tdeltaCoeffs();
385
386 if (this->tsnGradScheme_().corrected())
387 {
388 // Calculate sn gradient
389 tmp<SurfaceFieldType> tfaceGrad
390 (
391 new SurfaceFieldType
392 (
394 (
395 "snGradCorr("+vf.name()+')',
396 vf.instance(),
397 mesh,
400 ),
401 mesh,
402 vf.dimensions()
403 )
404 );
405
406 {
407 // Calculate gradient
408 tmp<GradFieldType> tgGrad
409 (
411 (
412 mesh,
413 mesh.gradScheme("grad(" + vf.name() + ')')
414 )().grad(vf, "grad(" + vf.name() + ')')
415 );
416 const auto& gGrad = tgGrad();
417
418 // Doing a dotinterpolate with nonOrthCorrectionVectors
419 const auto dotInterpolate = [&]
420 (
421 const vector& area,
422 const scalar lambda,
423
424 const GradType& ownVal,
425 const GradType& neiVal,
426
427 const vector& dotVector,
428
429 Type& result
430 )
431 {
432 result = dotVector&(lambda*(ownVal - neiVal) + neiVal);
433 };
434
436 (
437 mesh.surfaceInterpolation::weights(), // linear interpolation
438 gGrad, // volume field
439 mesh.nonOrthCorrectionVectors(),// surface multiplier
440 dotInterpolate,
441 tfaceGrad.ref()
442 );
443 }
444 const auto& faceGrad = tfaceGrad();
445
446 const auto snGrad = [&]
447 (
448 const vector& Sf,
449
450 const scalar weight,
451 const scalar ownGamma,
452 const scalar neiGamma,
453
454 const scalar dc,
455 const Type& ownVal,
456 const Type& neiVal,
457
458 const Type& correction
459 ) -> Type
460 {
461 const auto snGrad(dc*(neiVal-ownVal) + correction);
462 const scalar faceGamma(weight*(ownGamma-neiGamma)+neiGamma);
463 return mag(Sf)*faceGamma*snGrad;
464 };
465
467 (
468 weights, // gamma weights
469 gamma,
470
471 deltaCoeffs,
472 vf,
473
474 faceGrad, // face-based addition
475
476 snGrad,
477
478 result,
479 false // avoid boundary evaluation until volume division
480 );
481 }
482 else
483 {
484 const auto snGrad = [&]
485 (
486 const vector& Sf,
487
488 const scalar weight,
489 const scalar ownGamma,
490 const scalar neiGamma,
491
492 const scalar dc,
493 const Type& ownVal,
494 const Type& neiVal
495 ) -> Type
496 {
497 const auto snGrad(dc*(neiVal-ownVal));
498 const scalar faceGamma(weight*(ownGamma-neiGamma)+neiGamma);
499 return mag(Sf)*faceGamma*snGrad;
500 };
501
503 (
504 weights, // gamma weights
505 gamma,
506
507 deltaCoeffs,
508 vf,
509
510 snGrad,
511
512 result,
513 false // avoid boundary evaluation until volume division
514 );
515 }
516
517 result.primitiveFieldRef() /= mesh.V();
518 result.correctBoundaryConditions();
519
520 return tresult;
521}
523
524template<>
526<
528>
530(
532)
533{
534 typedef scalar Type;
535
538 typedef typename outerProduct<vector, Type>::type GradType;
540
541 const fvMesh& mesh = vf.mesh();
542
543 tmp<FieldType> tresult
544 (
545 new FieldType
546 (
548 (
549 "laplacian(" + vf.name() + ')',
550 vf.instance(),
551 mesh,
554 ),
555 mesh,
558 )
559 );
560 FieldType& result = tresult.ref();
561
563 << "fusedGaussLaplacianScheme<scalar, GType>::fvcLaplacian on "
564 << vf.name()
565 << " to generate " << result.name() << endl;
566
567
568 const auto tdeltaCoeffs(this->tsnGradScheme_().deltaCoeffs(vf));
569 const auto& deltaCoeffs = tdeltaCoeffs();
570
571
572 if (this->tsnGradScheme_().corrected())
573 {
574 // Calculate sn gradient
575 tmp<SurfaceFieldType> tfaceGrad
576 (
577 new SurfaceFieldType
578 (
580 (
581 "snGradCorr("+vf.name()+')',
582 vf.instance(),
583 mesh,
586 ),
587 mesh,
588 vf.dimensions()
589 )
590 );
591
592 {
593 // Calculate gradient
594 tmp<GradFieldType> tgGrad
595 (
597 (
598 mesh,
599 mesh.gradScheme("grad(" + vf.name() + ')')
600 )().grad(vf, "grad(" + vf.name() + ')')
601 );
602 const auto& gGrad = tgGrad();
603
604 // Doing a dotinterpolate with nonOrthCorrectionVectors
605 const auto dotInterpolate = [&]
606 (
607 const vector& area,
608 const scalar lambda,
609
610 const GradType& ownVal,
611 const GradType& neiVal,
612
613 const vector& dotVector,
614
615 Type& result
616 )
617 {
618 result = dotVector&(lambda*(ownVal - neiVal) + neiVal);
619 };
620
622 (
623 mesh.surfaceInterpolation::weights(), // linear interpolation
624 gGrad, // volume field
625 mesh.nonOrthCorrectionVectors(),// surface multiplier
626 dotInterpolate,
627 tfaceGrad.ref()
628 );
629 }
630 const auto& faceGrad = tfaceGrad();
631
632 const auto snGrad = [&]
633 (
634 const vector& Sf,
635 const scalar dc,
636 const Type& ownVal,
637 const Type& neiVal,
638 const Type& correction
639 ) -> Type
640 {
641 const auto snGrad(dc*(neiVal-ownVal) + correction);
642 return mag(Sf)*snGrad;
643 };
644
646 (
647 deltaCoeffs,
648 vf,
649 faceGrad, // face-based addition
650 snGrad,
651 result,
652 false // avoid boundary evaluation until volume division
653 );
654 }
655 else
656 {
657 const auto snGrad = [&]
658 (
659 const vector& Sf,
660 const scalar dc,
661 const Type& ownVal,
662 const Type& neiVal
663 ) -> Type
664 {
665 const auto snGrad(dc*(neiVal-ownVal));
666 return mag(Sf)*snGrad;
667 };
668
670 (
671 deltaCoeffs,
672 vf,
673 snGrad,
674 result,
675 false // avoid boundary evaluation until volume division
676 );
677 }
678
679 result.primitiveFieldRef() /= mesh.V();
680 result.correctBoundaryConditions();
681
682 return tresult;
683}
685
686template<>
688<
690>
692(
694)
695{
696 typedef vector Type;
697
700 typedef typename outerProduct<vector, Type>::type GradType;
702
703 const fvMesh& mesh = vf.mesh();
704
705 tmp<FieldType> tresult
706 (
707 new FieldType
708 (
710 (
711 "laplacian(" + vf.name() + ')',
712 vf.instance(),
713 mesh,
716 ),
717 mesh,
720 )
721 );
722 FieldType& result = tresult.ref();
723
725 << "fusedGaussLaplacianScheme<vector, GType>::fvcLaplacian on "
726 << vf.name()
727 << " to generate " << result.name() << endl;
728
729
730 const auto tdeltaCoeffs(this->tsnGradScheme_().deltaCoeffs(vf));
731 const auto& deltaCoeffs = tdeltaCoeffs();
732
733
734 if (this->tsnGradScheme_().corrected())
735 {
736 // Calculate sn gradient
737 tmp<SurfaceFieldType> tfaceGrad
738 (
739 new SurfaceFieldType
740 (
742 (
743 "snGradCorr("+vf.name()+')',
744 vf.instance(),
745 mesh,
748 ),
749 mesh,
750 vf.dimensions()
751 )
752 );
753
754 {
755 // Calculate gradient
756 tmp<GradFieldType> tgGrad
757 (
759 (
760 mesh,
761 mesh.gradScheme("grad(" + vf.name() + ')')
762 )().grad(vf, "grad(" + vf.name() + ')')
763 );
764 const auto& gGrad = tgGrad();
765
766 // Doing a dotinterpolate with nonOrthCorrectionVectors
767 const auto dotInterpolate = [&]
768 (
769 const vector& area,
770 const scalar lambda,
771
772 const GradType& ownVal,
773 const GradType& neiVal,
774
775 const vector& dotVector,
776
777 Type& result
778 )
779 {
780 result = dotVector&(lambda*(ownVal - neiVal) + neiVal);
781 };
782
784 (
785 mesh.surfaceInterpolation::weights(), // linear interpolation
786 gGrad, // volume field
787 mesh.nonOrthCorrectionVectors(),// surface multiplier
788 dotInterpolate,
789 tfaceGrad.ref()
790 );
791 }
792 const auto& faceGrad = tfaceGrad();
793
794 const auto snGrad = [&]
795 (
796 const vector& Sf,
797 const scalar dc,
798 const Type& ownVal,
799 const Type& neiVal,
800 const Type& correction
801 ) -> Type
802 {
803 const auto snGrad(dc*(neiVal-ownVal) + correction);
804 return mag(Sf)*snGrad;
805 };
806
808 (
809 deltaCoeffs,
810 vf,
811 faceGrad, // face-based addition
812 snGrad,
813 result,
814 false // avoid boundary evaluation until volume division
815 );
816 }
817 else
818 {
819 const auto snGrad = [&]
820 (
821 const vector& Sf,
822 const scalar dc,
823 const Type& ownVal,
824 const Type& neiVal
825 ) -> Type
826 {
827 const auto snGrad(dc*(neiVal-ownVal));
828 return mag(Sf)*snGrad;
829 };
830
832 (
833 deltaCoeffs,
834 vf,
835 snGrad,
836 result,
837 false // avoid boundary evaluation until volume division
838 );
839 }
840
841 result.primitiveFieldRef() /= mesh.V();
842 result.correctBoundaryConditions();
844 return tresult;
846
847
853
854
855// ************************************************************************* //
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
Generic GeometricField class.
@ 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
const fileName & instance() const noexcept
Read access to instance path component.
Definition IOobjectI.H:289
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcLaplacian(const GeometricField< Type, fvPatchField, volMesh > &)
static tmp< gradScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new gradScheme created on freestore.
Definition gradScheme.C:30
tmp< snGradScheme< Type > > tsnGradScheme_
tmp< surfaceInterpolationScheme< GType > > tinterpGammaScheme_
const fvMesh & mesh() const
Return mesh reference.
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank)>::type type
Definition products.H:118
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
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const scalar gamma
Definition EEqn.H:9
dynamicFvMesh & mesh
#define declareFvmLaplacianScalarGamma(Type)
Expression templates for fvMatrix.
#define makeFvLaplacianScheme(SS)
#define DebugPout
Report an information message using Foam::Pout.
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.
void surfaceSnSum(const surfaceScalarField &deltaCoeffs, const GeometricField< Type, fvPatchField, volMesh > &vf, const CellToFaceOp &cop, GeometricField< ResultType, fvPatchField, volMesh > &result, const bool doCorrectBoundaryConditions)
sum of snGrad
const dimensionSet dimArea(sqr(dimLength))
Tensor< scalar > tensor
Definition symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
Vector< scalar > vector
Definition vector.H:57
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars, i.e. SphericalTensor<scalar>.
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition symmTensor.H:55
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)