76 const Cmpt txx,
const Cmpt txy,
const Cmpt txz,
77 const Cmpt tyy,
const Cmpt tyz,
82 this->
v_[
YY] = tyy; this->
v_[
YZ] = tyz;
118template<Foam::direction Idx>
121 if (Idx == 0)
return x();
122 else if (Idx == 1)
return y();
123 else if (Idx == 2)
return z();
125 static_assert(Idx < 3,
"Invalid row access");
135 case 0:
return x();
break;
136 case 1:
return y();
break;
137 case 2:
return z();
break;
148template<Foam::direction Idx>
153 this->v_[XX] = v.
x(); this->v_[XY] = v.
y(); this->v_[XZ] = v.
z();
157 this->v_[XY] = v.
x(); this->v_[YY] = v.
y(); this->v_[YZ] = v.
z();
161 this->v_[XZ] = v.
x(); this->v_[YZ] = v.
y(); this->v_[ZZ] = v.
z();
164 static_assert(Idx < 3,
"Invalid row access");
171 const Vector<Cmpt>&
x,
172 const Vector<Cmpt>&
y,
173 const Vector<Cmpt>& z
178 this->
v_[
ZZ] =
z.z();
191 case 0: row<0>(v);
break;
192 case 1: row<1>(v);
break;
193 case 2: row<2>(v);
break;
261 xx()*yy()*zz() + xy()*yz()*xz()
275 return (yy()*zz() - yz()*yz());
280 return (xx()*zz() - xz()*xz());
295 yy()*zz() - yz()*yz(), xz()*yz() - xy()*zz(), xy()*yz() - xz()*yy(),
330 return SymmTensor<Cmpt>
355 const Cmpt detval = this->det2D(excludeCmpt);
357 return this->
adjunct2D(excludeCmpt)/detval;
365 const Cmpt detval = this->
det();
368 if (
mag(detval) < VSMALL)
371 <<
"Tensor not properly invertible, determinant:"
372 << detval <<
" tensor:" << *
this <<
nl
394 const scalar threshold = SMALL * (magSqr_xx + magSqr_yy + magSqr_zz);
396 const bool small_xx = (magSqr_xx < threshold);
397 const bool small_yy = (magSqr_yy < threshold);
398 const bool small_zz = (magSqr_zz < threshold);
400 if (small_xx || small_yy || small_zz)
408 const Cmpt detval = work.
det();
410 if (
mag(detval) < ROOTVSMALL)
426 const Cmpt detval = this->
det();
428 if (
mag(detval) < ROOTVSMALL)
431 return SymmTensor<Cmpt>(
Zero);
443 this->v_[XX] = st.
ii(); this->v_[XY] =
Zero; this->v_[XZ] =
Zero;
444 this->v_[YY] = st.
ii(); this->v_[YZ] =
Zero;
445 this->v_[ZZ] = st.
ii();
460 return st.
xx() + st.
yy() + st.
zz();
485inline SymmTensor<Cmpt>
twoSymm(
const SymmTensor<Cmpt>& st)
493inline SymmTensor<Cmpt>
devSymm(
const SymmTensor<Cmpt>& st)
501inline SymmTensor<Cmpt>
devTwoSymm(
const SymmTensor<Cmpt>& st)
519 return st - 2*
sph(st);
525inline Cmpt
det(
const SymmTensor<Cmpt>& st)
533inline SymmTensor<Cmpt>
cof(
const SymmTensor<Cmpt>& st)
541inline SymmTensor<Cmpt>
inv(
const SymmTensor<Cmpt>& st,
const Cmpt detval)
544 if (
mag(detval) < VSMALL)
547 <<
"SymmTensor not properly invertible, determinant:"
548 << detval <<
" tensor:" << st <<
nl
553 return st.adjunct()/detval;
567inline Cmpt
invariantI(
const SymmTensor<Cmpt>& st)
579 st.xx()*st.yy() + st.yy()*st.zz() + st.xx()*st.zz()
580 - st.xy()*st.xy() - st.yz()*st.yz() - st.xz()*st.xz()
595inline SymmTensor<Cmpt>
627inline SymmTensor<Cmpt>
sqr(
const Vector<Cmpt>& v)
629 return SymmTensor<Cmpt>
631 v.x()*v.x(), v.x()*v.y(), v.x()*v.z(),
632 v.y()*v.y(), v.y()*v.z(),
647 const scalar onet = (1-t);
651 onet*a.
xx() + t*
b.xx(),
652 onet*a.
xy() + t*
b.xy(),
653 onet*a.
xz() + t*
b.xz(),
654 onet*a.
yy() + t*
b.yy(),
655 onet*a.
yz() + t*
b.yz(),
665inline SymmTensor<Cmpt>
666operator+(
const SphericalTensor<Cmpt>& spt1,
const SymmTensor<Cmpt>& st2)
668 return SymmTensor<Cmpt>
670 spt1.ii() + st2.xx(), st2.xy(), st2.xz(),
671 spt1.ii() + st2.yy(), st2.yz(),
679inline SymmTensor<Cmpt>
680operator+(
const SymmTensor<Cmpt>& st1,
const SphericalTensor<Cmpt>& spt2)
682 return SymmTensor<Cmpt>
684 st1.xx() + spt2.ii(), st1.xy(), st1.xz(),
685 st1.yy() + spt2.ii(), st1.yz(),
698 spt1.ii() - st2.xx(), -st2.xy(), -st2.xz(),
699 spt1.ii() - st2.yy(), -st2.yz(),
707inline SymmTensor<Cmpt>
708operator-(
const SymmTensor<Cmpt>& st1,
const SphericalTensor<Cmpt>& spt2)
710 return SymmTensor<Cmpt>
712 st1.xx() - spt2.ii(), st1.xy(), st1.xz(),
713 st1.yy() - spt2.ii(), st1.yz(),
721inline Vector<Cmpt>
operator*(
const SymmTensor<Cmpt>& st)
723 return Vector<Cmpt>(st.yz(), -st.xz(), st.xy());
734 st.xx()/
s, st.xy()/
s, st.xz()/
s,
735 st.yy()/
s, st.yz()/
s,
744operator&(
const SymmTensor<Cmpt>& st1,
const SymmTensor<Cmpt>& st2)
748 st1.xx()*st2.xx() + st1.xy()*st2.xy() + st1.xz()*st2.xz(),
749 st1.xx()*st2.xy() + st1.xy()*st2.yy() + st1.xz()*st2.yz(),
750 st1.xx()*st2.xz() + st1.xy()*st2.yz() + st1.xz()*st2.zz(),
752 st1.xy()*st2.xx() + st1.yy()*st2.xy() + st1.yz()*st2.xz(),
753 st1.xy()*st2.xy() + st1.yy()*st2.yy() + st1.yz()*st2.yz(),
754 st1.xy()*st2.xz() + st1.yy()*st2.yz() + st1.yz()*st2.zz(),
756 st1.xz()*st2.xx() + st1.yz()*st2.xy() + st1.zz()*st2.xz(),
757 st1.xz()*st2.xy() + st1.yz()*st2.yy() + st1.zz()*st2.yz(),
758 st1.xz()*st2.xz() + st1.yz()*st2.yz() + st1.zz()*st2.zz()
770 spt1.
ii()*st2.
xx(), spt1.
ii()*st2.
xy(), spt1.
ii()*st2.
xz(),
771 spt1.
ii()*st2.
yy(), spt1.
ii()*st2.
yz(),
779inline SymmTensor<Cmpt>
780operator&(
const SymmTensor<Cmpt>& st1,
const SphericalTensor<Cmpt>& spt2)
782 return SymmTensor<Cmpt>
784 st1.xx()*spt2.ii(), st1.xy()*spt2.ii(), st1.xz()*spt2.ii(),
785 st1.yy()*spt2.ii(), st1.yz()*spt2.ii(),
798 st.
xx()*v.
x() + st.
xy()*v.
y() + st.
xz()*v.
z(),
799 st.
xy()*v.
x() + st.
yy()*v.
y() + st.
yz()*v.
z(),
800 st.
xz()*v.
x() + st.
yz()*v.
y() + st.
zz()*v.
z()
808operator&(
const Vector<Cmpt>& v,
const SymmTensor<Cmpt>& st)
812 v.x()*st.xx() + v.y()*st.xy() + v.z()*st.xz(),
813 v.x()*st.xy() + v.y()*st.yy() + v.z()*st.yz(),
814 v.x()*st.xz() + v.y()*st.yz() + v.z()*st.zz()
822operator&&(
const SymmTensor<Cmpt>& st1,
const SymmTensor<Cmpt>& st2)
826 st1.xx()*st2.xx() + 2*st1.xy()*st2.xy() + 2*st1.xz()*st2.xz()
827 + st1.yy()*st2.yy() + 2*st1.yz()*st2.yz()
836operator&&(
const SphericalTensor<Cmpt>& spt1,
const SymmTensor<Cmpt>& st2)
838 return (spt1.ii()*st2.xx() + spt1.ii()*st2.yy() + spt1.ii()*st2.zz());
847 return (st1.
xx()*spt2.
ii() + st1.
yy()*spt2.
ii() + st1.
zz()*spt2.
ii());
854class outerProduct<SymmTensor<Cmpt>, Cmpt>
858 typedef SymmTensor<Cmpt>
type;
870class innerProduct<SymmTensor<Cmpt>, SymmTensor<Cmpt>>
874 typedef Tensor<Cmpt>
type;
886class innerProduct<Vector<Cmpt>, SymmTensor<Cmpt>>
890 typedef Vector<Cmpt>
type;
903class typeOfSum<SymmTensor<Cmpt>, SphericalTensor<Cmpt>>
911class innerProduct<SphericalTensor<Cmpt>, SymmTensor<Cmpt>>
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
A templated (3 x 3) diagonal tensor of objects of <T>, effectively containing 1 element,...
const Cmpt & ii() const noexcept
A templated (3 x 3) symmetric tensor of objects of <T>, effectively containing 6 elements,...
scalar diagSqr() const
The L2-norm squared of the diagonal.
bool is_identity(const scalar tol=ROOTVSMALL) const
Is identity tensor?
SymmTensor< Cmpt > adjunct() const
Return adjunct matrix (transpose of cofactor matrix).
SymmTensor< Cmpt > inv2D(const direction excludeCmpt) const
Return inverse of 2D tensor (by excluding given direction).
const Cmpt & yy() const noexcept
void addDiag(const Vector< Cmpt > &v)
Add to the diagonal.
void rows(const Vector< Cmpt > &x, const Vector< Cmpt > &y, const Vector< Cmpt > &z)
Set row values.
Cmpt det2D(const direction excludeCmpt) const
The 2D determinant by excluding given direction.
Vector< Cmpt > row() const
Extract vector for given row: compile-time check of index.
void subtractDiag(const Vector< Cmpt > &v)
Subtract from the diagonal.
SymmTensor< Cmpt > cof() const
Return cofactor matrix (transpose of adjunct matrix).
Vector< Cmpt > z() const
Extract vector for row 2.
SymmTensor< Cmpt > adjunct2D(const direction excludeCmpt) const
Return 2D adjunct matrix by excluding given direction.
Cmpt det() const
The determinate.
Vector< Cmpt > y() const
Extract vector for row 1.
const Cmpt & xx() const noexcept
SymmTensor< Cmpt > safeInv() const
Return inverse, with (ad hoc) failsafe handling of 2D tensors.
const Cmpt & yz() const noexcept
Vector< Cmpt > diag() const
Extract the diagonal as a vector.
SymmTensor< Cmpt > inv() const
Return inverse.
const Cmpt & zz() const noexcept
const Cmpt & xy() const noexcept
SymmTensor & operator=(const SymmTensor &)=default
Copy assignment.
SymmTensor()=default
Default construct.
const Cmpt & xz() const noexcept
Vector< Cmpt > x() const
Extract vector for row 0.
A templated (3 x 3) tensor of objects of <T> derived from MatrixSpace.
static const SymmTensor< Cmpt > zero
VectorSpace< SymmTensor< Cmpt >, Cmpt, Ncmpts > vsType
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
const Cmpt & x() const noexcept
Access to the vector x component.
const Cmpt & z() const noexcept
Access to the vector z component.
const Cmpt & y() const noexcept
Access to the vector y component.
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) -2 >::type type
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank)>::type type
A traits class, which is primarily used for primitives and vector-space.
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
tmp< faMatrix< Type > > operator-(const faMatrix< Type > &)
Unary negation.
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar det(const dimensionedSphericalTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< faMatrix< Type > > operator+(const faMatrix< Type > &, const faMatrix< Type > &)
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Cmpt invariantII(const SymmTensor2D< Cmpt > &st)
Return the 2nd invariant of a SymmTensor2D.
Cmpt invariantI(const SymmTensor2D< Cmpt > &st)
Return the 1st invariant of a SymmTensor2D.
tmp< faMatrix< Type > > operator*(const areaScalarField::Internal &, const faMatrix< Type > &)
SymmTensor< Cmpt > devSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of the symmetric part of a SymmTensor.
Cmpt invariantIII(const SymmTensor< Cmpt > &st)
Return the 3rd invariant of a SymmTensor.
dimensionedScalar operator/(const scalar s1, const dimensionedScalar &ds2)
dimensionedSymmTensor cof(const dimensionedSymmTensor &dt)
tmp< GeometricField< Type, faPatchField, areaMesh > > operator&(const faMatrix< Type > &, const DimensionedField< Type, areaMesh > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
SphericalTensor< Cmpt > sph(const DiagTensor< Cmpt > &dt)
Return the spherical part of a DiagTensor as a SphericalTensor.
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
dimensioned< typename scalarProduct< Type1, Type2 >::type > operator&&(const dimensioned< Type1 > &, const dimensioned< Type2 > &)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
constexpr char nl
The newline '\n' character (0x0a).