9 #ifndef OPENVDB_MATH_CONJGRADIENT_HAS_BEEN_INCLUDED 10 #define OPENVDB_MATH_CONJGRADIENT_HAS_BEEN_INCLUDED 17 #include <tbb/parallel_for.h> 18 #include <tbb/parallel_reduce.h> 37 template<
typename ValueType>
class Vector;
55 template<
typename ValueType>
63 s.
absoluteError = std::numeric_limits<ValueType>::epsilon() * 100.0;
90 template<
typename PositiveDefMatrix>
93 const PositiveDefMatrix& A,
97 const State& termination = terminationDefaults<typename PositiveDefMatrix::ValueType>());
122 template<
typename PositiveDefMatrix,
typename Interrupter>
125 const PositiveDefMatrix& A,
129 Interrupter& interrupter,
130 const State& termination = terminationDefaults<typename PositiveDefMatrix::ValueType>());
151 ~Vector() { mSize = 0;
delete[] mData; mData =
nullptr; }
161 bool empty()
const {
return (mSize == 0); }
168 void swap(
Vector& other) { std::swap(mData, other.mData); std::swap(mSize, other.mSize); }
175 template<
typename Scalar>
void scale(
const Scalar& s);
192 template<
typename OtherValueType>
197 std::string str()
const;
209 inline T*
data() {
return mData; }
210 inline const T*
data()
const {
return mData; }
216 template<
typename Scalar>
struct ScaleOp;
217 struct DeterministicDotProductOp;
219 template<
typename OtherValueType>
struct EqOp;
236 template<
typename ValueType_, SizeType STENCIL_SIZE>
244 class ConstValueIter;
276 ConstRow getConstRow(
SizeType row)
const;
279 RowEditor getRowEditor(
SizeType row);
283 template<
typename Scalar>
void scale(
const Scalar& s);
284 template<
typename Scalar>
291 template<
typename VecValueType>
298 template<
typename VecValueType>
299 void vectorMultiply(
const VecValueType* inVec, VecValueType* resultVec)
const;
303 template<
typename OtherValueType>
311 std::string str()
const;
319 struct ConstRowData {
321 mVals(v), mCols(c), mSize(s) {}
326 template<
typename DataType_ = RowData>
330 using DataType = DataType_;
332 static SizeType capacity() {
return STENCIL_SIZE; }
334 RowBase(
const DataType& data): mData(data) {}
336 bool empty()
const {
return (mData.mSize == 0); }
337 const SizeType& size()
const {
return mData.mSize; }
343 ConstValueIter cbegin()
const;
347 template<
typename OtherDataType>
348 bool eq(
const RowBase<OtherDataType>& other,
354 template<
typename VecValueType>
355 VecValueType dot(
const VecValueType* inVec,
SizeType vecSize)
const;
358 template<
typename VecValueType>
362 std::string str()
const;
365 friend class ConstValueIter;
379 using ConstRowBase = RowBase<ConstRowData>;
388 if (mData.mSize == 0)
return SparseStencilMatrix::sZeroValue;
389 return mData.mVals[mCursor];
396 operator bool()
const {
return (mCursor < mData.mSize); }
402 ConstValueIter(
const RowData& d): mData(d.mVals, d.mCols, d.mSize), mCursor(0) {}
405 const ConstRowData mData;
433 template<
typename Scalar>
void scale(
const Scalar&);
434 template<
typename Scalar>
445 template<
typename VecValueType>
struct VecMultOp;
446 template<
typename Scalar>
struct RowScaleOp;
450 template<
typename OtherValueType>
struct EqOp;
453 std::unique_ptr<ValueType[]> mValueArray;
454 std::unique_ptr<SizeType[]> mColumnIdxArray;
455 std::unique_ptr<SizeType[]> mRowSizeArray;
492 CopyOp(
const T* from_, T* to_): from(from_), to(to_) {}
495 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) to[n] = from[n];
507 FillOp(T* data_,
const T& val_): data(data_), val(val_) {}
510 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) data[n] = val;
522 LinearOp(
const T& a_,
const T* x_,
const T* y_, T* out_): a(a_), x(x_), y(y_), out(out_) {}
526 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) out[n] = x[n] + y[n];
528 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) out[n] = -x[n] + y[n];
530 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) out[n] = a * x[n] + y[n];
547 os << (state.
success ?
"succeeded with " :
"")
572 if (mSize != other.mSize) {
575 mData =
new T[mSize];
591 if (mData)
delete[] mData;
607 template<
typename Scalar>
610 ScaleOp(T* data_,
const Scalar& s_):
data(data_), s(s_) {}
612 void operator()(
const SizeRange& range)
const {
613 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n)
data[n] *= s;
622 template<
typename Scalar>
626 tbb::parallel_for(
SizeRange(0, mSize), ScaleOp<Scalar>(mData, s));
635 a(a_), b(b_), binCount(binCount_), arraySize(arraySize_), reducetmp(reducetmp_) {}
639 const SizeType binSize = arraySize / binCount;
642 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) {
644 const SizeType end = (n == binCount-1) ? arraySize : begin + binSize;
647 T sum = zeroVal<T>();
648 for (
SizeType i = begin; i < end; ++i) { sum += a[i] * b[i]; }
666 assert(this->
size() == other.
size());
668 const T* aData = this->
data();
669 const T* bData = other.
data();
673 T result = zeroVal<T>();
675 if (arraySize < 1024) {
679 for (
SizeType n = 0; n < arraySize; ++n) {
680 result += aData[n] * bData[n];
692 tbb::parallel_for(
SizeRange(0, binCount),
695 for (
SizeType n = 0; n < binCount; ++n) {
696 result += partialSums[n];
711 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) {
726 T result = tbb::parallel_reduce(
SizeRange(0, this->
size()), zeroVal<T>(),
740 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) {
741 if (!std::isfinite(
data[n]))
return false;
756 bool finite = tbb::parallel_reduce(
SizeRange(0, this->
size()),
true,
758 [](
bool finite1,
bool finite2) {
return (finite1 && finite2); });
764 template<
typename OtherValueType>
767 EqOp(
const T* a_,
const OtherValueType* b_, T e): a(a_), b(b_), eps(e) {}
769 bool operator()(
const SizeRange& range,
bool equal)
const 772 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) {
780 const OtherValueType* b;
786 template<
typename OtherValueType>
790 if (this->
size() != other.
size())
return false;
791 bool equal = tbb::parallel_reduce(
SizeRange(0, this->
size()),
true,
792 EqOp<OtherValueType>(this->
data(), other.
data(), eps),
793 [](
bool eq1,
bool eq2) {
return (eq1 && eq2); });
802 std::ostringstream ostr;
806 ostr << sep << (*this)[n];
817 template<
typename ValueType, SizeType STENCIL_SIZE>
821 template<
typename ValueType, SizeType STENCIL_SIZE>
825 , mValueArray(new
ValueType[mNumRows * STENCIL_SIZE])
826 , mColumnIdxArray(new
SizeType[mNumRows * STENCIL_SIZE])
827 , mRowSizeArray(new
SizeType[mNumRows])
830 tbb::parallel_for(
SizeRange(0, mNumRows),
835 template<
typename ValueType, SizeType STENCIL_SIZE>
839 from(&from_), to(&to_) {}
843 const ValueType* fromVal = from->mValueArray.get();
844 const SizeType* fromCol = from->mColumnIdxArray.get();
845 ValueType* toVal = to->mValueArray.get();
846 SizeType* toCol = to->mColumnIdxArray.get();
847 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) {
848 toVal[n] = fromVal[n];
849 toCol[n] = fromCol[n];
857 template<
typename ValueType, SizeType STENCIL_SIZE>
860 : mNumRows(other.mNumRows)
861 , mValueArray(new
ValueType[mNumRows * STENCIL_SIZE])
862 , mColumnIdxArray(new
SizeType[mNumRows * STENCIL_SIZE])
863 , mRowSizeArray(new
SizeType[mNumRows])
871 tbb::parallel_for(
SizeRange(0, mNumRows),
876 template<
typename ValueType, SizeType STENCIL_SIZE>
881 assert(row < mNumRows);
886 template<
typename ValueType, SizeType STENCIL_SIZE>
890 assert(row < mNumRows);
895 template<
typename ValueType, SizeType STENCIL_SIZE>
903 template<
typename ValueType, SizeType STENCIL_SIZE>
904 template<
typename Scalar>
911 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) {
922 template<
typename ValueType, SizeType STENCIL_SIZE>
923 template<
typename Scalar>
928 tbb::parallel_for(
SizeRange(0, mNumRows), RowScaleOp<Scalar>(*
this, s));
932 template<
typename ValueType, SizeType STENCIL_SIZE>
933 template<
typename VecValueType>
937 mat(&m), in(i), out(o) {}
941 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) {
943 out[n] = row.dot(in, mat->numRows());
948 const VecValueType* in;
953 template<
typename ValueType, SizeType STENCIL_SIZE>
954 template<
typename VecValueType>
959 if (inVec.
size() != mNumRows) {
961 << mNumRows <<
"x" << mNumRows <<
" vs. " << inVec.
size() <<
")");
963 if (resultVec.
size() != mNumRows) {
965 << mNumRows <<
"x" << mNumRows <<
" vs. " << resultVec.
size() <<
")");
972 template<
typename ValueType, SizeType STENCIL_SIZE>
973 template<
typename VecValueType>
976 const VecValueType* inVec, VecValueType* resultVec)
const 979 tbb::parallel_for(
SizeRange(0, mNumRows),
980 VecMultOp<VecValueType>(*
this, inVec, resultVec));
984 template<
typename ValueType, SizeType STENCIL_SIZE>
985 template<
typename OtherValueType>
990 a(&a_), b(&b_), eps(e) {}
995 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) {
996 if (!a->getConstRow(n).eq(b->getConstRow(n), eps))
return false;
1004 const ValueType eps;
1008 template<
typename ValueType, SizeType STENCIL_SIZE>
1009 template<
typename OtherValueType>
1016 EqOp<OtherValueType>(*
this, other, eps),
1017 [](
bool eq1,
bool eq2) {
return (eq1 && eq2); });
1022 template<
typename ValueType, SizeType STENCIL_SIZE>
1030 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) {
1031 const ConstRow row = mat->getConstRow(n);
1033 if (!std::isfinite(*it))
return false;
1044 template<
typename ValueType, SizeType STENCIL_SIZE>
1050 IsFiniteOp(*
this), [](
bool finite1,
bool finite2) {
return (finite1&&finite2); });
1055 template<
typename ValueType, SizeType STENCIL_SIZE>
1059 std::ostringstream ostr;
1061 ostr << n <<
": " << this->
getConstRow(n).str() <<
"\n";
1067 template<
typename ValueType, SizeType STENCIL_SIZE>
1071 assert(i < mNumRows);
1072 const SizeType head = i * STENCIL_SIZE;
1073 return RowEditor(&mValueArray[head], &mColumnIdxArray[head], mRowSizeArray[i], mNumRows);
1077 template<
typename ValueType, SizeType STENCIL_SIZE>
1081 assert(i < mNumRows);
1082 const SizeType head = i * STENCIL_SIZE;
1083 return ConstRow(&mValueArray[head], &mColumnIdxArray[head], mRowSizeArray[i]);
1087 template<
typename ValueType, SizeType STENCIL_SIZE>
1088 template<
typename DataType>
1092 if (this->empty())
return mData.mSize;
1096 const SizeType* colPtr = std::lower_bound(mData.mCols, mData.mCols + mData.mSize, columnIdx);
1098 return static_cast<SizeType>(colPtr - mData.mCols);
1102 template<
typename ValueType, SizeType STENCIL_SIZE>
1103 template<
typename DataType>
1104 inline const ValueType&
1106 SizeType columnIdx,
bool& active)
const 1109 SizeType idx = this->find(columnIdx);
1110 if (idx < this->
size() && this->column(idx) == columnIdx) {
1112 return this->
value(idx);
1117 template<
typename ValueType, SizeType STENCIL_SIZE>
1118 template<
typename DataType>
1119 inline const ValueType&
1122 SizeType idx = this->find(columnIdx);
1123 if (idx < this->
size() && this->column(idx) == columnIdx) {
1124 return this->
value(idx);
1130 template<
typename ValueType, SizeType STENCIL_SIZE>
1131 template<
typename DataType>
1139 template<
typename ValueType, SizeType STENCIL_SIZE>
1140 template<
typename DataType>
1141 template<
typename OtherDataType>
1144 const RowBase<OtherDataType>& other, ValueType eps)
const 1146 if (this->
size() != other.size())
return false;
1147 for (
ConstValueIter it = cbegin(), oit = other.cbegin(); it || oit; ++it, ++oit) {
1148 if (it.column() != oit.column())
return false;
1155 template<
typename ValueType, SizeType STENCIL_SIZE>
1156 template<
typename DataType>
1157 template<
typename VecValueType>
1160 const VecValueType* inVec,
SizeType vecSize)
const 1162 VecValueType result = zeroVal<VecValueType>();
1164 result +=
static_cast<VecValueType
>(this->
value(idx) * inVec[this->column(idx)]);
1169 template<
typename ValueType, SizeType STENCIL_SIZE>
1170 template<
typename DataType>
1171 template<
typename VecValueType>
1176 return dot(inVec.
data(), inVec.
size());
1180 template<
typename ValueType, SizeType STENCIL_SIZE>
1181 template<
typename DataType>
1185 std::ostringstream ostr;
1188 ostr << sep <<
"(" << this->column(n) <<
", " << this->
value(n) <<
")";
1195 template<
typename ValueType, SizeType STENCIL_SIZE>
1198 const ValueType* valueHead,
const SizeType* columnHead,
const SizeType& rowSize):
1199 ConstRowBase(ConstRowData(const_cast<ValueType*>(valueHead),
1205 template<
typename ValueType, SizeType STENCIL_SIZE>
1209 RowBase<>(RowData(valueHead, columnHead, rowSize)), mNumColumns(colSize)
1214 template<
typename ValueType, SizeType STENCIL_SIZE>
1219 RowBase<>::mData.mSize = 0;
1223 template<
typename ValueType, SizeType STENCIL_SIZE>
1226 SizeType column,
const ValueType& value)
1228 assert(column < mNumColumns);
1230 RowData& data = RowBase<>::mData;
1234 SizeType offset = this->find(column);
1236 if (offset < data.mSize && data.mCols[offset] == column) {
1238 data.mVals[offset] = value;
1243 assert(data.mSize < this->capacity());
1245 if (offset >= data.mSize) {
1247 data.mVals[data.mSize] = value;
1248 data.mCols[data.mSize] = column;
1251 for (
SizeType i = data.mSize; i > offset; --i) {
1252 data.mVals[i] = data.mVals[i - 1];
1253 data.mCols[i] = data.mCols[i - 1];
1255 data.mVals[offset] = value;
1256 data.mCols[offset] = column;
1264 template<
typename ValueType, SizeType STENCIL_SIZE>
1265 template<
typename Scalar>
1269 for (
int idx = 0, N = this->size(); idx < N; ++idx) {
1270 RowBase<>::mData.mVals[idx] *= s;
1279 template<
typename MatrixType>
1287 using ValueType =
typename MatrixType::ValueType;
1295 tbb::parallel_for(
SizeRange(0, A.numRows()), InitOp(A, mDiag.data()));
1302 const SizeType size = mDiag.size();
1305 assert(r.
size() == size);
1307 tbb::parallel_for(
SizeRange(0, size), ApplyOp(mDiag.data(), r.
data(), z.
data()));
1317 InitOp(
const MatrixType& m, ValueType* v): mat(&m), vec(v) {}
1319 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) {
1320 const ValueType val = mat->
getValue(n, n);
1322 vec[n] =
static_cast<ValueType
>(1.0 / val);
1325 const MatrixType* mat; ValueType* vec;
1331 ApplyOp(
const ValueType* x_,
const ValueType* y_, ValueType* out_):
1332 x(x_), y(y_), out(out_) {}
1334 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) out[n] = x[n] * y[n];
1336 const ValueType *x, *y; ValueType* out;
1345 template<
typename MatrixType>
1349 struct CopyToLowerOp;
1353 using ValueType =
typename MatrixType::ValueType;
1363 , mLowerTriangular(matrix.
numRows())
1364 , mUpperTriangular(matrix.
numRows())
1371 tbb::parallel_for(
SizeRange(0, numRows), CopyToLowerOp(matrix, mLowerTriangular));
1391 mPassedCompatibilityCondition =
true;
1396 ValueType diagonalValue = crow_k.
getValue(k);
1399 if (diagonalValue < 1.e-5) {
1400 mPassedCompatibilityCondition =
false;
1404 diagonalValue =
Sqrt(diagonalValue);
1407 row_k.setValue(k, diagonalValue);
1410 typename MatrixType::ConstRow srcRow = matrix.getConstRow(k);
1411 typename MatrixType::ConstValueIter citer = srcRow.cbegin();
1412 for ( ; citer; ++citer) {
1414 if (ii < k+1)
continue;
1418 row_ii.setValue(k, *citer / diagonalValue);
1423 for ( ; citer; ++citer) {
1425 if (j < k+1)
continue;
1428 ValueType a_jk = row_j.
getValue(k);
1432 typename MatrixType::ConstRow mask = matrix.getConstRow(j);
1433 typename MatrixType::ConstValueIter maskIter = mask.cbegin();
1434 for ( ; maskIter; ++maskIter) {
1436 if (i < j)
continue;
1439 ValueType a_ij = crow_i.
getValue(j);
1440 ValueType a_ik = crow_i.
getValue(k);
1442 a_ij -= a_ik * a_jk;
1444 row_i.setValue(j, a_ij);
1450 tbb::parallel_for(
SizeRange(0, numRows),
1451 TransposeOp(matrix, mLowerTriangular, mUpperTriangular));
1456 bool isValid()
const override {
return mPassedCompatibilityCondition; }
1460 if (!mPassedCompatibilityCondition) {
1466 SizeType size = mLowerTriangular.numRows();
1468 zVec.
fill(zeroVal<ValueType>());
1469 ValueType* zData = zVec.
data();
1471 if (size == 0)
return;
1473 assert(rVec.
size() == size);
1474 assert(zVec.
size() == size);
1477 mTempVec.fill(zeroVal<ValueType>());
1478 ValueType* tmpData = mTempVec.data();
1479 const ValueType* rData = rVec.
data();
1482 for (
SizeType i = 0; i < size; ++i) {
1483 typename TriangularMatrix::ConstRow row = mLowerTriangular.getConstRow(i);
1484 ValueType diagonal = row.
getValue(i);
1485 ValueType dot = row.dot(mTempVec);
1486 tmpData[i] = (rData[i] - dot) / diagonal;
1487 if (!std::isfinite(tmpData[i])) {
1494 for (
SizeType ii = 0; ii < size; ++ii) {
1496 typename TriangularMatrix::ConstRow row = mUpperTriangular.getConstRow(i);
1497 ValueType diagonal = row.
getValue(i);
1498 ValueType dot = row.dot(zVec);
1499 zData[i] = (tmpData[i] - dot) / diagonal;
1500 if (!std::isfinite(zData[i])) {
1511 struct CopyToLowerOp
1513 CopyToLowerOp(
const MatrixType& m,
TriangularMatrix& l): mat(&m), lower(&l) {}
1515 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) {
1516 typename TriangularMatrix::RowEditor outRow = lower->getRowEditor(n);
1518 typename MatrixType::ConstRow inRow = mat->getConstRow(n);
1519 for (
typename MatrixType::ConstValueIter it = inRow.cbegin(); it; ++it) {
1520 if (it.column() > n)
continue;
1521 outRow.setValue(it.column(), *it);
1532 mat(&m), lower(&l), upper(&u) {}
1534 for (
SizeType n = range.begin(), N = range.end(); n < N; ++n) {
1535 typename TriangularMatrix::RowEditor outRow = upper->getRowEditor(n);
1538 typename MatrixType::ConstRow inRow = mat->getConstRow(n);
1539 for (
typename MatrixType::ConstValueIter it = inRow.cbegin(); it; ++it) {
1540 const SizeType column = it.column();
1541 if (column < n)
continue;
1542 outRow.setValue(column, lower->getValue(column, n));
1552 bool mPassedCompatibilityCondition;
1559 namespace internal {
1562 template<
typename T>
1564 axpy(
const T& a,
const T* xVec,
const T* yVec, T* resultVec,
SizeType size)
1570 template<
typename T>
1574 assert(xVec.
size() == yVec.
size());
1575 assert(xVec.
size() == result.
size());
1581 template<
typename MatrixOperator,
typename VecValueType>
1584 const VecValueType* b, VecValueType* r)
1587 A.vectorMultiply(x, r);
1593 template<
typename MatrixOperator,
typename T>
1599 assert(x.
size() == A.numRows());
1610 template<
typename PositiveDefMatrix>
1613 const PositiveDefMatrix& Amat,
1617 const State& termination)
1620 return solve(Amat, bVec, xVec, precond, interrupter, termination);
1624 template<
typename PositiveDefMatrix,
typename Interrupter>
1627 const PositiveDefMatrix& Amat,
1631 Interrupter& interrupter,
1632 const State& termination)
1634 using ValueType =
typename PositiveDefMatrix::ValueType;
1643 const SizeType size = Amat.numRows();
1648 if (size != bVec.
size()) {
1650 << size <<
"x" << size <<
" vs. " << bVec.
size() <<
")");
1652 if (size != xVec.
size()) {
1654 << size <<
"x" << size <<
" vs. " << xVec.
size() <<
")");
1663 const ValueType tmp = bVec.
infNorm();
1671 assert(rVec.isFinite());
1683 ValueType rDotZPrev(1);
1690 for ( ; iteration < termination.
iterations; ++iteration) {
1692 if (interrupter.wasInterrupted()) {
1702 precond.
apply(rVec, zVec);
1705 const ValueType rDotZ = rVec.dot(zVec);
1706 assert(std::isfinite(rDotZ));
1708 if (0 == iteration) {
1712 const ValueType beta = rDotZ / rDotZPrev;
1718 Amat.vectorMultiply(pVec, qVec);
1721 const ValueType pAp = pVec.dot(qVec);
1722 assert(std::isfinite(pAp));
1724 const ValueType alpha = rDotZ / pAp;
1734 l2Error = rVec.l2Norm();
1735 minL2Error =
Min(l2Error, minL2Error);
1740 if (l2Error > 2 * minL2Error) {
1771 #endif // OPENVDB_MATH_CONJGRADIENT_HAS_BEEN_INCLUDED JacobiPreconditioner(const MatrixType &A)
Definition: ConjGradient.h:1292
void operator()(const SizeRange &range) const
Definition: ConjGradient.h:509
void scale(const Scalar &)
Scale all of the entries in this row.
Definition: ConjGradient.h:1267
const T * data
Definition: ConjGradient.h:747
const T * data
Definition: ConjGradient.h:717
Definition: Exceptions.h:56
ValueType infNorm() const
Return the infinity norm of this vector.
Definition: ConjGradient.h:723
bool empty() const
Return true if this vector has no elements.
Definition: ConjGradient.h:161
ConstValueIter & operator++()
Definition: ConjGradient.h:395
Index32 SizeType
Definition: ConjGradient.h:33
void computeResidual(const MatrixOperator &A, const Vector< T > &x, const Vector< T > &b, Vector< T > &r)
Compute r = b − Ax.
Definition: ConjGradient.h:1595
SharedPtr< Preconditioner > Ptr
Definition: ConjGradient.h:468
bool isExactlyEqual(const T0 &a, const T1 &b)
Return true if a is exactly equal to b.
Definition: Math.h:444
Definition: ConjGradient.h:505
IsFiniteOp(const T *data_)
Definition: ConjGradient.h:735
bool isFinite() const
Return true if all values along the diagonal are finite.
Definition: ConjGradient.h:1311
ConstRow(const ValueType *valueHead, const SizeType *columnHead, const SizeType &rowSize)
Definition: ConjGradient.h:1197
#define OPENVDB_THROW(exception, message)
Definition: Exceptions.h:74
bool isFinite(const float x)
Return true if x is finite.
Definition: Math.h:376
typename TriangularMatrix::ConstRow TriangleConstRow
Definition: ConjGradient.h:1358
RowEditor(ValueType *valueHead, SizeType *columnHead, SizeType &rowSize, SizeType colSize)
Definition: ConjGradient.h:1207
General-purpose arithmetic and comparison routines, most of which accept arbitrary value types (or at...
Iterator over the stored values in a row of this matrix.
Definition: ConjGradient.h:383
#define OPENVDB_LOG_DEBUG_RUNTIME(message)
Log a debugging message in both debug and optimized builds.
Definition: logging.h:270
~Vector()
Definition: ConjGradient.h:151
Base class for interrupters.
Definition: NullInterrupter.h:25
LinearOp(const T &a_, const T *x_, const T *y_, T *out_)
Definition: ConjGradient.h:522
ValueType ValueType
Definition: ConjGradient.h:141
Definition: ConjGradient.h:490
Vector(SizeType n, const ValueType &val)
Construct a vector of n elements and initialize each element to the given value.
Definition: ConjGradient.h:149
double absoluteError
Definition: ConjGradient.h:50
IsFiniteOp(const SparseStencilMatrix &m)
Definition: ConjGradient.h:1025
const SparseStencilMatrix * mat
Definition: ConjGradient.h:1040
const TriangularMatrix & upperMatrix() const
Definition: ConjGradient.h:1507
Vector(SizeType n)
Construct a vector of n elements, with uninitialized values.
Definition: ConjGradient.h:147
const T val
Definition: ConjGradient.h:514
Definition: ConjGradient.h:631
void operator()(const SizeRange &range) const
Definition: ConjGradient.h:524
void vectorMultiply(const Vector< VecValueType > &inVec, Vector< VecValueType > &resultVec) const
Multiply this matrix by inVec and return the result in resultVec.
Definition: ConjGradient.h:956
void resize(SizeType n)
Reset this vector to have n elements, with uninitialized values.
Definition: ConjGradient.h:588
const SizeType binCount
Definition: ConjGradient.h:657
bool isApproxZero(const Type &x)
Return true if x is equal to zero to within the default floating-point comparison tolerance...
Definition: Math.h:350
Preconditioner(const SparseStencilMatrix< T, STENCIL_SIZE > &)
Definition: ConjGradient.h:470
void computeResidual(const MatrixOperator &A, const VecValueType *x, const VecValueType *b, VecValueType *r)
Compute r = b − Ax.
Definition: ConjGradient.h:1583
State terminationDefaults()
Return default termination conditions for a conjugate gradient solver.
Definition: ConjGradient.h:57
Vector & operator=(const Vector &)
Deep copy the given vector.
Definition: ConjGradient.h:568
Tolerance for floating-point comparison.
Definition: Math.h:147
void operator()(const SizeRange &range) const
Definition: ConjGradient.h:637
std::shared_ptr< T > SharedPtr
Definition: Types.h:114
Definition: ConjGradient.h:836
const T * data() const
Return a pointer to this vector's elements.
Definition: ConjGradient.h:210
SparseStencilMatrix * to
Definition: ConjGradient.h:853
Definition: ConjGradient.h:520
void fill(const ValueType &value)
Set all elements of this vector to value.
Definition: ConjGradient.h:600
Definition: ConjGradient.h:733
void operator()(const SizeRange &range) const
Definition: ConjGradient.h:841
InfNormOp(const T *data_)
Definition: ConjGradient.h:707
void axpy(const T &a, const T *xVec, const T *yVec, T *resultVec, SizeType size)
Compute ax + y.
Definition: ConjGradient.h:1564
Lightweight, variable-length vector.
Definition: ConjGradient.h:37
int iterations
Definition: ConjGradient.h:48
Vector()
Construct an empty vector.
Definition: ConjGradient.h:145
SharedPtr< JacobiPreconditioner > Ptr
Definition: ConjGradient.h:1290
double relativeError
Definition: ConjGradient.h:49
Information about the state of a conjugate gradient solution.
Definition: ConjGradient.h:46
T operator()(const SizeRange &range, T maxValue) const
Definition: ConjGradient.h:709
SizeType size() const
Return the number of rows in this matrix.
Definition: ConjGradient.h:259
SizeType numRows() const
Return the number of rows in this matrix.
Definition: ConjGradient.h:258
Definition: ConjGradient.h:705
const T * y
Definition: ConjGradient.h:534
std::string str() const
Return a string representation of this matrix.
Definition: ConjGradient.h:1057
const ValueType & operator*() const
Definition: ConjGradient.h:386
const T * a
Definition: ConjGradient.h:655
CopyOp(const T *from_, T *to_)
Definition: ConjGradient.h:492
float Sqrt(float x)
Return the square root of a floating-point value.
Definition: Math.h:764
RowEditor & operator*=(const Scalar &s)
Scale all of the entries in this row.
Definition: ConjGradient.h:435
void apply(const Vector< ValueType > &rVec, Vector< ValueType > &zVec) override
Apply this preconditioner to a residue vector: z = M−1r
Definition: ConjGradient.h:1458
void operator()(const SizeRange &range) const
Definition: ConjGradient.h:494
DeterministicDotProductOp(const T *a_, const T *b_, const SizeType binCount_, const SizeType arraySize_, T *reducetmp_)
Definition: ConjGradient.h:633
Sparse, square matrix representing a 3D stencil operator of size STENCIL_SIZE.
Definition: ConjGradient.h:39
T * data()
Return a pointer to this vector's elements.
Definition: ConjGradient.h:209
MatrixType::ValueType ValueType
Definition: ConjGradient.h:467
T * out
Definition: ConjGradient.h:535
State solve(const PositiveDefMatrix &A, const Vector< typename PositiveDefMatrix::ValueType > &b, Vector< typename PositiveDefMatrix::ValueType > &x, Preconditioner< typename PositiveDefMatrix::ValueType > &preconditioner, const State &termination=terminationDefaults< typename PositiveDefMatrix::ValueType >())
Solve Ax = b via the preconditioned conjugate gradient method.
Definition: ConjGradient.h:1612
ValueType dot(const Vector &) const
Return the dot product of this vector with the given vector, which must be the same size...
Definition: ConjGradient.h:664
Definition: Exceptions.h:13
const T * b
Definition: ConjGradient.h:656
const TriangularMatrix & lowerMatrix() const
Definition: ConjGradient.h:1506
Base class for conjugate gradient preconditioners.
Definition: ConjGradient.h:41
SizeType setValue(SizeType column, const ValueType &value)
Set the value of the entry in the specified column.
Definition: ConjGradient.h:1225
SharedPtr< IncompleteCholeskyPreconditioner > Ptr
Definition: ConjGradient.h:1356
ValueT value
Definition: GridBuilder.h:1287
bool isFinite() const
Return true if every element of this matrix has a finite value.
Definition: ConjGradient.h:1046
Definition: Exceptions.h:63
static const ValueType sZeroValue
Definition: ConjGradient.h:246
Vector & operator*=(const Scalar &s)
Multiply each element of this vector by s.
Definition: ConjGradient.h:176
OtherValueType ValueType
Definition: ConjGradient.h:240
void apply(const Vector< ValueType > &r, Vector< ValueType > &z) override
Apply this preconditioner to a residue vector: z = M−1r
Definition: ConjGradient.h:1300
const Type & Max(const Type &a, const Type &b)
Return the maximum of two values.
Definition: Math.h:598
T & at(SizeType i)
Return the value of this vector's ith element.
Definition: ConjGradient.h:201
void scale(const Scalar &s)
Multiply all elements in the matrix by s;.
Definition: ConjGradient.h:925
bool operator()(const SizeRange &range, bool finite) const
Definition: ConjGradient.h:737
const Type & Min(const Type &a, const Type &b)
Return the minimum of two values.
Definition: Math.h:659
MatType scale(const Vec3< typename MatType::value_type > &s)
Return a matrix that scales by s.
Definition: Mat.h:637
Definition: ConjGradient.h:1023
std::string str() const
Return a string representation of this vector.
Definition: ConjGradient.h:800
virtual void apply(const Vector< T > &r, Vector< T > &z)=0
Apply this preconditioner to a residue vector: z = M−1r
SparseStencilMatrix & operator*=(const Scalar &s)
Multiply all elements in the matrix by s;.
Definition: ConjGradient.h:285
Preconditioner using incomplete Cholesky factorization.
Definition: ConjGradient.h:43
const T & operator[](SizeType i) const
Return the value of this vector's ith element.
Definition: ConjGradient.h:204
void reset()
Definition: ConjGradient.h:398
const SizeType arraySize
Definition: ConjGradient.h:658
IncompleteCholeskyPreconditioner(const MatrixType &matrix)
Definition: ConjGradient.h:1361
Vector< ValueType > VectorType
Definition: ConjGradient.h:241
bool eq(const SparseStencilMatrix< OtherValueType, STENCIL_SIZE > &other, ValueType eps=Tolerance< ValueType >::value()) const
Return true if this matrix is equivalent to the given matrix to within the specified tolerance...
Definition: ConjGradient.h:1011
bool isApproxEqual(const Type &a, const Type &b, const Type &tolerance)
Return true if a is equal to b to within the given tolerance.
Definition: Math.h:407
Diagonal preconditioner.
Definition: ConjGradient.h:42
ConstRow getConstRow(SizeType row) const
Return a read-only view onto the given row of this matrix.
Definition: ConjGradient.h:1079
void scale(const Scalar &s)
Multiply each element of this vector by s.
Definition: ConjGradient.h:624
SharedPtr< SparseStencilMatrix > Ptr
Definition: ConjGradient.h:242
T * reducetmp
Definition: ConjGradient.h:659
State solve(const PositiveDefMatrix &A, const Vector< typename PositiveDefMatrix::ValueType > &b, Vector< typename PositiveDefMatrix::ValueType > &x, Preconditioner< typename PositiveDefMatrix::ValueType > &preconditioner, Interrupter &interrupter, const State &termination=terminationDefaults< typename PositiveDefMatrix::ValueType >())
Solve Ax = b via the preconditioned conjugate gradient method.
Definition: ConjGradient.h:1626
Read/write accessor to a row of this matrix.
Definition: ConjGradient.h:419
RowEditor getRowEditor(SizeType row)
Return a read/write view onto the given row of this matrix.
Definition: ConjGradient.h:1069
#define OPENVDB_LOG_WARN(message)
Log a warning message of the form 'someVar << "some text" << ...'.
Definition: logging.h:256
Read-only accessor to a row of this matrix.
Definition: ConjGradient.h:411
T * data
Definition: ConjGradient.h:513
SparseStencilMatrix(SizeType n)
Construct an n x n matrix with at most STENCIL_SIZE nonzero elements per row.
Definition: ConjGradient.h:823
typename TriangularMatrix::RowEditor TriangleRowEditor
Definition: ConjGradient.h:1359
MatrixCopyOp(const SparseStencilMatrix &from_, SparseStencilMatrix &to_)
Definition: ConjGradient.h:838
const T & at(SizeType i) const
Return the value of this vector's ith element.
Definition: ConjGradient.h:202
bool isFinite() const
Return true if every element of this vector has a finite value.
Definition: ConjGradient.h:753
SizeType column() const
Definition: ConjGradient.h:392
std::ostream & operator<<(std::ostream &os, const State &state)
Definition: ConjGradient.h:545
void swap(Vector &other)
Swap internal storage with another vector, which need not be the same size.
Definition: ConjGradient.h:168
SizeType size() const
Return the number of elements in this vector.
Definition: ConjGradient.h:159
uint32_t Index32
Definition: Types.h:52
Coord Abs(const Coord &xyz)
Definition: Coord.h:514
void axpy(const T &a, const Vector< T > &xVec, const Vector< T > &yVec, Vector< T > &result)
Compute ax + y.
Definition: ConjGradient.h:1572
bool isZero(const Type &x)
Return true if x is exactly equal to zero.
Definition: Math.h:338
#define OPENVDB_VERSION_NAME
The version namespace name for this library version.
Definition: version.h.in:116
bool operator()(const SizeRange &range, bool finite) const
Definition: ConjGradient.h:1027
bool isValid() const override
Definition: ConjGradient.h:1456
void increment()
Definition: ConjGradient.h:394
T * to
Definition: ConjGradient.h:499
const T * constData() const
Return a pointer to this vector's elements.
Definition: ConjGradient.h:211
T & operator[](SizeType i)
Return the value of this vector's ith element.
Definition: ConjGradient.h:203
const ValueType & getValue(SizeType row, SizeType col) const
Return the value at the given coordinates.
Definition: ConjGradient.h:888
const T * from
Definition: ConjGradient.h:498
const ValueType & operator()(SizeType row, SizeType col) const
Return the value at the given coordinates.
Definition: ConjGradient.h:897
SharedPtr< Vector > Ptr
Definition: ConjGradient.h:142
virtual bool isValid() const
Definition: ConjGradient.h:473
FillOp(T *data_, const T &val_)
Definition: ConjGradient.h:507
#define OPENVDB_USE_VERSION_NAMESPACE
Definition: version.h.in:202
void clear()
Set the number of entries in this row to zero.
Definition: ConjGradient.h:1216
bool success
Definition: ConjGradient.h:47
tbb::blocked_range< SizeType > SizeRange
Definition: ConjGradient.h:35
ValueType l2Norm() const
Return the L2 norm of this vector.
Definition: ConjGradient.h:185
void setValue(SizeType row, SizeType col, const ValueType &)
Set the value at the given coordinates.
Definition: ConjGradient.h:878
bool eq(const Vector< OtherValueType > &other, ValueType eps=Tolerance< ValueType >::value()) const
Return true if this vector is equivalent to the given vector to within the specified tolerance...
Definition: ConjGradient.h:788