37 #ifndef VIGRA_LINEAR_SOLVE_HXX 
   38 #define VIGRA_LINEAR_SOLVE_HXX 
   42 #include "mathutil.hxx" 
   44 #include "singular_value_decomposition.hxx" 
   55 template <
class T, 
class C1>
 
   56 T determinantByLUDecomposition(MultiArrayView<2, T, C1> 
const & a)
 
   61     vigra_precondition(n == m,
 
   62         "determinantByLUDecomposition(): square matrix required.");
 
   63     vigra_precondition(NumericTraits<T>::isIntegral::value == 
false,
 
   64         "determinantByLUDecomposition(): Input matrix must not be an integral type.");
 
   76             LU(i,j) = LU(i,j) -= s;
 
   98 template <
class T, 
class C1>
 
   99 typename NumericTraits<T>::Promote
 
  100 determinantByMinors(MultiArrayView<2, T, C1> 
const & mat)
 
  102     typedef typename NumericTraits<T>::Promote PromoteType;
 
  107             "determinantByMinors(): square matrix required.");
 
  109             NumericTraits<PromoteType>::isSigned::value,
 
  110             "determinantByMinors(): promote type must be signed.");
 
  117         Matrix<T> minor_mat(Shape2(m-1, n-1));
 
  118         PromoteType det = NumericTraits<PromoteType>::zero();
 
  129             const PromoteType 
sign = 1 - 2 * (i % 2);
 
  130             det += sign * mat(i, 0) * determinantByMinors(minor_mat);
 
  139 T givensCoefficients(T a, T b, T & c, T & s)
 
  167 bool givensRotationMatrix(T a, T b, Matrix<T> & gTranspose)
 
  171     givensCoefficients(a, b, gTranspose(0,0), gTranspose(0,1));
 
  172     gTranspose(1,1) = gTranspose(0,0);
 
  173     gTranspose(1,0) = -gTranspose(0,1);
 
  181 givensReflectionMatrix(T a, T b, Matrix<T> & g)
 
  185     givensCoefficients(a, b, g(0,0), g(0,1));
 
  192 template <
class T, 
class C1, 
class C2>
 
  194 qrGivensStepImpl(
MultiArrayIndex i, MultiArrayView<2, T, C1> r, MultiArrayView<2, T, C2> rhs)
 
  196     typedef typename Matrix<T>::difference_type Shape;
 
  201     vigra_precondition(m == 
rowCount(rhs),
 
  202                        "qrGivensStepImpl(): Matrix shape mismatch.");
 
  204     Matrix<T> givens(2,2);
 
  205     for(
int k=m-1; k>
static_cast<int>(i); --k)
 
  207         if(!givensReflectionMatrix(r(k-1,i), r(k,i), givens))
 
  210         r(k-1,i) = givens(0,0)*r(k-1,i) + givens(0,1)*r(k,i);
 
  213         r.subarray(Shape(k-1,i+1), Shape(k+1,n)) = givens*r.subarray(Shape(k-1,i+1), Shape(k+1,n));
 
  214         rhs.subarray(Shape(k-1,0), Shape(k+1,rhsCount)) = givens*rhs.subarray(Shape(k-1,0), Shape(k+1,rhsCount));
 
  216     return r(i,i) != 0.0;
 
  220 template <
class T, 
class C1, 
class C2, 
class Permutation>
 
  223                                   MultiArrayView<2, T, C1> &r, MultiArrayView<2, T, C2> &rhs, Permutation & permutation)
 
  225     typedef typename Matrix<T>::difference_type Shape;
 
  230     vigra_precondition(i < n && j < n,
 
  231                        "upperTriangularCyclicShiftColumns(): Shift indices out of range.");
 
  232     vigra_precondition(m == 
rowCount(rhs),
 
  233                        "upperTriangularCyclicShiftColumns(): Matrix shape mismatch.");
 
  245         permutation[k] = permutation[k+1];
 
  250     Matrix<T> givens(2,2);
 
  253         if(!givensReflectionMatrix(r(k,k), r(k+1,k), givens))
 
  256         r(k,k) = givens(0,0)*r(k,k) + givens(0,1)*r(k+1,k);
 
  259         r.subarray(Shape(k,k+1), Shape(k+2,n)) = givens*r.subarray(Shape(k,k+1), Shape(k+2,n));
 
  260         rhs.subarray(Shape(k,0), Shape(k+2,rhsCount)) = givens*rhs.subarray(Shape(k,0), Shape(k+2,rhsCount));
 
  265 template <
