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,
94 const Vector<typename PositiveDefMatrix::ValueType>&
b,
95 Vector<typename PositiveDefMatrix::ValueType>&
x,
96 Preconditioner<typename PositiveDefMatrix::ValueType>& preconditioner,
97 const State& termination = terminationDefaults<typename PositiveDefMatrix::ValueType>());
122 template<
typename PositiveDefMatrix,
typename Interrupter>
125 const PositiveDefMatrix& A,
126 const Vector<typename PositiveDefMatrix::ValueType>& b,
127 Vector<typename PositiveDefMatrix::ValueType>& x,
128 Preconditioner<typename PositiveDefMatrix::ValueType>& preconditioner,
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); }
175 template<
typename Scalar>
void scale(
const Scalar&
s);
192 template<
typename OtherValueType>
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;
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>
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); }
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>
359 VecValueType
dot(
const Vector<VecValueType>& inVec)
const;
365 friend class ConstValueIter;
379 using ConstRowBase = RowBase<ConstRowData>;
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;
522 LinearOp(
const T& a_,
const T* x_,
const T* y_, T* out_):
a(a_), x(x_),
y(y_),
out(out_) {}
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_) {}
622 template<
typename Scalar>
635 a(a_), b(b_), binCount(binCount_), arraySize(arraySize_), reducetmp(reducetmp_) {}
639 const SizeType binSize = arraySize / binCount;
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();
675 if (arraySize < 1024) {
680 result += aData[
n] * bData[
n];
696 result += partialSums[
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) {}
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])
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();
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])
876 template<
typename ValueType, SizeType STENCIL_SIZE>
881 assert(row < mNumRows);
882 this->getRowEditor(row).setValue(col, val);
886 template<
typename ValueType, SizeType STENCIL_SIZE>
887 inline const ValueType&
890 assert(row < mNumRows);
891 return this->getConstRow(row).getValue(col);
895 template<
typename ValueType, SizeType STENCIL_SIZE>
896 inline const ValueType&
899 return this->getValue(row,col);
903 template<
typename ValueType, SizeType STENCIL_SIZE>
904 template<
typename Scalar>
912 RowEditor
row = mat->getRowEditor(
n);
922 template<
typename ValueType, SizeType STENCIL_SIZE>
923 template<
typename Scalar>
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) {
942 ConstRow
row = mat->getConstRow(
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) {
960 OPENVDB_THROW(ArithmeticError,
"matrix and input vector have incompatible sizes ("
961 << mNumRows <<
"x" << mNumRows <<
" vs. " << inVec.size() <<
")");
963 if (resultVec.size() != mNumRows) {
964 OPENVDB_THROW(ArithmeticError,
"matrix and result vector have incompatible sizes ("
965 << mNumRows <<
"x" << mNumRows <<
" vs. " << resultVec.size() <<
")");
968 vectorMultiply(inVec.data(), resultVec.data());
972 template<
typename ValueType, SizeType STENCIL_SIZE>
973 template<
typename VecValueType>
976 const VecValueType* inVec, VecValueType* resultVec)
const
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;
1003 const SparseStencilMatrix<OtherValueType, STENCIL_SIZE>*
b;
1008 template<
typename ValueType, SizeType STENCIL_SIZE>
1009 template<
typename OtherValueType>
1014 if (this->numRows() != other.
numRows())
return false;
1015 bool equal = tbb::parallel_reduce(
SizeRange(0, this->numRows()),
true,
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) {
1033 if (!std::isfinite(*it))
return false;
1044 template<
typename ValueType, SizeType STENCIL_SIZE>
1049 bool finite = tbb::parallel_reduce(
SizeRange(0, this->numRows()),
true,
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
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&
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>
1132 inline typename SparseStencilMatrix<ValueType, STENCIL_SIZE>::ConstValueIter
1133 SparseStencilMatrix<ValueType, STENCIL_SIZE>::RowBase<DataType>::cbegin()
const
1135 return ConstValueIter(mData);
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>
1174 const Vector<VecValueType>& inVec)
const
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>
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>
1228 assert(column < mNumColumns);
1230 RowData&
data = RowBase<>::mData;
1236 if (offset < data.mSize && data.mCols[offset] == column) {
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];
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>
1304 assert(r.size() == z.size());
1305 assert(r.size() ==
size);
1317 InitOp(
const MatrixType&
m,
ValueType*
v): mat(&m), vec(v) {}
1319 for (
SizeType n = range.begin(), N = range.end();
n <
N; ++
n) {
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];
1345 template<
typename MatrixType>
1346 class IncompleteCholeskyPreconditioner:
public Preconditioner<typename MatrixType::ValueType>
1349 struct CopyToLowerOp;
1363 , mLowerTriangular(matrix.
numRows())
1364 , mUpperTriangular(matrix.
numRows())
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;
1432 typename MatrixType::ConstRow
mask = matrix.getConstRow(j);
1433 typename MatrixType::ConstValueIter maskIter = mask.cbegin();
1434 for ( ; maskIter; ++maskIter) {
1436 if (i < j)
continue;
1442 a_ij -= a_ik * a_jk;
1444 row_i.setValue(j, a_ij);
1451 TransposeOp(matrix, mLowerTriangular, mUpperTriangular));
1456 bool isValid()
const override {
return mPassedCompatibilityCondition; }
1460 if (!mPassedCompatibilityCondition) {
1461 OPENVDB_THROW(ArithmeticError,
"invalid Cholesky decomposition");
1468 zVec.fill(zeroVal<ValueType>());
1471 if (size == 0)
return;
1473 assert(rVec.size() ==
size);
1474 assert(zVec.size() ==
size);
1477 mTempVec.fill(zeroVal<ValueType>());
1483 typename TriangularMatrix::ConstRow
row = mLowerTriangular.getConstRow(i);
1486 tmpData[i] = (rData[i] -
dot) / diagonal;
1487 if (!std::isfinite(tmpData[i])) {
1496 typename TriangularMatrix::ConstRow
row = mUpperTriangular.getConstRow(i);
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);
1525 const MatrixType* mat; TriangularMatrix*
lower;
1531 TransposeOp(
const MatrixType&
m,
const TriangularMatrix&
l, TriangularMatrix& 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) {
1541 if (column <
n)
continue;
1542 outRow.setValue(column,
lower->getValue(column,
n));
1546 const MatrixType* mat;
const TriangularMatrix*
lower; TriangularMatrix*
upper;
1549 TriangularMatrix mLowerTriangular;
1550 TriangularMatrix mUpperTriangular;
1551 Vector<ValueType> mTempVec;
1552 bool mPassedCompatibilityCondition;
1559 namespace internal {
1562 template<
typename T>
1570 template<
typename T>
1574 assert(xVec.size() == yVec.size());
1575 assert(xVec.size() == result.size());
1576 axpy(a, xVec.data(), yVec.data(), result.data(), xVec.size());
1581 template<
typename MatrixOperator,
typename VecValueType>
1584 const VecValueType* b, VecValueType*
r)
1587 A.vectorMultiply(x, r);
1593 template<
typename MatrixOperator,
typename T>
1597 assert(x.size() == b.size());
1598 assert(x.size() == r.size());
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;
1648 if (size != bVec.size()) {
1649 OPENVDB_THROW(ArithmeticError,
"A and b have incompatible sizes"
1650 << size <<
"x" << size <<
" vs. " << bVec.size() <<
")");
1652 if (size != xVec.size()) {
1653 OPENVDB_THROW(ArithmeticError,
"A and x have incompatible sizes"
1654 << size <<
"x" << size <<
" vs. " << xVec.size() <<
")");
1671 assert(rVec.isFinite());
1690 for ( ; iteration < termination.
iterations; ++iteration) {
1692 if (interrupter.wasInterrupted()) {
1693 OPENVDB_THROW(RuntimeError,
"conjugate gradient solver was interrupted");
1702 precond.apply(rVec, zVec);
1706 assert(std::isfinite(rDotZ));
1708 if (0 == iteration) {
1712 const ValueType beta = rDotZ / rDotZPrev;
1718 Amat.vectorMultiply(pVec, qVec);
1722 assert(std::isfinite(pAp));
1734 l2Error = rVec.l2Norm();
1735 minL2Error =
Min(l2Error, minL2Error);
1740 if (l2Error > 2 * minL2Error) {
1771 #endif // OPENVDB_MATH_CONJGRADIENT_HAS_BEEN_INCLUDED
void operator()(const SizeRange &range) const
SharedPtr< JacobiPreconditioner > Ptr
bool isFinite() const
Return true if every element of this matrix has a finite value.
void apply(const Vector< ValueType > &rVec, Vector< ValueType > &zVec) override
void parallel_for(int64_t start, int64_t end, std::function< void(int64_t index)> &&task, parallel_options opt=parallel_options(0, Split_Y, 1))
GLboolean GLboolean GLboolean b
std::string upper(string_view a)
Return an all-upper case version of a (locale-independent).
void scale(const Scalar &s)
Multiply each element of this vector by s.
RowEditor getRowEditor(SizeType row)
Return a read/write view onto the given row of this matrix.
void scale(const Scalar &)
Scale all of the entries in this row.
Lightweight, variable-length vector.
void resize(SizeType n)
Reset this vector to have n elements, with uninitialized values.
T * data()
Return a pointer to this vector's elements.
const ValueType & operator()(SizeType row, SizeType col) const
Return the value at the given coordinates.
#define OPENVDB_LOG_WARN(mesg)
void axpy(const T &a, const T *xVec, const T *yVec, T *resultVec, SizeType size)
Compute ax + y.
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...
const ValueType & getValue(SizeType row, SizeType col) const
Return the value at the given coordinates.
const TriangularMatrix & lowerMatrix() const
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...
bool isExactlyEqual(const T0 &a, const T1 &b)
Return true if a is exactly equal to b.
void apply(const Vector< ValueType > &r, Vector< ValueType > &z) override
void computeResidual(const MatrixOperator &A, const Vector< T > &x, const Vector< T > &b, Vector< T > &r)
Compute r = b − Ax.
void swap(UT::ArraySet< Key, MULTI, MAX_LOAD_FACTOR_256, Clearer, Hash, KeyEqual > &a, UT::ArraySet< Key, MULTI, MAX_LOAD_FACTOR_256, Clearer, Hash, KeyEqual > &b)
IncompleteCholeskyPreconditioner(const MatrixType &matrix)
void axpy(const T &a, const Vector< T > &xVec, const Vector< T > &yVec, Vector< T > &result)
Compute ax + y.
JacobiPreconditioner(const MatrixType &A)
ConstRow(const ValueType *valueHead, const SizeType *columnHead, const SizeType &rowSize)
Read-only accessor to a row of this matrix.
#define OPENVDB_USE_VERSION_NAMESPACE
**But if you need a or simply need to know when the task has note that the like this
Iterator over the stored values in a row of this matrix.
Information about the state of a conjugate gradient solution.
InfNormOp(const T *data_)
SizeType size() const
Return the number of rows in this matrix.
Dummy NOOP interrupter class defining interface.
std::string str() const
Return a string representation of this vector.
Tolerance for floating-point comparison.
SizeType numRows() const
Return the number of rows in this matrix.
const T * data() const
Return a pointer to this vector's elements.
const T & at(SizeType i) const
Return the value of this vector's ith element.
ValueType l2Norm() const
Return the L2 norm of this vector.
DeterministicDotProductOp(const T *a_, const T *b_, const SizeType binCount_, const SizeType arraySize_, T *reducetmp_)
Sparse, square matrix representing a 3D stencil operator of size STENCIL_SIZE.
tbb::blocked_range< SizeType > SizeRange
Read/write accessor to a row of this matrix.
std::shared_ptr< T > SharedPtr
typename TriangularMatrix::RowEditor TriangleRowEditor
T & operator[](SizeType i)
Return the value of this vector's ith element.
SYS_FORCE_INLINE const_iterator end() const
SharedPtr< Preconditioner > Ptr
bool operator()(const SizeRange &range, bool finite) const
T & at(SizeType i)
Return the value of this vector's ith element.
bool isApproxEqual(const Type &a, const Type &b, const Type &tolerance)
Return true if a is equal to b to within the given tolerance.
void operator()(const SizeRange &range) const
ImageBuf OIIO_API min(Image_or_Const A, Image_or_Const B, ROI roi={}, int nthreads=0)
std::ostream & operator<<(std::ostream &os, const State &state)
Preconditioner using incomplete Cholesky factorization.
Coord Abs(const Coord &xyz)
ValueType dot(const Vector &) const
Return the dot product of this vector with the given vector, which must be the same size...
SparseStencilMatrix(SizeType n)
Construct an n x n matrix with at most STENCIL_SIZE nonzero elements per row.
bool empty() const
Return true if this vector has no elements.
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.
FMT_CONSTEXPR bool find(Ptr first, Ptr last, T value, Ptr &out)
MatrixCopyOp(const SparseStencilMatrix &from_, SparseStencilMatrix &to_)
Vector(SizeType n)
Construct a vector of n elements, with uninitialized values.
General-purpose arithmetic and comparison routines, most of which accept arbitrary value types (or at...
const SparseStencilMatrix * mat
void clear()
Set the number of entries in this row to zero.
fpreal64 dot(const CE_VectorT< T > &a, const CE_VectorT< T > &b)
Base class for conjugate gradient preconditioners.
GLboolean GLboolean GLboolean GLboolean a
GLsizei const GLchar *const * string
const T * constData() const
Return a pointer to this vector's elements.
void computeResidual(const MatrixOperator &A, const VecValueType *x, const VecValueType *b, VecValueType *r)
Compute r = b − Ax.
ImageBuf OIIO_API max(Image_or_Const A, Image_or_Const B, ROI roi={}, int nthreads=0)
float Sqrt(float x)
Return the square root of a floating-point value.
static const ValueType sZeroValue
const TriangularMatrix & upperMatrix() const
GLdouble GLdouble GLdouble z
Vector & operator*=(const Scalar &s)
Multiply each element of this vector by s.
CopyOp(const T *from_, T *to_)
void swap(Vector &other)
Swap internal storage with another vector, which need not be the same size.
void operator()(const SizeRange &range) const
bool operator()(const SizeRange &range, bool finite) const
ValueType infNorm() const
Return the infinity norm of this vector.
Vector< ValueType > VectorType
const T & operator[](SizeType i) const
Return the value of this vector's ith element.
RowEditor(ValueType *valueHead, SizeType *columnHead, SizeType &rowSize, SizeType colSize)
Vector()
Construct an empty vector.
bool isApproxZero(const Type &x)
Return true if x is equal to zero to within the default floating-point comparison tolerance...
void setValue(SizeType row, SizeType col, const ValueType &)
Set the value at the given coordinates.
GLenum GLenum void void * column
Vector & operator=(const Vector &)
Deep copy the given vector.
FillOp(T *data_, const T &val_)
void fill(const ValueType &value)
Set all elements of this vector to value.
virtual ~Preconditioner()=default
SizeType size() const
Return the number of elements in this vector.
void scale(const Scalar &s)
Multiply all elements in the matrix by s;.
typename MatrixType::ValueType ValueType
SharedPtr< SparseStencilMatrix > Ptr
SharedPtr< IncompleteCholeskyPreconditioner > Ptr
#define OPENVDB_LOG_DEBUG_RUNTIME(mesg)
void vectorMultiply(const Vector< VecValueType > &inVec, Vector< VecValueType > &resultVec) const
Multiply this matrix by inVec and return the result in resultVec.
SparseStencilMatrix & operator*=(const Scalar &s)
Multiply all elements in the matrix by s;.
ConstRow getConstRow(SizeType row) const
Return a read-only view onto the given row of this matrix.
SizeType setValue(SizeType column, const ValueType &value)
Set the value of the entry in the specified column.
std::string lower(string_view a)
Return an all-upper case version of a (locale-independent).
virtual bool isValid() const
Preconditioner(const SparseStencilMatrix< T, STENCIL_SIZE > &)
RowEditor & operator*=(const Scalar &s)
Scale all of the entries in this row.
const Type & Min(const Type &a, const Type &b)
Return the minimum of two values.
GA_API const UT_StringHolder N
GLsizei const GLfloat * value
bool isFinite() const
Return true if all values along the diagonal are finite.
bool equal(T1 a, T2 b, T3 t)
const ValueType & operator*() const
T operator()(const SizeRange &range, T maxValue) const
bool isValid() const override
bool isFinite() const
Return true if every element of this vector has a finite value.
IsFiniteOp(const T *data_)
std::string str() const
Return a string representation of this matrix.
GLfloat GLfloat GLfloat alpha
void operator()(const SizeRange &range) const
#define OPENVDB_VERSION_NAME
The version namespace name for this library version.
virtual void apply(const Vector< T > &r, Vector< T > &z)=0
Apply this preconditioner to a residue vector: z = M−1r
typename TriangularMatrix::ConstRow TriangleConstRow
State terminationDefaults()
Return default termination conditions for a conjugate gradient solver.
const Type & Max(const Type &a, const Type &b)
Return the maximum of two values.
LinearOp(const T &a_, const T *x_, const T *y_, T *out_)
typename MatrixType::ValueType ValueType
bool isZero(const Type &x)
Return true if x is exactly equal to zero.
Vector(SizeType n, const ValueType &val)
Construct a vector of n elements and initialize each element to the given value.
#define OPENVDB_THROW(exception, message)
ConstValueIter & operator++()
IsFiniteOp(const SparseStencilMatrix &m)
void operator()(const SizeRange &range) const