OpenVDB  3.1.0
ConjGradient.h
Go to the documentation of this file.
1 //
3 // Copyright (c) 2012-2015 DreamWorks Animation LLC
4 //
5 // All rights reserved. This software is distributed under the
6 // Mozilla Public License 2.0 ( http://www.mozilla.org/MPL/2.0/ )
7 //
8 // Redistributions of source code must retain the above copyright
9 // and license notice and the following restrictions and disclaimer.
10 //
11 // * Neither the name of DreamWorks Animation nor the names of
12 // its contributors may be used to endorse or promote products derived
13 // from this software without specific prior written permission.
14 //
15 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
16 // "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
17 // LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
18 // A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
19 // OWNER OR CONTRIBUTORS BE LIABLE FOR ANY INDIRECT, INCIDENTAL,
20 // SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
21 // LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
22 // DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
23 // THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
24 // (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
25 // OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
26 // IN NO EVENT SHALL THE COPYRIGHT HOLDERS' AND CONTRIBUTORS' AGGREGATE
27 // LIABILITY FOR ALL CLAIMS REGARDLESS OF THEIR BASIS EXCEED US$250.00.
28 //
30 
35 
36 #ifndef OPENVDB_MATH_CONJGRADIENT_HAS_BEEN_INCLUDED
37 #define OPENVDB_MATH_CONJGRADIENT_HAS_BEEN_INCLUDED
38 
39 #include <openvdb/Exceptions.h>
40 #include <openvdb/Types.h>
41 #include <openvdb/util/logging.h>
43 #include "Math.h" // for Abs(), isZero(), Max(), Sqrt()
44 #include <boost/shared_ptr.hpp>
45 #include <boost/scoped_array.hpp>
46 #include <tbb/parallel_for.h>
47 #include <tbb/parallel_reduce.h>
48 #include <algorithm> // for std::lower_bound()
49 #include <cassert>
50 #include <cmath> // for std::isfinite()
51 #include <limits>
52 #include <sstream>
53 
54 
55 namespace openvdb {
57 namespace OPENVDB_VERSION_NAME {
58 namespace math {
59 namespace pcg {
60 
61 typedef Index32 SizeType;
62 
63 typedef tbb::blocked_range<SizeType> SizeRange;
64 
65 template<typename ValueType> class Vector;
66 
67 template<typename ValueType, SizeType STENCIL_SIZE> class SparseStencilMatrix;
68 
69 template<typename ValueType> class Preconditioner;
70 template<typename MatrixType> class JacobiPreconditioner;
71 template<typename MatrixType> class IncompleteCholeskyPreconditioner;
72 
74 struct State {
75  bool success;
77  double relativeError;
78  double absoluteError;
79 };
80 
81 
83 template<typename ValueType>
84 inline State
86 {
87  State s;
88  s.success = false;
89  s.iterations = 50;
90  s.relativeError = 1.0e-6;
91  s.absoluteError = std::numeric_limits<ValueType>::epsilon() * 100.0;
92  return s;
93 }
94 
95 
97 
98 
118 template<typename PositiveDefMatrix>
119 inline State
120 solve(
121  const PositiveDefMatrix& A,
122  const Vector<typename PositiveDefMatrix::ValueType>& b,
123  Vector<typename PositiveDefMatrix::ValueType>& x,
124  Preconditioner<typename PositiveDefMatrix::ValueType>& preconditioner,
125  const State& termination = terminationDefaults<typename PositiveDefMatrix::ValueType>());
126 
127 
150 template<typename PositiveDefMatrix, typename Interrupter>
151 inline State
152 solve(
153  const PositiveDefMatrix& A,
154  const Vector<typename PositiveDefMatrix::ValueType>& b,
155  Vector<typename PositiveDefMatrix::ValueType>& x,
156  Preconditioner<typename PositiveDefMatrix::ValueType>& preconditioner,
157  Interrupter& interrupter,
158  const State& termination = terminationDefaults<typename PositiveDefMatrix::ValueType>());
159 
160 
162 
163 
165 template<typename T>
166 class Vector
167 {
168 public:
169  typedef T ValueType;
170  typedef boost::shared_ptr<Vector> Ptr;
171 
173  Vector(): mData(NULL), mSize(0) {}
175  Vector(SizeType n): mData(new T[n]), mSize(n) {}
177  Vector(SizeType n, const ValueType& val): mData(new T[n]), mSize(n) { this->fill(val); }
178 
179  ~Vector() { mSize = 0; delete[] mData; mData = NULL; }
180 
182  Vector(const Vector&);
184  Vector& operator=(const Vector&);
185 
187  SizeType size() const { return mSize; }
189  bool empty() const { return (mSize == 0); }
190 
193  void resize(SizeType n);
194 
196  void swap(Vector& other) { std::swap(mData, other.mData); std::swap(mSize, other.mSize); }
197 
199  void fill(const ValueType& value);
200 
202  template<typename Scalar> void scale(const Scalar& s);
204  template<typename Scalar> Vector& operator*=(const Scalar& s) { this->scale(s); return *this; }
206 
208  ValueType dot(const Vector&) const;
209 
211  ValueType infNorm() const;
213  ValueType l2Norm() const { return Sqrt(this->dot(*this)); }
214 
216  bool isFinite() const;
217 
220  template<typename OtherValueType>
221  bool eq(const Vector<OtherValueType>& other,
222  ValueType eps = Tolerance<ValueType>::value()) const;
223 
225  std::string str() const;
226 
228  inline T& at(SizeType i) { return mData[i]; }
230  inline const T& at(SizeType i) const { return mData[i]; }
231  inline T& operator[](SizeType i) { return this->at(i); }
232  inline const T& operator[](SizeType i) const { return this->at(i); }
234 
236  inline T* data() { return mData; }
238  inline const T* data() const { return mData; }
239  inline const T* constData() const { return mData; }
241 
242 private:
243  // Functor for use with tbb::parallel_for()
244  template<typename Scalar> struct ScaleOp;
245  struct DeterministicDotProductOp;
246  // Functors for use with tbb::parallel_reduce()
247  template<typename OtherValueType> struct EqOp;
248  struct InfNormOp;
249  struct IsFiniteOp;
250 
251  T* mData;
252  SizeType mSize;
253 };
254 
257 //typedef Vector<double> LinearVector;
258 
259 
261 
262 
265 template<typename ValueType_, SizeType STENCIL_SIZE>
267 {
268 public:
269  typedef ValueType_ ValueType;
271  typedef boost::shared_ptr<SparseStencilMatrix> Ptr;
272 
273  class ConstValueIter;
274  class ConstRow;
275  class RowEditor;
276 
277  static const ValueType sZeroValue;
278 
280  SparseStencilMatrix(SizeType n);
281 
284 
286  SizeType numRows() const { return mNumRows; }
288  SizeType size() const { return mNumRows; }
290 
294  void setValue(SizeType row, SizeType col, const ValueType&);
295 
297  const ValueType& getValue(SizeType row, SizeType col) const;
301  const ValueType& operator()(SizeType row, SizeType col) const;
303 
305  ConstRow getConstRow(SizeType row) const;
306 
308  RowEditor getRowEditor(SizeType row);
309 
311  template<typename Scalar> void scale(const Scalar& s);
313  template<typename Scalar>
314  SparseStencilMatrix& operator*=(const Scalar& s) { this->scale(s); return *this; }
316 
320  template<typename VecValueType>
321  void vectorMultiply(const Vector<VecValueType>& inVec, Vector<VecValueType>& resultVec) const;
322 
327  template<typename VecValueType>
328  void vectorMultiply(const VecValueType* inVec, VecValueType* resultVec) const;
329 
332  template<typename OtherValueType>
334  ValueType eps = Tolerance<ValueType>::value()) const;
335 
337  bool isFinite() const;
338 
340  std::string str() const;
341 
342 private:
343  struct RowData {
344  RowData(ValueType* v, SizeType* c, SizeType& s): mVals(v), mCols(c), mSize(s) {}
345  ValueType* mVals; SizeType* mCols; SizeType& mSize;
346  };
347 
348  struct ConstRowData {
349  ConstRowData(const ValueType* v, const SizeType* c, const SizeType& s):
350  mVals(v), mCols(c), mSize(s) {}
351  const ValueType* mVals; const SizeType* mCols; const SizeType& mSize;
352  };
353 
355  template<typename DataType_ = RowData>
356  class RowBase
357  {
358  public:
359  typedef DataType_ DataType;
360 
361  static SizeType capacity() { return STENCIL_SIZE; }
362 
363  RowBase(const DataType& data): mData(data) {}
364 
365  bool empty() const { return (mData.mSize == 0); }
366  const SizeType& size() const { return mData.mSize; }
367 
368  const ValueType& getValue(SizeType columnIdx, bool& active) const;
369  const ValueType& getValue(SizeType columnIdx) const;
370 
372  ConstValueIter cbegin() const;
373 
376  template<typename OtherDataType>
377  bool eq(const RowBase<OtherDataType>& other,
378  ValueType eps = Tolerance<ValueType>::value()) const;
379 
383  template<typename VecValueType>
384  VecValueType dot(const VecValueType* inVec, SizeType vecSize) const;
385 
387  template<typename VecValueType>
388  VecValueType dot(const Vector<VecValueType>& inVec) const;
389 
391  std::string str() const;
392 
393  protected:
394  friend class ConstValueIter;
395 
396  const ValueType& value(SizeType i) const { return mData.mVals[i]; }
397  SizeType column(SizeType i) const { return mData.mCols[i]; }
398 
403  SizeType find(SizeType columnIdx) const;
404 
405  DataType mData;
406  };
407 
408  typedef RowBase<ConstRowData> ConstRowBase;
409 
410 public:
413  {
414  public:
415  const ValueType& operator*() const
416  {
417  if (mData.mSize == 0) return SparseStencilMatrix::sZeroValue;
418  return mData.mVals[mCursor];
419  }
420 
421  SizeType column() const { return mData.mCols[mCursor]; }
422 
423  void increment() { mCursor++; }
424  ConstValueIter& operator++() { increment(); return *this; }
425  operator bool() const { return (mCursor < mData.mSize); }
426 
427  void reset() { mCursor = 0; }
428 
429  private:
430  friend class SparseStencilMatrix;
431  ConstValueIter(const RowData& d): mData(d.mVals, d.mCols, d.mSize), mCursor(0) {}
432  ConstValueIter(const ConstRowData& d): mData(d), mCursor(0) {}
433 
434  const ConstRowData mData;
435  SizeType mCursor;
436  };
437 
438 
440  class ConstRow: public ConstRowBase
441  {
442  public:
443  ConstRow(const ValueType* valueHead, const SizeType* columnHead, const SizeType& rowSize);
444  }; // class ConstRow
445 
446 
448  class RowEditor: public RowBase<>
449  {
450  public:
451  RowEditor(ValueType* valueHead, SizeType* columnHead, SizeType& rowSize, SizeType colSize);
452 
454  void clear();
455 
458  SizeType setValue(SizeType column, const ValueType& value);
459 
461  template<typename Scalar> void scale(const Scalar&);
463  template<typename Scalar>
464  RowEditor& operator*=(const Scalar& s) { this->scale(s); return *this; }
466 
467  private:
468  const SizeType mNumColumns; // used only for bounds checking
469  }; // class RowEditor
470 
471 private:
472  // Functors for use with tbb::parallel_for()
473  struct MatrixCopyOp;
474  template<typename VecValueType> struct VecMultOp;
475  template<typename Scalar> struct RowScaleOp;
476 
477  // Functors for use with tbb::parallel_reduce()
478  struct IsFiniteOp;
479  template<typename OtherValueType> struct EqOp;
480 
481  const SizeType mNumRows;
482  boost::scoped_array<ValueType> mValueArray;
483  boost::scoped_array<SizeType> mColumnIdxArray;
484  boost::scoped_array<SizeType> mRowSizeArray;
485 }; // class SparseStencilMatrix
486 
487 
489 
490 
492 template<typename T>
493 class Preconditioner
494 {
495 public:
496  typedef T ValueType;
497  typedef boost::shared_ptr<Preconditioner> Ptr;
498 
499  template<SizeType STENCIL_SIZE> Preconditioner(const SparseStencilMatrix<T, STENCIL_SIZE>&) {}
500  virtual ~Preconditioner() {}
501 
502  virtual bool isValid() const { return true; }
503 
508  virtual void apply(const Vector<T>& r, Vector<T>& z) = 0;
509 };
510 
511 
513 
514 
515 namespace internal {
516 
517 // Functor for use with tbb::parallel_for() to copy data from one array to another
518 template<typename T>
519 struct CopyOp
520 {
521  CopyOp(const T* from_, T* to_): from(from_), to(to_) {}
522 
523  void operator()(const SizeRange& range) const {
524  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) to[n] = from[n];
525  }
526 
527  const T* from;
528  T* to;
529 };
530 
531 
532 // Functor for use with tbb::parallel_for() to fill an array with a constant value
533 template<typename T>
534 struct FillOp
535 {
536  FillOp(T* data_, const T& val_): data(data_), val(val_) {}
537 
538  void operator()(const SizeRange& range) const {
539  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) data[n] = val;
540  }
541 
542  T* data;
543  const T val;
544 };
545 
546 
547 // Functor for use with tbb::parallel_for() that computes a * x + y
548 template<typename T>
549 struct LinearOp
550 {
551  LinearOp(const T& a_, const T* x_, const T* y_, T* out_): a(a_), x(x_), y(y_), out(out_) {}
552 
553  void operator()(const SizeRange& range) const {
554  if (isExactlyEqual(a, T(1))) {
555  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) out[n] = x[n] + y[n];
556  } else if (isExactlyEqual(a, T(-1))) {
557  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) out[n] = -x[n] + y[n];
558  } else {
559  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) out[n] = a * x[n] + y[n];
560  }
561  }
562 
563  const T a, *x, *y;
564  T* out;
565 };
566 
567 } // namespace internal
568 
569 
571 
572 
573 inline std::ostream&
574 operator<<(std::ostream& os, const State& state)
575 {
576  os << (state.success ? "succeeded with " : "")
577  << "rel. err. " << state.relativeError << ", abs. err. " << state.absoluteError
578  << " after " << state.iterations << " iteration" << (state.iterations == 1 ? "" : "s");
579  return os;
580 }
581 
582 
584 
585 
586 template<typename T>
587 inline
588 Vector<T>::Vector(const Vector& other): mData(new T[other.mSize]), mSize(other.mSize)
589 {
590  tbb::parallel_for(SizeRange(0, mSize),
591  internal::CopyOp<T>(/*from=*/other.mData, /*to=*/mData));
592 }
593 
594 
595 template<typename T>
596 inline
598 {
599  // Update the internal storage to the correct size
600 
601  if (mSize != other.mSize) {
602  mSize = other.mSize;
603  delete[] mData;
604  mData = new T[mSize];
605  }
606 
607  // Deep copy the data
608  tbb::parallel_for(SizeRange(0, mSize),
609  internal::CopyOp<T>(/*from=*/other.mData, /*to=*/mData));
610 
611  return *this;
612 }
613 
614 
615 template<typename T>
616 inline void
617 Vector<T>::resize(SizeType n)
618 {
619  if (n != mSize) {
620  if (mData) delete[] mData;
621  mData = new T[n];
622  mSize = n;
623  }
624 }
625 
626 
627 template<typename T>
628 inline void
629 Vector<T>::fill(const ValueType& value)
630 {
631  tbb::parallel_for(SizeRange(0, mSize), internal::FillOp<T>(mData, value));
632 }
633 
634 
635 template<typename T>
636 template<typename Scalar>
637 struct Vector<T>::ScaleOp
638 {
639  ScaleOp(T* data_, const Scalar& s_): data(data_), s(s_) {}
640 
641  void operator()(const SizeRange& range) const {
642  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) data[n] *= s;
643  }
644 
645  T* data;
646  const Scalar s;
647 };
648 
649 
650 template<typename T>
651 template<typename Scalar>
652 inline void
653 Vector<T>::scale(const Scalar& s)
654 {
655  tbb::parallel_for(SizeRange(0, mSize), ScaleOp<Scalar>(mData, s));
656 }
657 
658 
659 template<typename T>
661 {
662  DeterministicDotProductOp(const T* a_, const T* b_,
663  const SizeType binCount_, const SizeType arraySize_, T* reducetmp_):
664  a(a_), b(b_), binCount(binCount_), arraySize(arraySize_), reducetmp(reducetmp_) {}
665 
666  void operator()(const SizeRange& range) const
667  {
668 
669  const SizeType binSize = arraySize / binCount;
670 
671  // Iterate over bins (array segments)
672  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) {
673  const SizeType begin = n * binSize;
674  const SizeType end = (n == binCount-1) ? arraySize : begin + binSize;
675 
676  // Compute the partial sum for this array segment
677  T sum = zeroVal<T>();
678  for (SizeType i = begin; i < end; ++i) {
679 
680  sum += a[i] * b[i];
681  }
682  // Store the partial sum
683  reducetmp[n] = sum;
684  }
685  }
686 
687 
688  const T* a;
689  const T* b;
690  const SizeType binCount;
691  const SizeType arraySize;
693 };
694 
695 template<typename T>
696 inline T
697 Vector<T>::dot(const Vector<T>& other) const
698 {
699  assert(this->size() == other.size());
700 
701  const T* aData = this->data();
702  const T* bData = other.data();
703 
704  SizeType arraySize = this->size();
705 
706  T result = zeroVal<T>();
707 
708  if (arraySize < 1024) {
709 
710  // Compute the dot product in serial for small arrays
711 
712  for (SizeType n = 0; n < arraySize; ++n) {
713  result += aData[n] * bData[n];
714  }
715 
716  } else {
717 
718  // Compute the dot product by segmenting the arrays into
719  // a predetermined number of sub arrays in parallel and
720  // accumulate the finial result in series.
721 
722  const SizeType binCount = 100;
723  T partialSums[100];
724 
725  tbb::parallel_for(SizeRange(0, binCount),
726  DeterministicDotProductOp(aData, bData, binCount, arraySize, partialSums));
727 
728  for (SizeType n = 0; n < binCount; ++n) {
729  result += partialSums[n];
730  }
731  }
732 
733  return result;
734 }
735 
736 
737 template<typename T>
738 struct Vector<T>::InfNormOp
739 {
740  InfNormOp(const T* data_): data(data_) {}
741 
742  T operator()(const SizeRange& range, T maxValue) const
743  {
744  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) {
745  maxValue = Max(maxValue, Abs(data[n]));
746  }
747  return maxValue;
748  }
749 
750  static T join(T max1, T max2) { return Max(max1, max2); }
751 
752  const T* data;
753 };
754 
755 
756 template<typename T>
757 inline T
759 {
760  // Parallelize over the elements of this vector.
761  T result = tbb::parallel_reduce(SizeRange(0, this->size()), /*seed=*/zeroVal<T>(),
762  InfNormOp(this->data()), InfNormOp::join);
763  return result;
764 }
765 
766 
767 template<typename T>
768 struct Vector<T>::IsFiniteOp
769 {
770  IsFiniteOp(const T* data_): data(data_) {}
771 
772  bool operator()(const SizeRange& range, bool finite) const
773  {
774  if (finite) {
775  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) {
776  if (!std::isfinite(data[n])) return false;
777  }
778  }
779  return finite;
780  }
781 
782  static bool join(bool finite1, bool finite2) { return (finite1 && finite2); }
783 
784  const T* data;
785 };
786 
787 
788 template<typename T>
789 inline bool
791 {
792  // Parallelize over the elements of this vector.
793  bool finite = tbb::parallel_reduce(SizeRange(0, this->size()), /*seed=*/true,
794  IsFiniteOp(this->data()), IsFiniteOp::join);
795  return finite;
796 }
797 
798 
799 template<typename T>
800 template<typename OtherValueType>
801 struct Vector<T>::EqOp
802 {
803  EqOp(const T* a_, const OtherValueType* b_, T e): a(a_), b(b_), eps(e) {}
804 
805  bool operator()(const SizeRange& range, bool equal) const
806  {
807  if (equal) {
808  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) {
809  if (!isApproxEqual(a[n], b[n], eps)) return false;
810  }
811  }
812  return equal;
813  }
814 
815  static bool join(bool eq1, bool eq2) { return (eq1 && eq2); }
816 
817  const T* a;
818  const OtherValueType* b;
819  const T eps;
820 };
821 
822 
823 template<typename T>
824 template<typename OtherValueType>
825 inline bool
826 Vector<T>::eq(const Vector<OtherValueType>& other, ValueType eps) const
827 {
828  if (this->size() != other.size()) return false;
829  bool equal = tbb::parallel_reduce(SizeRange(0, this->size()), /*seed=*/true,
830  EqOp<OtherValueType>(this->data(), other.data(), eps), EqOp<OtherValueType>::join);
831  return equal;
832 }
833 
834 
835 template<typename T>
836 inline std::string
838 {
839  std::ostringstream ostr;
840  ostr << "[";
841  std::string sep;
842  for (SizeType n = 0, N = this->size(); n < N; ++n) {
843  ostr << sep << (*this)[n];
844  sep = ", ";
845  }
846  ostr << "]";
847  return ostr.str();
848 }
849 
850 
852 
853 
854 template<typename ValueType, SizeType STENCIL_SIZE>
855 const ValueType SparseStencilMatrix<ValueType, STENCIL_SIZE>::sZeroValue = zeroVal<ValueType>();
856 
857 
858 template<typename ValueType, SizeType STENCIL_SIZE>
859 inline
861  : mNumRows(numRows)
862  , mValueArray(new ValueType[mNumRows * STENCIL_SIZE])
863  , mColumnIdxArray(new SizeType[mNumRows * STENCIL_SIZE])
864  , mRowSizeArray(new SizeType[mNumRows])
865 {
866  // Initialize the matrix to a null state by setting the size of each row to zero.
867  tbb::parallel_for(SizeRange(0, mNumRows),
868  internal::FillOp<SizeType>(mRowSizeArray.get(), /*value=*/0));
869 }
870 
871 
872 template<typename ValueType, SizeType STENCIL_SIZE>
873 struct SparseStencilMatrix<ValueType, STENCIL_SIZE>::MatrixCopyOp
874 {
876  from(&from_), to(&to_) {}
877 
878  void operator()(const SizeRange& range) const
879  {
880  const ValueType* fromVal = from->mValueArray.get();
881  const SizeType* fromCol = from->mColumnIdxArray.get();
882  ValueType* toVal = to->mValueArray.get();
883  SizeType* toCol = to->mColumnIdxArray.get();
884  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) {
885  toVal[n] = fromVal[n];
886  toCol[n] = fromCol[n];
887  }
888  }
889 
891 };
892 
893 
894 template<typename ValueType, SizeType STENCIL_SIZE>
895 inline
897  : mNumRows(other.mNumRows)
898  , mValueArray(new ValueType[mNumRows * STENCIL_SIZE])
899  , mColumnIdxArray(new SizeType[mNumRows * STENCIL_SIZE])
900  , mRowSizeArray(new SizeType[mNumRows])
901 {
902  SizeType size = mNumRows * STENCIL_SIZE;
903 
904  // Copy the value and column index arrays from the other matrix to this matrix.
905  tbb::parallel_for(SizeRange(0, size), MatrixCopyOp(/*from=*/other, /*to=*/*this));
906 
907  // Copy the row size array from the other matrix to this matrix.
908  tbb::parallel_for(SizeRange(0, mNumRows),
909  internal::CopyOp<SizeType>(/*from=*/other.mRowSizeArray.get(), /*to=*/mRowSizeArray.get()));
910 }
911 
912 
913 template<typename ValueType, SizeType STENCIL_SIZE>
914 inline void
916  const ValueType& val)
917 {
918  assert(row < mNumRows);
919  this->getRowEditor(row).setValue(col, val);
920 }
921 
922 
923 template<typename ValueType, SizeType STENCIL_SIZE>
924 inline const ValueType&
926 {
927  assert(row < mNumRows);
928  return this->getConstRow(row).getValue(col);
929 }
930 
931 
932 template<typename ValueType, SizeType STENCIL_SIZE>
933 inline const ValueType&
935 {
936  return this->getValue(row,col);
937 }
938 
939 
940 template<typename ValueType, SizeType STENCIL_SIZE>
941 template<typename Scalar>
942 struct SparseStencilMatrix<ValueType, STENCIL_SIZE>::RowScaleOp
943 {
944  RowScaleOp(SparseStencilMatrix& m, const Scalar& s_): mat(&m), s(s_) {}
945 
946  void operator()(const SizeRange& range) const
947  {
948  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) {
949  RowEditor row = mat->getRowEditor(n);
950  row.scale(s);
951  }
952  }
953 
954  SparseStencilMatrix* mat;
955  const Scalar s;
956 };
957 
958 
959 template<typename ValueType, SizeType STENCIL_SIZE>
960 template<typename Scalar>
961 inline void
963 {
964  // Parallelize over the rows in the matrix.
965  tbb::parallel_for(SizeRange(0, mNumRows), RowScaleOp<Scalar>(*this, s));
966 }
967 
968 
969 template<typename ValueType, SizeType STENCIL_SIZE>
970 template<typename VecValueType>
971 struct SparseStencilMatrix<ValueType, STENCIL_SIZE>::VecMultOp
972 {
973  VecMultOp(const SparseStencilMatrix& m, const VecValueType* i, VecValueType* o):
974  mat(&m), in(i), out(o) {}
975 
976  void operator()(const SizeRange& range) const
977  {
978  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) {
979  ConstRow row = mat->getConstRow(n);
980  out[n] = row.dot(in, mat->numRows());
981  }
982  }
983 
984  const SparseStencilMatrix* mat;
985  const VecValueType* in;
986  VecValueType* out;
987 };
988 
989 
990 template<typename ValueType, SizeType STENCIL_SIZE>
991 template<typename VecValueType>
992 inline void
994  const Vector<VecValueType>& inVec, Vector<VecValueType>& resultVec) const
995 {
996  if (inVec.size() != mNumRows) {
997  OPENVDB_THROW(ArithmeticError, "matrix and input vector have incompatible sizes ("
998  << mNumRows << "x" << mNumRows << " vs. " << inVec.size() << ")");
999  }
1000  if (resultVec.size() != mNumRows) {
1001  OPENVDB_THROW(ArithmeticError, "matrix and result vector have incompatible sizes ("
1002  << mNumRows << "x" << mNumRows << " vs. " << resultVec.size() << ")");
1003  }
1004 
1005  vectorMultiply(inVec.data(), resultVec.data());
1006 }
1007 
1008 
1009 template<typename ValueType, SizeType STENCIL_SIZE>
1010 template<typename VecValueType>
1011 inline void
1013  const VecValueType* inVec, VecValueType* resultVec) const
1014 {
1015  // Parallelize over the rows in the matrix.
1016  tbb::parallel_for(SizeRange(0, mNumRows),
1017  VecMultOp<VecValueType>(*this, inVec, resultVec));
1018 }
1019 
1020 
1021 template<typename ValueType, SizeType STENCIL_SIZE>
1022 template<typename OtherValueType>
1023 struct SparseStencilMatrix<ValueType, STENCIL_SIZE>::EqOp
1024 {
1025  EqOp(const SparseStencilMatrix& a_,
1026  const SparseStencilMatrix<OtherValueType, STENCIL_SIZE>& b_, ValueType e):
1027  a(&a_), b(&b_), eps(e) {}
1028 
1029  bool operator()(const SizeRange& range, bool equal) const
1030  {
1031  if (equal) {
1032  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) {
1033  if (!a->getConstRow(n).eq(b->getConstRow(n), eps)) return false;
1034  }
1035  }
1036  return equal;
1037  }
1038 
1039  static bool join(bool eq1, bool eq2) { return (eq1 && eq2); }
1040 
1041  const SparseStencilMatrix* a;
1042  const SparseStencilMatrix<OtherValueType, STENCIL_SIZE>* b;
1043  const ValueType eps;
1044 };
1045 
1046 
1047 template<typename ValueType, SizeType STENCIL_SIZE>
1048 template<typename OtherValueType>
1049 inline bool
1051  const SparseStencilMatrix<OtherValueType, STENCIL_SIZE>& other, ValueType eps) const
1052 {
1053  if (this->numRows() != other.numRows()) return false;
1054  bool equal = tbb::parallel_reduce(SizeRange(0, this->numRows()), /*seed=*/true,
1055  EqOp<OtherValueType>(*this, other, eps), EqOp<OtherValueType>::join);
1056  return equal;
1057 }
1058 
1059 
1060 template<typename ValueType, SizeType STENCIL_SIZE>
1061 struct SparseStencilMatrix<ValueType, STENCIL_SIZE>::IsFiniteOp
1062 {
1063  IsFiniteOp(const SparseStencilMatrix& m): mat(&m) {}
1064 
1065  bool operator()(const SizeRange& range, bool finite) const
1066  {
1067  if (finite) {
1068  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) {
1069  const ConstRow row = mat->getConstRow(n);
1070  for (ConstValueIter it = row.cbegin(); it; ++it) {
1071  if (!std::isfinite(*it)) return false;
1072  }
1073  }
1074  }
1075  return finite;
1076  }
1077 
1078  static bool join(bool finite1, bool finite2) { return (finite1 && finite2); }
1079 
1081 };
1082 
1083 
1084 template<typename ValueType, SizeType STENCIL_SIZE>
1085 inline bool
1087 {
1088  // Parallelize over the rows of this matrix.
1089  bool finite = tbb::parallel_reduce(SizeRange(0, this->numRows()), /*seed=*/true,
1090  IsFiniteOp(*this), IsFiniteOp::join);
1091  return finite;
1092 }
1093 
1094 
1095 template<typename ValueType, SizeType STENCIL_SIZE>
1096 inline std::string
1098 {
1099  std::ostringstream ostr;
1100  for (SizeType n = 0, N = this->size(); n < N; ++n) {
1101  ostr << n << ": " << this->getConstRow(n).str() << "\n";
1102  }
1103  return ostr.str();
1104 }
1105 
1106 
1107 template<typename ValueType, SizeType STENCIL_SIZE>
1110 {
1111  assert(i < mNumRows);
1112  const SizeType head = i * STENCIL_SIZE;
1113  return RowEditor(&mValueArray[head], &mColumnIdxArray[head], mRowSizeArray[i], mNumRows);
1114 }
1115 
1116 
1117 template<typename ValueType, SizeType STENCIL_SIZE>
1120 {
1121  assert(i < mNumRows);
1122  const SizeType head = i * STENCIL_SIZE; // index for this row into main storage
1123  return ConstRow(&mValueArray[head], &mColumnIdxArray[head], mRowSizeArray[i]);
1124 }
1125 
1126 
1127 template<typename ValueType, SizeType STENCIL_SIZE>
1128 template<typename DataType>
1129 inline SizeType
1131 {
1132  if (this->empty()) return mData.mSize;
1133 
1134  // Get a pointer to the first column index that is equal to or greater than the given index.
1135  // (This assumes that the data is sorted by column.)
1136  const SizeType* colPtr = std::lower_bound(mData.mCols, mData.mCols + mData.mSize, columnIdx);
1137  // Return the offset of the pointer from the beginning of the array.
1138  return static_cast<SizeType>(colPtr - mData.mCols);
1139 }
1140 
1141 
1142 template<typename ValueType, SizeType STENCIL_SIZE>
1143 template<typename DataType>
1144 inline const ValueType&
1146  SizeType columnIdx, bool& active) const
1147 {
1148  active = false;
1149  SizeType idx = this->find(columnIdx);
1150  if (idx < this->size() && this->column(idx) == columnIdx) {
1151  active = true;
1152  return this->value(idx);
1153  }
1155 }
1156 
1157 template<typename ValueType, SizeType STENCIL_SIZE>
1158 template<typename DataType>
1159 inline const ValueType&
1161 {
1162  SizeType idx = this->find(columnIdx);
1163  if (idx < this->size() && this->column(idx) == columnIdx) {
1164  return this->value(idx);
1165  }
1167 }
1168 
1169 
1170 template<typename ValueType, SizeType STENCIL_SIZE>
1171 template<typename DataType>
1172 inline typename SparseStencilMatrix<ValueType, STENCIL_SIZE>::ConstValueIter
1173 SparseStencilMatrix<ValueType, STENCIL_SIZE>::RowBase<DataType>::cbegin() const
1174 {
1175  return ConstValueIter(mData);
1176 }
1177 
1178 
1179 template<typename ValueType, SizeType STENCIL_SIZE>
1180 template<typename DataType>
1181 template<typename OtherDataType>
1182 inline bool
1184  const RowBase<OtherDataType>& other, ValueType eps) const
1185 {
1186  if (this->size() != other.size()) return false;
1187  for (ConstValueIter it = cbegin(), oit = other.cbegin(); it || oit; ++it, ++oit) {
1188  if (it.column() != oit.column()) return false;
1189  if (!isApproxEqual(*it, *oit, eps)) return false;
1190  }
1191  return true;
1192 }
1193 
1194 
1195 template<typename ValueType, SizeType STENCIL_SIZE>
1196 template<typename DataType>
1197 template<typename VecValueType>
1198 inline VecValueType
1199 SparseStencilMatrix<ValueType, STENCIL_SIZE>::RowBase<DataType>::dot(
1200  const VecValueType* inVec, SizeType vecSize) const
1201 {
1202  VecValueType result = zeroVal<VecValueType>();
1203  for (SizeType idx = 0, N = std::min(vecSize, this->size()); idx < N; ++idx) {
1204  result += static_cast<VecValueType>(this->value(idx) * inVec[this->column(idx)]);
1205  }
1206  return result;
1207 }
1208 
1209 template<typename ValueType, SizeType STENCIL_SIZE>
1210 template<typename DataType>
1211 template<typename VecValueType>
1212 inline VecValueType
1213 SparseStencilMatrix<ValueType, STENCIL_SIZE>::RowBase<DataType>::dot(
1214  const Vector<VecValueType>& inVec) const
1215 {
1216  return dot(inVec.data(), inVec.size());
1217 }
1218 
1219 
1220 template<typename ValueType, SizeType STENCIL_SIZE>
1221 template<typename DataType>
1222 inline std::string
1224 {
1225  std::ostringstream ostr;
1226  std::string sep;
1227  for (SizeType n = 0, N = this->size(); n < N; ++n) {
1228  ostr << sep << "(" << this->column(n) << ", " << this->value(n) << ")";
1229  sep = ", ";
1230  }
1231  return ostr.str();
1232 }
1233 
1234 
1235 template<typename ValueType, SizeType STENCIL_SIZE>
1236 inline
1238  const ValueType* valueHead, const SizeType* columnHead, const SizeType& rowSize):
1239  ConstRowBase(ConstRowData(const_cast<ValueType*>(valueHead),
1240  const_cast<SizeType*>(columnHead), const_cast<SizeType&>(rowSize)))
1241 {
1242 }
1243 
1244 
1245 template<typename ValueType, SizeType STENCIL_SIZE>
1246 inline
1248  ValueType* valueHead, SizeType* columnHead, SizeType& rowSize, SizeType colSize):
1249  RowBase<>(RowData(valueHead, columnHead, rowSize)), mNumColumns(colSize)
1250 {
1251 }
1252 
1253 
1254 template<typename ValueType, SizeType STENCIL_SIZE>
1255 inline void
1257 {
1258  // Note: since mSize is a reference, this modifies the underlying matrix.
1259  RowBase<>::mData.mSize = 0;
1260 }
1261 
1262 
1263 template<typename ValueType, SizeType STENCIL_SIZE>
1264 inline SizeType
1266  SizeType column, const ValueType& value)
1267 {
1268  assert(column < mNumColumns);
1269 
1270  RowData& data = RowBase<>::mData;
1271 
1272  // Get the offset of the first column index that is equal to or greater than
1273  // the column to be modified.
1274  SizeType offset = this->find(column);
1275 
1276  if (offset < data.mSize && data.mCols[offset] == column) {
1277  // If the column already exists, just update its value.
1278  data.mVals[offset] = value;
1279  return data.mSize;
1280  }
1281 
1282  // Check that it is safe to add a new column.
1283  assert(data.mSize < this->capacity());
1284 
1285  if (offset >= data.mSize) {
1286  // The new column's index is larger than any existing index. Append the new column.
1287  data.mVals[data.mSize] = value;
1288  data.mCols[data.mSize] = column;
1289  } else {
1290  // Insert the new column at the computed offset after shifting subsequent columns.
1291  for (SizeType i = data.mSize; i > offset; --i) {
1292  data.mVals[i] = data.mVals[i - 1];
1293  data.mCols[i] = data.mCols[i - 1];
1294  }
1295  data.mVals[offset] = value;
1296  data.mCols[offset] = column;
1297  }
1298  ++data.mSize;
1299 
1300  return data.mSize;
1301 }
1302 
1303 
1304 template<typename ValueType, SizeType STENCIL_SIZE>
1305 template<typename Scalar>
1306 inline void
1308 {
1309  for (int idx = 0, N = this->size(); idx < N; ++idx) {
1310  RowBase<>::mData.mVals[idx] *= s;
1311  }
1312 }
1313 
1314 
1316 
1317 
1319 template<typename MatrixType>
1320 class JacobiPreconditioner: public Preconditioner<typename MatrixType::ValueType>
1321 {
1322 private:
1323  struct InitOp;
1324  struct ApplyOp;
1325 
1326 public:
1327  typedef typename MatrixType::ValueType ValueType;
1330  typedef boost::shared_ptr<JacobiPreconditioner> Ptr;
1331 
1332  JacobiPreconditioner(const MatrixType& A): BaseType(A), mDiag(A.numRows())
1333  {
1334  // Initialize vector mDiag with the values from the matrix diagonal.
1335  tbb::parallel_for(SizeRange(0, A.numRows()), InitOp(A, mDiag.data()));
1336  }
1337 
1339 
1340  virtual void apply(const Vector<ValueType>& r, Vector<ValueType>& z)
1341  {
1342  const SizeType size = mDiag.size();
1343 
1344  assert(r.size() == z.size());
1345  assert(r.size() == size);
1346 
1347  tbb::parallel_for(SizeRange(0, size), ApplyOp(mDiag.data(), r.data(), z.data()));
1348  }
1349 
1351  bool isFinite() const { return mDiag.isFinite(); }
1352 
1353 private:
1354  // Functor for use with tbb::parallel_for()
1355  struct InitOp
1356  {
1357  InitOp(const MatrixType& m, ValueType* v): mat(&m), vec(v) {}
1358  void operator()(const SizeRange& range) const {
1359  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) {
1360  const ValueType val = mat->getValue(n, n);
1361  assert(!isApproxZero(val, ValueType(0.0001)));
1362  vec[n] = static_cast<ValueType>(1.0 / val);
1363  }
1364  }
1365  const MatrixType* mat; ValueType* vec;
1366  };
1367 
1368  // Functor for use with tbb::parallel_reduce()
1369  struct ApplyOp
1370  {
1371  ApplyOp(const ValueType* x_, const ValueType* y_, ValueType* out_):
1372  x(x_), y(y_), out(out_) {}
1373  void operator()(const SizeRange& range) const {
1374  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) out[n] = x[n] * y[n];
1375  }
1376  const ValueType *x, *y; ValueType* out;
1377  };
1378 
1379  // The Jacobi preconditioner is a diagonal matrix
1380  VectorType mDiag;
1381 }; // class JacobiPreconditioner
1382 
1383 
1385 template<typename MatrixType>
1386 class IncompleteCholeskyPreconditioner: public Preconditioner<typename MatrixType::ValueType>
1387 {
1388 private:
1389  struct CopyToLowerOp;
1390  struct TransposeOp;
1391 
1392 public:
1393  typedef typename MatrixType::ValueType ValueType;
1396  typedef boost::shared_ptr<IncompleteCholeskyPreconditioner> Ptr;
1398  typedef typename TriangularMatrix::ConstRow TriangleConstRow;
1399  typedef typename TriangularMatrix::RowEditor TriangleRowEditor;
1400 
1401  IncompleteCholeskyPreconditioner(const MatrixType& matrix)
1402  : BaseType(matrix)
1403  , mLowerTriangular(matrix.numRows())
1404  , mUpperTriangular(matrix.numRows())
1405  , mTempVec(matrix.numRows())
1406  {
1407  // Size of matrix
1408  const SizeType numRows = mLowerTriangular.numRows();
1409 
1410  // Copy the upper triangular part to the lower triangular part.
1411  tbb::parallel_for(SizeRange(0, numRows), CopyToLowerOp(matrix, mLowerTriangular));
1412 
1413  // Build the Incomplete Cholesky Matrix
1414  //
1415  // Algorithm:
1416  //
1417  // for (k = 0; k < size; ++k) {
1418  // A(k,k) = sqrt(A(k,k));
1419  // for (i = k +1, i < size; ++i) {
1420  // if (A(i,k) == 0) continue;
1421  // A(i,k) = A(i,k) / A(k,k);
1422  // }
1423  // for (j = k+1; j < size; ++j) {
1424  // for (i = j; i < size; ++i) {
1425  // if (A(i,j) == 0) continue;
1426  // A(i,j) -= A(i,k)*A(j,k);
1427  // }
1428  // }
1429  // }
1430 
1431  mPassedCompatibilityCondition = true;
1432 
1433  for (SizeType k = 0; k < numRows; ++k) {
1434 
1435  TriangleConstRow crow_k = mLowerTriangular.getConstRow(k);
1436  ValueType diagonalValue = crow_k.getValue(k);
1437 
1438  // Test if the matrix build has failed.
1439  if (diagonalValue < 1.e-5) {
1440  mPassedCompatibilityCondition = false;
1441  break;
1442  }
1443 
1444  diagonalValue = Sqrt(diagonalValue);
1445 
1446  TriangleRowEditor row_k = mLowerTriangular.getRowEditor(k);
1447  row_k.setValue(k, diagonalValue);
1448 
1449  // Exploit the fact that the matrix is symmetric.
1450  typename MatrixType::ConstRow srcRow = matrix.getConstRow(k);
1451  typename MatrixType::ConstValueIter citer = srcRow.cbegin();
1452  for ( ; citer; ++citer) {
1453  SizeType ii = citer.column();
1454  if (ii < k+1) continue; // look above diagonal
1455 
1456  TriangleRowEditor row_ii = mLowerTriangular.getRowEditor(ii);
1457 
1458  row_ii.setValue(k, *citer / diagonalValue);
1459  }
1460 
1461  // for (j = k+1; j < size; ++j) replaced by row iter below
1462  citer.reset(); // k,j entries
1463  for ( ; citer; ++citer) {
1464  SizeType j = citer.column();
1465  if (j < k+1) continue;
1466 
1467  TriangleConstRow row_j = mLowerTriangular.getConstRow(j);
1468  ValueType a_jk = row_j.getValue(k); // a_jk is non zero if a_kj is non zero
1469 
1470  // Entry (i,j) is non-zero if matrix(j,i) is nonzero
1471 
1472  typename MatrixType::ConstRow mask = matrix.getConstRow(j);
1473  typename MatrixType::ConstValueIter maskIter = mask.cbegin();
1474  for ( ; maskIter; ++maskIter) {
1475  SizeType i = maskIter.column();
1476  if (i < j) continue;
1477 
1478  TriangleConstRow crow_i = mLowerTriangular.getConstRow(i);
1479  ValueType a_ij = crow_i.getValue(j);
1480  ValueType a_ik = crow_i.getValue(k);
1481  TriangleRowEditor row_i = mLowerTriangular.getRowEditor(i);
1482  a_ij -= a_ik * a_jk;
1483 
1484  row_i.setValue(j, a_ij);
1485  }
1486  }
1487  }
1488 
1489  // Build the transpose of the IC matrix: mUpperTriangular
1490  tbb::parallel_for(SizeRange(0, numRows),
1491  TransposeOp(matrix, mLowerTriangular, mUpperTriangular));
1492  }
1493 
1495 
1496  virtual bool isValid() const { return mPassedCompatibilityCondition; }
1497 
1498  virtual void apply(const Vector<ValueType>& rVec, Vector<ValueType>& zVec)
1499  {
1500  if (!mPassedCompatibilityCondition) {
1501  OPENVDB_THROW(ArithmeticError, "invalid Cholesky decomposition");
1502  }
1503 
1504  // Solve mUpperTriangular * mLowerTriangular * rVec = zVec;
1505 
1506  SizeType size = mLowerTriangular.numRows();
1507 
1508  zVec.fill(zeroVal<ValueType>());
1509  ValueType* zData = zVec.data();
1510 
1511  if (size == 0) return;
1512 
1513  assert(rVec.size() == size);
1514  assert(zVec.size() == size);
1515 
1516  // Allocate a temp vector
1517  mTempVec.fill(zeroVal<ValueType>());
1518  ValueType* tmpData = mTempVec.data();
1519  const ValueType* rData = rVec.data();
1520 
1521  // Solve mLowerTriangular * tmp = rVec;
1522  for (SizeType i = 0; i < size; ++i) {
1523  typename TriangularMatrix::ConstRow row = mLowerTriangular.getConstRow(i);
1524  ValueType diagonal = row.getValue(i);
1525  ValueType dot = row.dot(mTempVec);
1526  tmpData[i] = (rData[i] - dot) / diagonal;
1527  if (!std::isfinite(tmpData[i])) {
1528  OPENVDB_LOG_DEBUG_RUNTIME("1 diagonal was " << diagonal);
1529  OPENVDB_LOG_DEBUG_RUNTIME("1a diagonal " << row.getValue(i));
1530  }
1531  }
1532 
1533  // Solve mUpperTriangular * zVec = tmp;
1534  for (SizeType ii = 0; ii < size; ++ii) {
1535  SizeType i = size - 1 - ii;
1536  typename TriangularMatrix::ConstRow row = mUpperTriangular.getConstRow(i);
1537  ValueType diagonal = row.getValue(i);
1538  ValueType dot = row.dot(zVec);
1539  zData[i] = (tmpData[i] - dot) / diagonal;
1540  if (!std::isfinite(zData[i])) {
1541  OPENVDB_LOG_DEBUG_RUNTIME("2 diagonal was " << diagonal);
1542  }
1543  }
1544  }
1545 
1546  const TriangularMatrix& lowerMatrix() const { return mLowerTriangular; }
1547  const TriangularMatrix& upperMatrix() const { return mUpperTriangular; }
1548 
1549 private:
1550  // Functor for use with tbb::parallel_for()
1551  struct CopyToLowerOp
1552  {
1553  CopyToLowerOp(const MatrixType& m, TriangularMatrix& l): mat(&m), lower(&l) {}
1554  void operator()(const SizeRange& range) const {
1555  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) {
1556  typename TriangularMatrix::RowEditor outRow = lower->getRowEditor(n);
1557  outRow.clear();
1558  typename MatrixType::ConstRow inRow = mat->getConstRow(n);
1559  for (typename MatrixType::ConstValueIter it = inRow.cbegin(); it; ++it) {
1560  if (it.column() > n) continue; // skip above diagonal
1561  outRow.setValue(it.column(), *it);
1562  }
1563  }
1564  }
1565  const MatrixType* mat; TriangularMatrix* lower;
1566  };
1567 
1568  // Functor for use with tbb::parallel_for()
1569  struct TransposeOp
1570  {
1571  TransposeOp(const MatrixType& m, const TriangularMatrix& l, TriangularMatrix& u):
1572  mat(&m), lower(&l), upper(&u) {}
1573  void operator()(const SizeRange& range) const {
1574  for (SizeType n = range.begin(), N = range.end(); n < N; ++n) {
1575  typename TriangularMatrix::RowEditor outRow = upper->getRowEditor(n);
1576  outRow.clear();
1577  // Use the fact that matrix is symmetric.
1578  typename MatrixType::ConstRow inRow = mat->getConstRow(n);
1579  for (typename MatrixType::ConstValueIter it = inRow.cbegin(); it; ++it) {
1580  const SizeType column = it.column();
1581  if (column < n) continue; // only set upper triangle
1582  outRow.setValue(column, lower->getValue(column, n));
1583  }
1584  }
1585  }
1586  const MatrixType* mat; const TriangularMatrix* lower; TriangularMatrix* upper;
1587  };
1588 
1589  TriangularMatrix mLowerTriangular;
1590  TriangularMatrix mUpperTriangular;
1591  Vector<ValueType> mTempVec;
1592  bool mPassedCompatibilityCondition;
1593 }; // class IncompleteCholeskyPreconditioner
1594 
1595 
1597 
1598 
1599 namespace internal {
1600 
1602 template<typename T>
1603 inline void
1604 axpy(const T& a, const T* xVec, const T* yVec, T* resultVec, SizeType size)
1605 {
1606  tbb::parallel_for(SizeRange(0, size), LinearOp<T>(a, xVec, yVec, resultVec));
1607 }
1608 
1610 template<typename T>
1611 inline void
1612 axpy(const T& a, const Vector<T>& xVec, const Vector<T>& yVec, Vector<T>& result)
1613 {
1614  assert(xVec.size() == yVec.size());
1615  assert(xVec.size() == result.size());
1616  axpy(a, xVec.data(), yVec.data(), result.data(), xVec.size());
1617 }
1618 
1619 
1621 template<typename MatrixOperator, typename VecValueType>
1622 inline void
1623 computeResidual(const MatrixOperator& A, const VecValueType* x,
1624  const VecValueType* b, VecValueType* r)
1625 {
1626  // Compute r = A * x.
1627  A.vectorMultiply(x, r);
1628  // Compute r = b - r.
1629  tbb::parallel_for(SizeRange(0, A.numRows()), LinearOp<VecValueType>(-1.0, r, b, r));
1630 }
1631 
1633 template<typename MatrixOperator, typename T>
1634 inline void
1635 computeResidual(const MatrixOperator& A, const Vector<T>& x, const Vector<T>& b, Vector<T>& r)
1636 {
1637  assert(x.size() == b.size());
1638  assert(x.size() == r.size());
1639  assert(x.size() == A.numRows());
1640 
1641  computeResidual(A, x.data(), b.data(), r.data());
1642 }
1643 
1644 } // namespace internal
1645 
1646 
1648 
1649 
1650 template<typename PositiveDefMatrix>
1651 inline State
1653  const PositiveDefMatrix& Amat,
1657  const State& termination)
1658 {
1659  util::NullInterrupter interrupter;
1660  return solve(Amat, bVec, xVec, precond, interrupter, termination);
1661 }
1662 
1663 
1664 template<typename PositiveDefMatrix, typename Interrupter>
1665 inline State
1667  const PositiveDefMatrix& Amat,
1671  Interrupter& interrupter,
1672  const State& termination)
1673 {
1674  typedef typename PositiveDefMatrix::ValueType ValueType;
1675  typedef Vector<ValueType> VectorType;
1676 
1677  State result;
1678  result.success = false;
1679  result.iterations = 0;
1680  result.relativeError = 0.0;
1681  result.absoluteError = 0.0;
1682 
1683  const SizeType size = Amat.numRows();
1684  if (size == 0) {
1685  OPENVDB_LOG_WARN("pcg::solve(): matrix has dimension zero");
1686  return result;
1687  }
1688  if (size != bVec.size()) {
1689  OPENVDB_THROW(ArithmeticError, "A and b have incompatible sizes"
1690  << size << "x" << size << " vs. " << bVec.size() << ")");
1691  }
1692  if (size != xVec.size()) {
1693  OPENVDB_THROW(ArithmeticError, "A and x have incompatible sizes"
1694  << size << "x" << size << " vs. " << xVec.size() << ")");
1695  }
1696 
1697  // Temp vectors
1698  VectorType zVec(size); // transformed residual (M^-1 r)
1699  VectorType pVec(size); // search direction
1700  VectorType qVec(size); // A * p
1701 
1702  // Compute norm of B (the source)
1703  const ValueType tmp = bVec.infNorm();
1704  const ValueType infNormOfB = isZero(tmp) ? 1.f : tmp;
1705 
1706  // Compute rVec: residual = b - Ax.
1707  VectorType rVec(size); // vector of residuals
1708 
1709  internal::computeResidual(Amat, xVec, bVec, rVec);
1710 
1711  assert(rVec.isFinite());
1712 
1713  // Normalize the residual norm with the source norm and look for early out.
1714  result.absoluteError = static_cast<double>(rVec.infNorm());
1715  result.relativeError = static_cast<double>(result.absoluteError / infNormOfB);
1716  if (result.relativeError <= termination.relativeError) {
1717  result.success = true;
1718  return result;
1719  }
1720 
1721  // Iterations of the CG solve
1722 
1723  ValueType rDotZPrev(1); // inner product of <z,r>
1724 
1725  // Keep track of the minimum error to monitor convergence.
1726  ValueType minL2Error = std::numeric_limits<ValueType>::max();
1727  ValueType l2Error;
1728 
1729  int iteration = 0;
1730  for ( ; iteration < termination.iterations; ++iteration) {
1731 
1732  if (interrupter.wasInterrupted()) {
1733  OPENVDB_THROW(RuntimeError, "conjugate gradient solver was interrupted");
1734  }
1735 
1736  OPENVDB_LOG_DEBUG_RUNTIME("pcg::solve() " << result);
1737 
1738  result.iterations = iteration + 1;
1739 
1740  // Apply preconditioner to residual
1741  // z_{k} = M^-1 r_{k}
1742  precond.apply(rVec, zVec);
1743 
1744  // <r,z>
1745  const ValueType rDotZ = rVec.dot(zVec);
1746  assert(std::isfinite(rDotZ));
1747 
1748  if (0 == iteration) {
1749  // Initialize
1750  pVec = zVec;
1751  } else {
1752  const ValueType beta = rDotZ / rDotZPrev;
1753  // p = beta * p + z
1754  internal::axpy(beta, pVec, zVec, /*result */pVec);
1755  }
1756 
1757  // q_{k} = A p_{k}
1758  Amat.vectorMultiply(pVec, qVec);
1759 
1760  // alpha = <r_{k-1}, z_{k-1}> / <p_{k},q_{k}>
1761  const ValueType pAp = pVec.dot(qVec);
1762  assert(std::isfinite(pAp));
1763 
1764  const ValueType alpha = rDotZ / pAp;
1765  rDotZPrev = rDotZ;
1766 
1767  // x_{k} = x_{k-1} + alpha * p_{k}
1768  internal::axpy(alpha, pVec, xVec, /*result=*/xVec);
1769 
1770  // r_{k} = r_{k-1} - alpha_{k-1} A p_{k}
1771  internal::axpy(-alpha, qVec, rVec, /*result=*/rVec);
1772 
1773  // update tolerances
1774  l2Error = rVec.l2Norm();
1775  minL2Error = Min(l2Error, minL2Error);
1776 
1777  result.absoluteError = static_cast<double>(rVec.infNorm());
1778  result.relativeError = static_cast<double>(result.absoluteError / infNormOfB);
1779 
1780  if (l2Error > 2 * minL2Error) {
1781  // The solution started to diverge.
1782  result.success = false;
1783  break;
1784  }
1785  if (!std::isfinite(result.absoluteError)) {
1786  // Total divergence of solution
1787  result.success = false;
1788  break;
1789  }
1790  if (result.absoluteError <= termination.absoluteError) {
1791  // Convergence
1792  result.success = true;
1793  break;
1794  }
1795  if (result.relativeError <= termination.relativeError) {
1796  // Convergence
1797  result.success = true;
1798  break;
1799  }
1800  }
1801  OPENVDB_LOG_DEBUG_RUNTIME("pcg::solve() " << result);
1802 
1803  return result;
1804 }
1805 
1806 } // namespace pcg
1807 } // namespace math
1808 } // namespace OPENVDB_VERSION_NAME
1809 } // namespace openvdb
1810 
1811 #endif // OPENVDB_MATH_CONJGRADIENT_HAS_BEEN_INCLUDED
1812 
1814 //
1815 // Copyright (c) 2012-2015 DreamWorks Animation LLC
1816 //
1817 // All rights reserved. This software is distributed under the
1818 // Mozilla Public License 2.0 ( http://www.mozilla.org/MPL/2.0/ )
1819 //
1820 // Redistributions of source code must retain the above copyright
1821 // and license notice and the following restrictions and disclaimer.
1822 //
1823 // * Neither the name of DreamWorks Animation nor the names of
1824 // its contributors may be used to endorse or promote products derived
1825 // from this software without specific prior written permission.
1826 //
1827 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
1828 // "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
1829 // LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
1830 // A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
1831 // OWNER OR CONTRIBUTORS BE LIABLE FOR ANY INDIRECT, INCIDENTAL,
1832 // SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
1833 // LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
1834 // DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
1835 // THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
1836 // (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
1837 // OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
1838 // IN NO EVENT SHALL THE COPYRIGHT HOLDERS' AND CONTRIBUTORS' AGGREGATE
1839 // LIABILITY FOR ALL CLAIMS REGARDLESS OF THEIR BASIS EXCEED US$250.00.
1840 //
const ValueType & operator*() const
Definition: ConjGradient.h:415
Preconditioner(const SparseStencilMatrix< T, STENCIL_SIZE > &)
Definition: ConjGradient.h:499
Vector(SizeType n, const ValueType &val)
Construct a vector of n elements and initialize each element to the given value.
Definition: ConjGradient.h:177
void axpy(const T &a, const Vector< T > &xVec, const Vector< T > &yVec, Vector< T > &result)
Compute ax + y.
Definition: ConjGradient.h:1612
LinearOp(const T &a_, const T *x_, const T *y_, T *out_)
Definition: ConjGradient.h:551
virtual ~JacobiPreconditioner()
Definition: ConjGradient.h:1338
JacobiPreconditioner(const MatrixType &A)
Definition: ConjGradient.h:1332
T * out
Definition: ConjGradient.h:564
Definition: ConjGradient.h:519
boost::shared_ptr< SparseStencilMatrix > Ptr
Definition: ConjGradient.h:271
State terminationDefaults()
Return default termination conditions for a conjugate gradient solver.
Definition: ConjGradient.h:85
Preconditioner< ValueType > BaseType
Definition: ConjGradient.h:1394
void scale(const Scalar &s)
Multiply all elements in the matrix by s;.
Definition: ConjGradient.h:962
Vector< float > VectorS
Definition: ConjGradient.h:255
virtual void apply(const Vector< ValueType > &rVec, Vector< ValueType > &zVec)
Apply this preconditioner to a residue vector: z = M−1r
Definition: ConjGradient.h:1498
Base class for conjugate gradient preconditioners.
Definition: ConjGradient.h:69
const ValueType & getValue(SizeType row, SizeType col) const
Return the value at the given coordinates.
Definition: ConjGradient.h:925
static bool join(bool finite1, bool finite2)
Definition: ConjGradient.h:782
virtual bool isValid() const
Definition: ConjGradient.h:502
TriangularMatrix::RowEditor TriangleRowEditor
Definition: ConjGradient.h:1399
General-purpose arithmetic and comparison routines, most of which accept arbitrary value types (or at...
SparseStencilMatrix & operator*=(const Scalar &s)
Multiply all elements in the matrix by s;.
Definition: ConjGradient.h:314
int32_t Abs(int32_t i)
Return the absolute value of the given quantity.
Definition: Math.h:293
#define OPENVDB_LOG_DEBUG_RUNTIME(message)
Log a debugging message in both debug and optimized builds.
Definition: logging.h:48
Definition: ConjGradient.h:534
static const ValueType sZeroValue
Definition: ConjGradient.h:275
Preconditioner< ValueType > BaseType
Definition: ConjGradient.h:1328
IncompleteCholeskyPreconditioner(const MatrixType &matrix)
Definition: ConjGradient.h:1401
MatrixType::ValueType ValueType
Definition: ConjGradient.h:1390
bool isExactlyEqual(const T0 &a, const T1 &b)
Return true if a is exactly equal to b.
Definition: Math.h:407
Vector< ValueType > VectorType
Definition: ConjGradient.h:1395
bool isFinite() const
Return true if every element of this vector has a finite value.
Definition: ConjGradient.h:790
void scale(const Scalar &)
Scale all of the entries in this row.
Definition: ConjGradient.h:1307
Read/write accessor to a row of this matrix.
Definition: ConjGradient.h:448
static bool join(bool finite1, bool finite2)
Definition: ConjGradient.h:1078
#define OPENVDB_THROW(exception, message)
Definition: Exceptions.h:97
const T * constData() const
Return a pointer to this vector's elements.
Definition: ConjGradient.h:239
SizeType column() const
Definition: ConjGradient.h:421
Preconditioner using incomplete Cholesky factorization.
Definition: ConjGradient.h:71
void computeResidual(const MatrixOperator &A, const VecValueType *x, const VecValueType *b, VecValueType *r)
Compute r = b − Ax.
Definition: ConjGradient.h:1623
MatrixCopyOp(const SparseStencilMatrix &from_, SparseStencilMatrix &to_)
Definition: ConjGradient.h:875
uint32_t Index32
Definition: Types.h:56
CopyOp(const T *from_, T *to_)
Definition: ConjGradient.h:521
const T & operator[](SizeType i) const
Return the value of this vector's ith element.
Definition: ConjGradient.h:232
Vector()
Construct an empty vector.
Definition: ConjGradient.h:173
Vector & operator=(const Vector &)
Deep copy the given vector.
Definition: ConjGradient.h:597
const TriangularMatrix & lowerMatrix() const
Definition: ConjGradient.h:1546
void clear()
Set the number of entries in this row to zero.
Definition: ConjGradient.h:1256
std::ostream & operator<<(std::ostream &os, const State &state)
Definition: ConjGradient.h:574
Definition: ConjGradient.h:738
const T * data
Definition: ConjGradient.h:752
bool operator()(const SizeRange &range, bool finite) const
Definition: ConjGradient.h:772
Diagonal preconditioner.
Definition: ConjGradient.h:70
void vectorMultiply(const Vector< VecValueType > &inVec, Vector< VecValueType > &resultVec) const
Multiply this matrix by inVec and return the result in resultVec.
Definition: ConjGradient.h:993
static T join(T max1, T max2)
Definition: ConjGradient.h:750
FillOp(T *data_, const T &val_)
Definition: ConjGradient.h:536
bool isFinite() const
Return true if every element of this matrix has a finite value.
Definition: ConjGradient.h:1086
bool isApproxEqual(const Type &a, const Type &b)
Return true if a is equal to b to within the default floating-point comparison tolerance.
Definition: Math.h:370
ValueType_ ValueType
Definition: ConjGradient.h:269
ValueType l2Norm() const
Return the L2 norm of this vector.
Definition: ConjGradient.h:213
Vector< ValueType > VectorType
Definition: ConjGradient.h:270
Index32 SizeType
Definition: ConjGradient.h:61
ConstRow getConstRow(SizeType row) const
Return a read-only view onto the given row of this matrix.
Definition: ConjGradient.h:1119
IsFiniteOp(const SparseStencilMatrix &m)
Definition: ConjGradient.h:1063
const boost::disable_if_c< VecTraits< T >::IsVec, T >::type & min(const T &a, const T &b)
Definition: Composite.h:105
ConstValueIter & operator++()
Definition: ConjGradient.h:424
Read-only accessor to a row of this matrix.
Definition: ConjGradient.h:440
T * data
Definition: ConjGradient.h:542
SparseStencilMatrix(SizeType n)
Construct an n x n matrix with at most STENCIL_SIZE nonzero elements per row.
Definition: ConjGradient.h:860
Vector & operator*=(const Scalar &s)
Multiply each element of this vector by s.
Definition: ConjGradient.h:204
SparseStencilMatrix * to
Definition: ConjGradient.h:890
DeterministicDotProductOp(const T *a_, const T *b_, const SizeType binCount_, const SizeType arraySize_, T *reducetmp_)
Definition: ConjGradient.h:662
T & operator[](SizeType i)
Return the value of this vector's ith element.
Definition: ConjGradient.h:231
bool isFinite() const
Return true if all values along the diagonal are finite.
Definition: ConjGradient.h:1351
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:697
ConstRow(const ValueType *valueHead, const SizeType *columnHead, const SizeType &rowSize)
Definition: ConjGradient.h:1237
Sparse, square matrix representing a 3D stencil operator of size STENCIL_SIZE.
Definition: ConjGradient.h:67
void operator()(const SizeRange &range) const
Definition: ConjGradient.h:523
void operator()(const SizeRange &range) const
Definition: ConjGradient.h:538
void operator()(const SizeRange &range) const
Definition: ConjGradient.h:553
boost::shared_ptr< IncompleteCholeskyPreconditioner > Ptr
Definition: ConjGradient.h:1396
const SizeType binCount
Definition: ConjGradient.h:690
const T & at(SizeType i) const
Return the value of this vector's ith element.
Definition: ConjGradient.h:230
SizeType setValue(SizeType column, const ValueType &value)
Set the value of the entry in the specified column.
Definition: ConjGradient.h:1265
#define OPENVDB_VERSION_NAME
Definition: version.h:43
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:1652
const SizeType arraySize
Definition: ConjGradient.h:691
void setValue(SizeType row, SizeType col, const ValueType &)
Set the value at the given coordinates.
Definition: ConjGradient.h:915
void axpy(const T &a, const T *xVec, const T *yVec, T *resultVec, SizeType size)
Compute ax + y.
Definition: ConjGradient.h:1604
Dummy NOOP interrupter class defining interface.
Definition: NullInterrupter.h:52
boost::shared_ptr< JacobiPreconditioner > Ptr
Definition: ConjGradient.h:1330
void operator()(const SizeRange &range) const
Definition: ConjGradient.h:666
Vector< double > VectorD
Definition: ConjGradient.h:256
Lightweight, variable-length vector.
Definition: ConjGradient.h:65
const T * from
Definition: ConjGradient.h:527
#define OPENVDB_LOG_WARN(message)
Log a warning message of the form 'someVar << "some text" << ...'.
Definition: logging.h:39
T * data()
Return a pointer to this vector's elements.
Definition: ConjGradient.h:237
Information about the state of a conjugate gradient solution.
Definition: ConjGradient.h:74
double relativeError
Definition: ConjGradient.h:77
virtual ~Preconditioner()
Definition: ConjGradient.h:500
bool isZero(const Type &x)
Return true if x is exactly equal to zero.
Definition: Math.h:324
Vector(SizeType n)
Construct a vector of n elements, with uninitialized values.
Definition: ConjGradient.h:175
RowEditor getRowEditor(SizeType row)
Return a read/write view onto the given row of this matrix.
Definition: ConjGradient.h:1109
const TriangularMatrix & upperMatrix() const
Definition: ConjGradient.h:1547
bool isFinite(const Type &x)
Return true if x is finite.
Definition: Math.h:363
Definition: Exceptions.h:86
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:1666
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:1050
virtual ~IncompleteCholeskyPreconditioner()
Definition: ConjGradient.h:1494
Definition: Exceptions.h:39
std::string str() const
Return a string representation of this vector.
Definition: ConjGradient.h:837
bool isApproxZero(const Type &x)
Return true if x is equal to zero to within the default floating-point comparison tolerance...
Definition: Math.h:336
const SparseStencilMatrix * mat
Definition: ConjGradient.h:1080
SizeType numRows() const
Return the number of rows in this matrix.
Definition: ConjGradient.h:287
const Type & Max(const Type &a, const Type &b)
Return the maximum of two values.
Definition: Math.h:561
RowEditor(ValueType *valueHead, SizeType *columnHead, SizeType &rowSize, SizeType colSize)
Definition: ConjGradient.h:1247
RowEditor & operator*=(const Scalar &s)
Scale all of the entries in this row.
Definition: ConjGradient.h:464
void computeResidual(const MatrixOperator &A, const Vector< T > &x, const Vector< T > &b, Vector< T > &r)
Compute r = b − Ax.
Definition: ConjGradient.h:1635
Iterator over the stored values in a row of this matrix.
Definition: ConjGradient.h:412
std::string str() const
Return a string representation of this matrix.
Definition: ConjGradient.h:1097
bool operator()(const SizeRange &range, bool finite) const
Definition: ConjGradient.h:1065
T ValueType
Definition: ConjGradient.h:496
void fill(const ValueType &value)
Set all elements of this vector to value.
Definition: ConjGradient.h:629
T * to
Definition: ConjGradient.h:528
MatrixType::ValueType ValueType
Definition: ConjGradient.h:1324
void resize(SizeType n)
Reset this vector to have n elements, with uninitialized values.
Definition: ConjGradient.h:617
void operator()(const SizeRange &range) const
Definition: ConjGradient.h:878
const ValueType & operator()(SizeType row, SizeType col) const
Return the value at the given coordinates.
Definition: ConjGradient.h:934
bool success
Definition: ConjGradient.h:75
ValueType infNorm() const
Return the infinity norm of this vector.
Definition: ConjGradient.h:758
const boost::disable_if_c< VecTraits< T >::IsVec, T >::type & max(const T &a, const T &b)
Definition: Composite.h:109
boost::shared_ptr< Vector > Ptr
Definition: ConjGradient.h:170
Tolerance for floating-point comparison.
Definition: Math.h:125
TriangularMatrix::ConstRow TriangleConstRow
Definition: ConjGradient.h:1398
boost::shared_ptr< Preconditioner > Ptr
Definition: ConjGradient.h:497
virtual void apply(const Vector< ValueType > &r, Vector< ValueType > &z)
Apply this preconditioner to a residue vector: z = M−1r
Definition: ConjGradient.h:1340
int iterations
Definition: ConjGradient.h:76
Definition: Exceptions.h:78
const T * data
Definition: ConjGradient.h:784
virtual bool isValid() const
Definition: ConjGradient.h:1496
bool empty() const
Return true if this vector has no elements.
Definition: ConjGradient.h:189
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:826
Vector< ValueType > VectorType
Definition: ConjGradient.h:1329
const T * y
Definition: ConjGradient.h:563
SizeType size() const
Return the number of rows in this matrix.
Definition: ConjGradient.h:288
IsFiniteOp(const T *data_)
Definition: ConjGradient.h:770
#define OPENVDB_USE_VERSION_NAMESPACE
Definition: version.h:71
float Sqrt(float x)
Return the square root of a floating-point value.
Definition: Math.h:709
double absoluteError
Definition: ConjGradient.h:78
tbb::blocked_range< SizeType > SizeRange
Definition: ConjGradient.h:63
T ValueType
Definition: ConjGradient.h:169
void swap(Vector &other)
Swap internal storage with another vector, which need not be the same size.
Definition: ConjGradient.h:196
MatType scale(const Vec3< typename MatType::value_type > &s)
Return a matrix that scales by s.
Definition: Mat.h:594
void scale(const Scalar &s)
Multiply each element of this vector by s.
Definition: ConjGradient.h:653
InfNormOp(const T *data_)
Definition: ConjGradient.h:740
~Vector()
Definition: ConjGradient.h:179
SizeType size() const
Return the number of elements in this vector.
Definition: ConjGradient.h:187
SparseStencilMatrix< ValueType, 4 > TriangularMatrix
Definition: ConjGradient.h:1397
const Type & Min(const Type &a, const Type &b)
Return the minimum of two values.
Definition: Math.h:622
const T * data() const
Return a pointer to this vector's elements.
Definition: ConjGradient.h:238
T operator()(const SizeRange &range, T maxValue) const
Definition: ConjGradient.h:742
const T val
Definition: ConjGradient.h:543