class T, 
class C1, 
class C2, 
class Permutation>
 
  268                            MultiArrayView<2, T, C1> &r, MultiArrayView<2, T, C2> &rhs, Permutation & permutation)
 
  270     typedef typename Matrix<T>::difference_type Shape;
 
  275     vigra_precondition(i < n && j < n,
 
  276                        "upperTriangularSwapColumns(): Swap indices out of range.");
 
  277     vigra_precondition(m == 
rowCount(rhs),
 
  278                        "upperTriangularSwapColumns(): Matrix shape mismatch.");
 
  286     std::swap(permutation[i], permutation[j]);
 
  288     Matrix<T> givens(2,2);
 
  289     for(
int k=m-1; k>
static_cast<int>(i); --k)
 
  291         if(!givensReflectionMatrix(r(k-1,i), r(k,i), givens))
 
  294         r(k-1,i) = givens(0,0)*r(k-1,i) + givens(0,1)*r(k,i);
 
  297         r.subarray(Shape(k-1,i+1), Shape(k+1,n)) = givens*r.subarray(Shape(k-1,i+1), Shape(k+1,n));
 
  298         rhs.subarray(Shape(k-1,0), Shape(k+1,rhsCount)) = givens*rhs.subarray(Shape(k-1,0), Shape(k+1,rhsCount));
 
  303         if(!givensReflectionMatrix(r(k,k), r(k+1,k), givens))
 
  306         r(k,k) = givens(0,0)*r(k,k) + givens(0,1)*r(k+1,k);
 
  309         r.subarray(Shape(k,k+1), Shape(k+2,n)) = givens*r.subarray(Shape(k,k+1), Shape(k+2,n));
 
  310         rhs.subarray(Shape(k,0), Shape(k+2,rhsCount)) = givens*rhs.subarray(Shape(k,0), Shape(k+2,rhsCount));
 
  315 template <
class T, 
class C1, 
class C2, 
class U>
 
  316 bool householderVector(MultiArrayView<2, T, C1> 
const & v, MultiArrayView<2, T, C2> & u, U & vnorm)
 
  318     vnorm = (v(0,0) > 0.0)
 
  323     if(f == NumericTraits<U>::zero())
 
  325         u.init(NumericTraits<T>::zero());
 
  330         u(0,0) = (v(0,0) - vnorm) / f;
 
  338 template <
class T, 
class C1, 
class C2, 
class C3>
 
  341                       MultiArrayView<2, T, C2> & rhs, MultiArrayView<2, T, C3> & householderMatrix)
 
  343     typedef typename Matrix<T>::difference_type Shape;
 
  349     vigra_precondition(i < n && i < m,
 
  350         "qrHouseholderStepImpl(): Index i out of range.");
 
  354     bool nontrivial = householderVector(
columnVector(r, Shape(i,i), m), u, vnorm);
 
  357     columnVector(r, Shape(i+1,i), m).init(NumericTraits<T>::zero());
 
  369     return r(i,i) != 0.0;
 
  372 template <
class T, 
class C1, 
class C2>
 
  374 qrColumnHouseholderStep(
MultiArrayIndex i, MultiArrayView<2, T, C1> &r, MultiArrayView<2, T, C2> &rhs)
 
  376     Matrix<T> dontStoreHouseholderVectors; 
 
  377     return qrHouseholderStepImpl(i, r, rhs, dontStoreHouseholderVectors);
 
  380 template <
