30 template<
class MatrixType>
39 template<
class MatrixType>
61 for (
label bi=
k; bi<m; bi++)
76 for (
label bi=
k+1; bi<m; bi++)
79 rSumSqrUk +=
sqr(uk[bi]);
81 rSumSqrUk = 2/rSumSqrUk;
84 for (
label bi=
k; bi<m; bi++)
86 for (
label bj=
k; bj<m; bj++)
88 Qk(bi, bj) = -rSumSqrUk*uk[bi]*uk[bj];
94 for (
label bi=
k; bi<m; bi++)
101 Rk(bi, bk) += Qk(bi, bj)*R_(bj, bk);
111 for (
label bi=
k; bi<m; bi++)
115 R_(bi, bj) = -Rk(bi, bj);
117 for (
label bj=
k; bj<m; bj++)
119 Q_(bi, bj) = -Q_(bi, bj);
126 for (
label bi=
k; bi<m; bi++)
130 R_(bi, bj) = Rk(bi, bj);
136 for (
label bi=0; bi<m; bi++)
138 for (
label bk=
k; bk<m; bk++)
141 for (
label bj=
k; bj<m; bj++)
143 Rk(bi, bk) += Q_(bi, bj)*Qk(bj, bk);
147 for (
label bi=0; bi<m; bi++)
151 Q_(bi, bj) = Rk(bi, bj);
158 template<
class MatrixType>
166 for (
int i=
n-1; i>=0; i--)
170 for (
label j=i+1; j<
n; j++)
172 sum -=
x[j]*R_(i, j);
175 if (
mag(R_(i, i)) < small)
178 <<
"Back-substitution failed due to small diagonal"
187 template<
class MatrixType>
194 const label m = Q_.m();
197 for (
label i=0; i<m; i++)
200 for (
label j=0; j<m; j++)
202 x[i] += Q_(j, i)*source[j];
210 template<
class MatrixType>
226 template<
class MatrixType>
235 for (
label i=0; i<m; i++)
237 for (
label j=0; j<m; j++)
242 inv.block(m, 1, 0, i) =
x;
Pre-declare SubField and related Field type.
label m() const
Return the number of rows.
Class templated on matrix type to perform the QR decomposition using Householder reflections on a squ...
void decompose(const MatrixType &M)
Decompose given matrix.
QRMatrix()
Construct null.
void solve(Field< cmptType > &x, const Field< cmptType > &source) const
Solve the linear system with the given source.
MatrixType::cmptType cmptType
QMatrixType inv() const
Return the inverse of a square matrix.
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
errorManip< error > abort(error &err)
static const Identity< scalar > I
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)