class T, 
class C1, 
class C2>
 
  382 qrRowHouseholderStep(
MultiArrayIndex i, MultiArrayView<2, T, C1> &r, MultiArrayView<2, T, C2> & householderMatrix)
 
  384     Matrix<T> dontTransformRHS; 
 
  385     MultiArrayView<2, T, StridedArrayTag> rt = 
transpose(r),
 
  387     return qrHouseholderStepImpl(i, rt, dontTransformRHS, ht);
 
  391 template <
class T, 
class C1, 
class C2, 
class SNType>
 
  393 incrementalMaxSingularValueApproximation(MultiArrayView<2, T, C1> 
const & newColumn, 
 
  394                                          MultiArrayView<2, T, C2> & z, SNType & v) 
 
  396     typedef typename Matrix<T>::difference_type Shape;
 
  407     z(n,0) = s*newColumn(n,0);
 
  411 template <
class T, 
class C1, 
class C2, 
class SNType>
 
  413 incrementalMinSingularValueApproximation(MultiArrayView<2, T, C1> 
const & newColumn, 
 
  414                                          MultiArrayView<2, T, C2> & z, SNType & v, 
double tolerance) 
 
  416     typedef typename Matrix<T>::difference_type Shape;
 
  426     T 
gamma = newColumn(n,0);
 
  440     z(n,0) = (s - c*yv) / gamma;
 
  441     v *= 
norm(gamma) / 
hypot(c*gamma, v*(s - c*yv));
 
  445 template <
class T, 
class C1, 
class C2, 
class C3>
 
  447 qrTransformToTriangularImpl(MultiArrayView<2, T, C1> & r, MultiArrayView<2, T, C2> & rhs, MultiArrayView<2, T, C3> & householder,
 
  448                             ArrayVector<MultiArrayIndex> & permutation, 
double epsilon)
 
  450     typedef typename Matrix<T>::difference_type Shape;
 
  451     typedef typename NormTraits<MultiArrayView<2, T, C1> >::NormType NormType;
 
  452     typedef typename NormTraits<MultiArrayView<2, T, C1> >::SquaredNormType SNType;
 
  458     vigra_precondition(m >= n,
 
  459         "qrTransformToTriangularImpl(): Coefficient matrix with at least as many rows as columns required.");
 
  462     bool transformRHS = rhsCount > 0;
 
  463     vigra_precondition(!transformRHS || m == 
rowCount(rhs),
 
  464                        "qrTransformToTriangularImpl(): RHS matrix shape mismatch.");
 
  466     bool storeHouseholderSteps = 
columnCount(householder) > 0;
 
  467     vigra_precondition(!storeHouseholderSteps || r.shape() == householder.shape(),
 
  468                        "qrTransformToTriangularImpl(): Householder matrix shape mismatch.");
 
  470     bool pivoting = permutation.size() > 0;
 
  471     vigra_precondition(!pivoting || n == static_cast<MultiArrayIndex>(permutation.size()),
 
  472                        "qrTransformToTriangularImpl(): Permutation array size mismatch.");
 
  477     Matrix<SNType> columnSquaredNorms;
 
  480         columnSquaredNorms.reshape(Shape(1,n));
 
  484         int pivot = 
argMax(columnSquaredNorms);
 
  488             std::swap(columnSquaredNorms[0], columnSquaredNorms[pivot]);
 
  489             std::swap(permutation[0], permutation[pivot]);
 
  493     qrHouseholderStepImpl(0, r, rhs, householder);
 
  496     NormType maxApproxSingularValue = 
norm(r(0,0)),
 
  497              minApproxSingularValue = maxApproxSingularValue;
 
  499     double tolerance = (epsilon == 0.0)
 
  500                           ? m*maxApproxSingularValue*NumericTraits<T>::epsilon()
 
  503     bool simpleSingularValueApproximation = (n < 4);
 
  504     Matrix<T> zmax, zmin;
 
  505     if(minApproxSingularValue <= tolerance)
 
  509         simpleSingularValueApproximation = 
true;
 
  511     if(!simpleSingularValueApproximation)
 
  513         zmax.reshape(Shape(m,1));
 
  514         zmin.reshape(Shape(m,1));
 
  516         zmin(0,0) = 1.0 / r(0,0);
 
  526             if(pivot != static_cast<int>(k))
 
  529                 std::swap(columnSquaredNorms[k], columnSquaredNorms[pivot]);
 
  530                 std::swap(permutation[k], permutation[pivot]);
 
  534         qrHouseholderStepImpl(k, r, rhs, householder);
 
  536         if(simpleSingularValueApproximation)
 
  538             NormType nv = 
norm(r(k,k));        
 
  539             maxApproxSingularValue = std::max(nv, maxApproxSingularValue);
 
  540             minApproxSingularValue = std::min(nv, minApproxSingularValue);
 
  544             incrementalMaxSingularValueApproximation(
columnVector(r, Shape(0,k),k+1), zmax, maxApproxSingularValue);
 
  545             incrementalMinSingularValueApproximation(
columnVector(r, Shape(0,k),k+1), zmin, minApproxSingularValue, tolerance);
 
  549         Matrix<T> u(k+1,k+1), s(k+1, 1), v(k+1,k+1);
 
  551         std::cerr << 
"estimate, svd " << k << 
": " << minApproxSingularValue << 
" " << s(k,0) << 
"\n";
 
  555             tolerance = m*maxApproxSingularValue*NumericTraits<T>::epsilon();
 
  557         if(minApproxSingularValue > tolerance)
 
  562     return static_cast<unsigned int>(rank);
 
  565 template <
class T, 
class C1, 
class C2>
 
  567 qrTransformToUpperTriangular(MultiArrayView<2, T, C1> & r, MultiArrayView<2, T, C2> & rhs, 
 
  568                              ArrayVector<MultiArrayIndex> & permutation, 
double epsilon = 0.0)
 
  570     Matrix<T> dontStoreHouseholderVectors; 
 
  571     return qrTransformToTriangularImpl(r, rhs, dontStoreHouseholderVectors, permutation, epsilon);
 
  575 template <
class T, 
class C1, 
class C2, 
class C3>
 
  577 qrTransformToLowerTriangular(MultiArrayView<2, T, C1> & r, MultiArrayView<2, T, C2> & rhs, MultiArrayView<2, T, C3> & householderMatrix, 
 
  578                       double epsilon = 0.0)
 
  580     ArrayVector<MultiArrayIndex> permutation(static_cast<unsigned int>(
rowCount(rhs)));
 
  581     for(
MultiArrayIndex k=0; k<static_cast<MultiArrayIndex>(permutation.size()); ++k)
 
  583     Matrix<T> dontTransformRHS; 
 
  584     MultiArrayView<2, T, StridedArrayTag> rt = 
transpose(r),
 
  586     unsigned int rank = qrTransformToTriangularImpl(rt, dontTransformRHS, ht, permutation, epsilon);
 
  589     Matrix<T> tempRHS(rhs);
 
  590     for(
MultiArrayIndex k=0; k<static_cast<MultiArrayIndex>(permutation.size()); ++k)
 
  596 template <
class T, 
class C1, 
class C2>
 
  598 qrTransformToUpperTriangular(MultiArrayView<2, T, C1> & r, MultiArrayView<2, T, C2> & rhs, 
 
  599                       double epsilon = 0.0)
 
  601     ArrayVector<MultiArrayIndex> noPivoting; 
 
  603     return (qrTransformToUpperTriangular(r, rhs, noPivoting, epsilon) == 
 
  608 template <
class T, 
class C1, 
class C2>
 
  610 qrTransformToLowerTriangular(MultiArrayView<2, T, C1> & r, MultiArrayView<2, T, C2> & householder, 
 
  611                       double epsilon = 0.0)
 
  613     Matrix<T> noPivoting; 
 
  615     return (qrTransformToLowerTriangular(r, noPivoting, householder, epsilon) == 
 
  616            static_cast<unsigned int>(
rowCount(r)));
 
  620 template <
class T, 
class C1, 
class C2, 
class Permutation>
 
  621 void inverseRowPermutation(MultiArrayView<2, T, C1> &permuted, MultiArrayView<2, T, C2> &res,
 
  622                            Permutation 
const & permutation)
 
  626             res(permutation[l], k) = permuted(l,k);
 
  629 template <
class T, 
class C1, 
class C2>
 
  630 void applyHouseholderColumnReflections(MultiArrayView<2, T, C1> 
const &householder, MultiArrayView<2, T, C2> &res)
 
  632     typedef typename Matrix<T>::difference_type Shape;
 
  637     for(
int k = m-1; k >= 0; --k)
 
  639         MultiArrayView<2, T, C1> u = 
columnVector(householder, Shape(k,k), n);
 
  647 template <
class T, 
class C1, 
class C2, 
class C3>
 
  649 linearSolveQRReplace(MultiArrayView<2, T, C1> &A, MultiArrayView<2, T, C2> &b,
 
  650                      MultiArrayView<2, T, C3> & res, 
 
  651                      double epsilon = 0.0)
 
  653     typedef typename Matrix<T>::difference_type Shape;
 
  663            "linearSolveQR(): RHS and solution must have the same number of columns.");
 
  664     vigra_precondition(m == 
rowCount(b),
 
  665            "linearSolveQR(): Coefficient matrix and RHS must have the same number of rows.");
 
  666     vigra_precondition(n == 
rowCount(res),
 
  667            "linearSolveQR(): Mismatch between column count of coefficient matrix and row count of solution.");
 
  668     vigra_precondition(epsilon >= 0.0,
 
  669            "linearSolveQR(): 'epsilon' must be non-negative.");
 
  674         Matrix<T> householderMatrix(n, m);
 
  675         MultiArrayView<2, T, StridedArrayTag> ht = 
transpose(householderMatrix);
 
  676         rank = 
static_cast<MultiArrayIndex>(detail::qrTransformToLowerTriangular(A, b, ht, epsilon));
 
  677         res.subarray(Shape(rank,0), Shape(n, rhsCount)).init(NumericTraits<T>::zero());
 
  681             MultiArrayView<2, T, C1> Asub = A.subarray(ul, Shape(m,rank));
 
  682             detail::qrTransformToUpperTriangular(Asub, b, epsilon);
 
  684                                        b.subarray(ul, Shape(rank,rhsCount)), 
 
  685                                        res.subarray(ul, Shape(rank, rhsCount)));
 
  691                                        b.subarray(ul, Shape(rank, rhsCount)), 
 
  692                                        res.subarray(ul, Shape(rank, rhsCount)));
 
  694         detail::applyHouseholderColumnReflections(householderMatrix.subarray(ul, Shape(n, rank)), res);
 
  699         ArrayVector<MultiArrayIndex> permutation(static_cast<unsigned int>(n));
 
  703         rank = 
static_cast<MultiArrayIndex>(detail::qrTransformToUpperTriangular(A, b, permutation, epsilon));
 
  705         Matrix<T> permutedSolution(n, rhsCount);
 
  709             Matrix<T> householderMatrix(n, rank);
 
  710             MultiArrayView<2, T, StridedArrayTag> ht = 
transpose(householderMatrix);
 
  711             MultiArrayView<2, T, C1> Asub = A.subarray(ul, Shape(rank,n));
 
  712             detail::qrTransformToLowerTriangular(Asub, ht, epsilon);
 
  714                                        b.subarray(ul, Shape(rank, rhsCount)), 
 
  715                                        permutedSolution.subarray(ul, Shape(rank, rhsCount)));
 
  716             detail::applyHouseholderColumnReflections(householderMatrix, permutedSolution);
 
  722                                        b.subarray(ul, Shape(rank,rhsCount)), 
 
  725         detail::inverseRowPermutation(permutedSolution, res, permutation);
 
  727     return static_cast<unsigned int>(rank);
 
  730 template <
class T, 
class C1, 
class C2, 
class C3>
 
  731 unsigned int linearSolveQR(MultiArrayView<2, T, C1> 
const & A, MultiArrayView<2, T, C2> 
const & b,
 
  732                                   MultiArrayView<2, T, C3> & res)
 
  734     Matrix<T> r(A), rhs(b);
 
  735     return linearSolveQRReplace(r, rhs, res);
 
  759 template <
class T, 
class C1, 
class C2>
 
  767        "inverse(): shape of output matrix must be the transpose of the input matrix' shape.");
 
  776                                    transpose(q).subarray(Shape(0,0), Shape(m,n)), 
 
  785                                    transpose(q).subarray(Shape(0,0), Shape(n,m)), 
 
  813 template <
class T, 
class C>
 
  817     vigra_precondition(
inverse(v, ret),
 
  818         "inverse(): matrix is not invertible.");
 
  823 TemporaryMatrix<T> 
inverse(
const TemporaryMatrix<T> &v)
 
  827         vigra_precondition(
inverse(v, 
const_cast<TemporaryMatrix<T> &
>(v)),
 
  828             "inverse(): matrix is not invertible.");
 
  834         vigra_precondition(
inverse(v, ret),
 
  835             "inverse(): matrix is not invertible.");
 
  857 template <
class T, 
class C1>
 
  858 typename NumericTraits<T>::Promote
 
  861     typedef typename NumericTraits<T>::Promote PromoteType;
 
  863     vigra_precondition(
rowCount(a) == n,
 
  864                "determinant(): Square matrix required.");    
 
  868     if(method == 
"default")
 
  870         method = NumericTraits<T>::isIntegral::value ? 
"minor" : 
"lu";
 
  875         return a(0,0)*a(1,1) - a(0,1)*a(1,0);
 
  878         return detail::determinantByLUDecomposition(a);
 
  880     else if(method == 
"minor")
 
  882         return detail::determinantByMinors(a);
 
  884     else if(method == 
"cholesky")
 
  888            "determinant(): Cholesky method requires symmetric positive definite matrix.");
 
  896         vigra_precondition(
false, 
"determinant(): Unknown solution method.");
 
  898     return PromoteType();
 
  910 template <
class T, 
class C1>
 
  914     vigra_precondition(
rowCount(a) == n,
 
  915                "logDeterminant(): Square matrix required.");    
 
  918         vigra_precondition(a(0,0) > 0.0,
 
  919                    "logDeterminant(): Matrix not positive definite.");    
 
  924         T det = a(0,0)*a(1,1) - a(0,1)*a(1,0);
 
  925         vigra_precondition(det > 0.0,
 
  926                    "logDeterminant(): Matrix not positive definite.");    
 
  933                 "logDeterminant(): Matrix not positive definite.");  
 
  958 template <
class T, 
class C1, 
class C2>
 
  963     vigra_precondition(NumericTraits<T>::isIntegral::value == 
false,
 
  964                        "choleskyDecomposition(): Input matrix must not be an integral type.");
 
  965     vigra_precondition(
rowCount(A) == n,
 
  966                        "choleskyDecomposition(): Input matrix must be square.");
 
  968                        "choleskyDecomposition(): Output matrix must have same shape as input matrix.");
 
  970                        "choleskyDecomposition(): Input matrix must be symmetric.");
 
  980                s += L(k, i)*L(j, i);
 
  982             L(j, k) = s = (A(j, k) - s)/L(k, k);
 
 1015 template <
class T, 
class C1, 
class C2, 
class C3>
 
 1018                      double epsilon = 0.0)
 
 1024                        "qrDecomposition(): Matrix shape mismatch.");
 
 1026     q = identityMatrix<T>(m);
 
 1030     return (static_cast<MultiArrayIndex>(detail::qrTransformToUpperTriangular(r, tq, noPivoting, epsilon) == std::min(m,n)));
 
 1035 template <
class T, 
class C1, 
class C2, 
class C3>
 
 1067 template <
class T, 
class C1, 
class C2, 
class C3>
 
 1074         "linearSolveUpperTriangular(): square coefficient matrix required.");
 
 1076         "linearSolveUpperTriangular(): matrix shape mismatch.");
 
 1080         for(
int i=m-1; i>=0; --i)
 
 1082             if(r(i,i) == NumericTraits<T>::zero())
 
 1086                  sum -= r(i, j) * x(j, k);
 
 1087             x(i, k) = sum / r(i, i);
 
 1117 template <
class T, 
class C1, 
class C2, 
class C3>
 
 1123     vigra_precondition(m == 
rowCount(l),
 
 1124         "linearSolveLowerTriangular(): square coefficient matrix required.");
 
 1126         "linearSolveLowerTriangular(): matrix shape mismatch.");
 
 1132             if(l(i,i) == NumericTraits<T>::zero())
 
 1136                  sum -= l(i, j) * x(j, k);
 
 1137             x(i, k) = sum / l(i, i);
 
 1166 template <
class T, 
class C1, 
class C2, 
class C3>
 
 1247 template <
class T, 
class C1, 
class C2, 
class C3>
 
 1251                  std::string method = 
"QR")
 
 1256     vigra_precondition(n <= m,
 
 1257         "linearSolve(): Coefficient matrix A must have at least as many rows as columns.");
 
 1258     vigra_precondition(n == 
rowCount(res) && 
 
 1260         "linearSolve(): matrix shape mismatch.");
 
 1263     if(method == 
"cholesky")
 
 1266             "linearSolve(): Cholesky method requires square coefficient matrix.");
 
 1267         Matrix<T> L(A.
shape());
 
 1272     else if(method == 
"qr")
 
 1276     else if(method == 
"ne")
 
 1280     else if(method == 
"svd")
 
 1283         Matrix<T> u(A.
shape()), s(n, 1), v(n, n);
 
 1293                 t(k,l) = NumericTraits<T>::zero();
 
 1301         vigra_precondition(
false, 
"linearSolve(): Unknown solution method.");
 
 1306 template <
class T, 
class C1, 
int N>
 
 1307 bool linearSolve(MultiArrayView<2, T, C1> 
const & A, 
 
 1308                  TinyVector<T, N> 
const & b,
 
 1309                  TinyVector<T, N> & res, 
 
 1310                  std::string method = 
"QR")
 
 1313     return linearSolve(A, MultiArrayView<2, T>(shape, b.data()), MultiArrayView<2, T>(shape, res.data()), method);
 
 1333 #endif // VIGRA_LINEAR_SOLVE_HXX 
bool qrDecomposition(MultiArrayView< 2, T, C1 > const &a, MultiArrayView< 2, T, C2 > &q, MultiArrayView< 2, T, C3 > &r, double epsilon=0.0)
Definition: linear_solve.hxx:1016
FixedPoint16< 2, OverflowHandling > atan2(FixedPoint16< IntBits, OverflowHandling > y, FixedPoint16< IntBits, OverflowHandling > x)
Arctangent. Accuracy better than 1/3 degree (9 significant bits). 
Definition: fixedpoint.hxx:1654
MultiArrayView< 2, T, C > columnVector(MultiArrayView< 2, T, C > const &m, MultiArrayIndex d)
Definition: matrix.hxx:727
std::string tolower(std::string s)
Definition: utilities.hxx:96
MultiArrayIndex rowCount(const MultiArrayView< 2, T, C > &x)
Definition: matrix.hxx:671
void transpose(const MultiArrayView< 2, T, C1 > &v, MultiArrayView< 2, T, C2 > &r)
Definition: matrix.hxx:965
linalg::TemporaryMatrix< T > sin(MultiArrayView< 2, T, C > const &v)
bool linearSolveLowerTriangular(const MultiArrayView< 2, T, C1 > &l, const MultiArrayView< 2, T, C2 > &b, MultiArrayView< 2, T, C3 > x)
Definition: linear_solve.hxx:1118
const difference_type & shape() const 
Definition: multi_array.hxx:1648
linalg::TemporaryMatrix< T > sign(MultiArrayView< 2, T, C > const &v)
bool reverseElimination(const MultiArrayView< 2, T, C1 > &r, const MultiArrayView< 2, T, C2 > &b, MultiArrayView< 2, T, C3 > x)
Definition: linear_solve.hxx:1037
bool isSymmetric(const MultiArrayView< 2, T, C > &v)
Definition: matrix.hxx:779
void choleskySolve(MultiArrayView< 2, T, C1 > const &L, MultiArrayView< 2, T, C2 > const &b, MultiArrayView< 2, T, C3 > &x)
Definition: linear_solve.hxx:1168
std::ptrdiff_t MultiArrayIndex
Definition: multi_fwd.hxx:60
FixedPoint16< IntBits, OverflowHandling > hypot(FixedPoint16< IntBits, OverflowHandling > v1, FixedPoint16< IntBits, OverflowHandling > v2)
Length of hypotenuse. 
Definition: fixedpoint.hxx:1640
Definition: multi_fwd.hxx:63
NumericTraits< V >::Promote sum(TinyVectorBase< V, SIZE, D1, D2 > const &l)
sum of the vector's elements 
Definition: tinyvector.hxx:2073
linalg::TemporaryMatrix< T > sq(MultiArrayView< 2, T, C > const &v)
Matrix< T, ALLLOC >::NormType norm(const Matrix< T, ALLLOC > &a)
doxygen_overloaded_function(template<...> void separableConvolveBlockwise) template< unsigned int N
Separated convolution on ChunkedArrays. 
Matrix< T, ALLLOC >::SquaredNormType squaredNorm(const Matrix< T, ALLLOC > &a)
TinyVector< MultiArrayIndex, N > type
Definition: multi_shape.hxx:272
int argMax(MultiArrayView< 2, T, C > const &a)
Find the index of the maximum element in a matrix. 
Definition: matrix.hxx:1986
NormTraits< T >::SquaredNormType dot(const MultiArrayView< 2, T, C1 > &x, const MultiArrayView< 2, T, C2 > &y)
Definition: matrix.hxx:1342
double gamma(double x)
The gamma function. 
Definition: mathutil.hxx:1587
linalg::TemporaryMatrix< T > log(MultiArrayView< 2, T, C > const &v)
MultiArrayView< 2, T, C > rowVector(MultiArrayView< 2, T, C > const &m, MultiArrayIndex d)
Definition: matrix.hxx:697
MultiArrayIndex columnCount(const MultiArrayView< 2, T, C > &x)
Definition: matrix.hxx:684
bool choleskyDecomposition(MultiArrayView< 2, T, C1 > const &A, MultiArrayView< 2, T, C2 > &L)
Definition: linear_solve.hxx:959
T logDeterminant(MultiArrayView< 2, T, C1 > const &a)
Definition: linear_solve.hxx:911
bool linearSolveUpperTriangular(const MultiArrayView< 2, T, C1 > &r, const MultiArrayView< 2, T, C2 > &b, MultiArrayView< 2, T, C3 > x)
Definition: linear_solve.hxx:1068
unsigned int singularValueDecomposition(MultiArrayView< 2, T, C1 > const &A, MultiArrayView< 2, T, C2 > &U, MultiArrayView< 2, T, C3 > &S, MultiArrayView< 2, T, C4 > &V)
Definition: singular_value_decomposition.hxx:75
bool inverse(const MultiArrayView< 2, T, C1 > &v, MultiArrayView< 2, T, C2 > &res)
Definition: linear_solve.hxx:760
linalg::TemporaryMatrix< T > abs(MultiArrayView< 2, T, C > const &v)
linalg::TemporaryMatrix< T > cos(MultiArrayView< 2, T, C > const &v)
SquareRootTraits< FixedPoint< IntBits, FracBits > >::SquareRootResult sqrt(FixedPoint< IntBits, FracBits > v)
square root. 
Definition: fixedpoint.hxx:616
NumericTraits< T >::Promote determinant(MultiArrayView< 2, T, C1 > const &a, std::string method="default")
Definition: linear_solve.hxx:859