/usr/include/dune/istl/matrix.hh is in libdune-istl-dev 2.5.1-1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 | // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_ISTL_MATRIX_HH
#define DUNE_ISTL_MATRIX_HH
/** \file
\brief A dynamic dense block matrix class
*/
#include <memory>
#include <cmath>
#include <dune/common/ftraits.hh>
#include <dune/istl/istlexception.hh>
#include <dune/istl/bvector.hh>
namespace Dune {
namespace MatrixImp
{
/**
\brief A Vector of blocks with different blocksizes.
This class started as a copy of VariableBlockVector, which used to be used for the internal memory managerment
of the 'Matrix' class. However, that mechanism stopped working when I started using the RandomAccessIteratorFacade
in VariableBlockVector (308dd85483108f8baaa4051251e2c75e2a9aed32, to make VariableBlockVector pass a number of
tightened interface compliance tests), and I couldn't quite figure out how to fix that. However, using
VariableBlockVector in Matrix internally was a hack anyway, so I simply took the working version of VariableBlockVector
and copied it here under the new name of DenseMatrixBase. This is still hacky, but one step closer to an
elegant solution.
*/
template<class B, class A=std::allocator<B> >
class DenseMatrixBase : public block_vector_unmanaged<B,A>
// this derivation gives us all the blas level 1 and norms
// on the large array. However, access operators have to be
// overwritten.
{
public:
//===== type definitions and constants
//! export the type representing the field
typedef typename B::field_type field_type;
//! export the allocator type
typedef A allocator_type;
//! The size type for the index access
typedef typename A::size_type size_type;
/** \brief Type of the elements of the outer vector, i.e., dynamic vectors of B
*
* Note that this is *not* the type referred to by the iterators and random access operators,
* which return proxy objects.
*/
typedef BlockVector<B,A> value_type;
/** \brief Same as value_type, here for historical reasons
*/
typedef BlockVector<B,A> block_type;
// just a shorthand
typedef BlockVectorWindow<B,A> window_type;
typedef window_type reference;
typedef const window_type const_reference;
//===== constructors and such
/** constructor without arguments makes empty vector,
object cannot be used yet
*/
DenseMatrixBase () : block_vector_unmanaged<B,A>()
{
// nothing is known ...
rows_ = 0;
columns_ = 0;
}
/** make vector with given number of blocks each having a constant size,
object is fully usable then.
\param _nblocks Number of blocks
\param m Number of elements in each block
*/
DenseMatrixBase (size_type rows, size_type columns) : block_vector_unmanaged<B,A>()
{
// and we can allocate the big array in the base class
this->n = rows*columns;
columns_ = columns;
if (this->n>0)
{
this->p = allocator_.allocate(this->n);
new (this->p)B[this->n];
}
else
{
this->n = 0;
this->p = 0;
}
// we can allocate the windows now
rows_ = rows;
}
//! copy constructor, has copy semantics
DenseMatrixBase (const DenseMatrixBase& a)
{
// allocate the big array in the base class
this->n = a.n;
columns_ = a.columns_;
if (this->n>0)
{
// allocate and construct objects
this->p = allocator_.allocate(this->n);
new (this->p)B[this->n];
// copy data
for (size_type i=0; i<this->n; i++)
this->p[i]=a.p[i];
}
else
{
this->n = 0;
this->p = nullptr;
}
// we can allocate the windows now
rows_ = a.rows_;
}
//! free dynamic memory
~DenseMatrixBase ()
{
if (this->n>0) {
size_type i=this->n;
while (i)
this->p[--i].~B();
allocator_.deallocate(this->p,this->n);
}
}
//! same effect as constructor with same argument
void resize (size_type rows, size_type columns)
{
// deconstruct objects and deallocate memory if necessary
if (this->n>0) {
size_type i=this->n;
while (i)
this->p[--i].~B();
allocator_.deallocate(this->p,this->n);
}
// and we can allocate the big array in the base class
this->n = rows*columns;
if (this->n>0)
{
this->p = allocator_.allocate(this->n);
new (this->p)B[this->n];
}
else
{
this->n = 0;
this->p = nullptr;
}
// we can allocate the windows now
rows_ = rows;
columns_ = columns;
}
//! assignment
DenseMatrixBase& operator= (const DenseMatrixBase& a)
{
if (&a!=this) // check if this and a are different objects
{
columns_ = a.columns_;
// reallocate arrays if necessary
// Note: still the block sizes may vary !
if (this->n!=a.n || rows_!=a.rows_)
{
// deconstruct objects and deallocate memory if necessary
if (this->n>0) {
size_type i=this->n;
while (i)
this->p[--i].~B();
allocator_.deallocate(this->p,this->n);
}
// allocate the big array in the base class
this->n = a.n;
if (this->n>0)
{
// allocate and construct objects
this->p = allocator_.allocate(this->n);
new (this->p)B[this->n];
}
else
{
this->n = 0;
this->p = nullptr;
}
// Copy number of rows
rows_ = a.rows_;
}
// and copy the data
for (size_type i=0; i<this->n; i++)
this->p[i]=a.p[i];
}
return *this;
}
//===== assignment from scalar
//! assign from scalar
DenseMatrixBase& operator= (const field_type& k)
{
(static_cast<block_vector_unmanaged<B,A>&>(*this)) = k;
return *this;
}
//===== access to components
// has to be overwritten from base class because it must
// return access to the windows
//! random access to blocks
reference operator[] (size_type i)
{
#ifdef DUNE_ISTL_WITH_CHECKING
if (i>=rows_) DUNE_THROW(ISTLError,"index out of range");
#endif
return window_type(this->p + i*columns_, columns_);
}
//! same for read only access
const_reference operator[] (size_type i) const
{
#ifdef DUNE_ISTL_WITH_CHECKING
if (i<0 || i>=rows_) DUNE_THROW(ISTLError,"index out of range");
#endif
return window_type(this->p + i*columns_, columns_);
}
// forward declaration
class ConstIterator;
//! Iterator class for sequential access
class Iterator
{
public:
//! constructor, no arguments
Iterator ()
: window_(nullptr,0)
{
i = 0;
}
Iterator (Iterator& other) = default;
Iterator (Iterator&& other) = default;
//! constructor
Iterator (B* data, size_type columns, size_type _i)
: i(_i),
window_(data + _i*columns, columns)
{}
/** \brief Move assignment */
Iterator& operator=(Iterator&& other)
{
i = other.i;
// Do NOT use window_.operator=, because that copies the window content, not just the window!
window_.set(other.window_.getsize(),other.window_.getptr());
return *this;
}
/** \brief Copy assignment */
Iterator& operator=(Iterator& other)
{
i = other.i;
// Do NOT use window_.operator=, because that copies the window content, not just the window!
window_.set(other.window_.getsize(),other.window_.getptr());
return *this;
}
//! prefix increment
Iterator& operator++()
{
++i;
window_.setptr(window_.getptr()+window_.getsize());
return *this;
}
//! prefix decrement
Iterator& operator--()
{
--i;
window_.setptr(window_.getptr()-window_.getsize());
return *this;
}
//! equality
bool operator== (const Iterator& it) const
{
return window_.getptr() == it.window_.getptr();
}
//! inequality
bool operator!= (const Iterator& it) const
{
return window_.getptr() != it.window_.getptr();
}
//! equality
bool operator== (const ConstIterator& it) const
{
return window_.getptr() == it.window_.getptr();
}
//! inequality
bool operator!= (const ConstIterator& it) const
{
return window_.getptr() != it.window_.getptr();
}
//! dereferencing
window_type& operator* () const
{
return window_;
}
//! arrow
window_type* operator-> () const
{
return &window_;
}
// return index corresponding to pointer
size_type index () const
{
return i;
}
friend class ConstIterator;
private:
size_type i;
mutable window_type window_;
};
//! begin Iterator
Iterator begin ()
{
return Iterator(this->p, columns_, 0);
}
//! end Iterator
Iterator end ()
{
return Iterator(this->p, columns_, rows_);
}
//! @returns an iterator that is positioned before
//! the end iterator of the vector, i.e. at the last entry.
Iterator beforeEnd ()
{
return Iterator(this->p, columns_, rows_-1);
}
//! @returns an iterator that is positioned before
//! the first entry of the vector.
Iterator beforeBegin () const
{
return Iterator(this->p, columns_, -1);
}
//! random access returning iterator (end if not contained)
Iterator find (size_type i)
{
return Iterator(this->p, columns_, std::min(i,rows_));
}
//! random access returning iterator (end if not contained)
ConstIterator find (size_type i) const
{
return ConstIterator(this->p, columns_, std::min(i,rows_));
}
//! ConstIterator class for sequential access
class ConstIterator
{
public:
//! constructor
ConstIterator ()
: window_(nullptr,0)
{
i = 0;
}
//! constructor from pointer
ConstIterator (const B* data, size_type columns, size_type _i)
: i(_i),
window_(const_cast<B*>(data + _i * columns), columns)
{}
//! constructor from non_const iterator
ConstIterator (const Iterator& it)
: i(it.i), window_(it.window_.getptr(),it.window_.getsize())
{}
ConstIterator& operator=(Iterator&& other)
{
i = other.i;
// Do NOT use window_.operator=, because that copies the window content, not just the window!
window_.set(other.window_.getsize(),other.window_.getptr());
return *this;
}
ConstIterator& operator=(Iterator& other)
{
i = other.i;
// Do NOT use window_.operator=, because that copies the window content, not just the window!
window_.set(other.window_.getsize(),other.window_.getptr());
return *this;
}
//! prefix increment
ConstIterator& operator++()
{
++i;
window_.setptr(window_.getptr()+window_.getsize());
return *this;
}
//! prefix decrement
ConstIterator& operator--()
{
--i;
window_.setptr(window_.getptr()-window_.getsize());
return *this;
}
//! equality
bool operator== (const ConstIterator& it) const
{
return window_.getptr() == it.window_.getptr();
}
//! inequality
bool operator!= (const ConstIterator& it) const
{
return window_.getptr() != it.window_.getptr();
}
//! equality
bool operator== (const Iterator& it) const
{
return window_.getptr() == it.window_.getptr();
}
//! inequality
bool operator!= (const Iterator& it) const
{
return window_.getptr() != it.window_.getptr();
}
//! dereferencing
const window_type& operator* () const
{
return window_;
}
//! arrow
const window_type* operator-> () const
{
return &window_;
}
// return index corresponding to pointer
size_type index () const
{
return i;
}
friend class Iterator;
private:
size_type i;
mutable window_type window_;
};
/** \brief Export the iterator type using std naming rules */
using iterator = Iterator;
/** \brief Export the const iterator type using std naming rules */
using const_iterator = ConstIterator;
//! begin ConstIterator
ConstIterator begin () const
{
return ConstIterator(this->p, columns_, 0);
}
//! end ConstIterator
ConstIterator end () const
{
return ConstIterator(this->p, columns_, rows_);
}
//! @returns an iterator that is positioned before
//! the end iterator of the vector. i.e. at the last element.
ConstIterator beforeEnd() const
{
return ConstIterator(this->p, columns_, rows_-1);
}
//! end ConstIterator
ConstIterator rend () const
{
return ConstIterator(this->p, columns_, -1);
}
//===== sizes
//! number of blocks in the vector (are of variable size here)
size_type N () const
{
return rows_;
}
private:
size_type rows_; // number of matrix rows
size_type columns_; // number of matrix columns
A allocator_;
};
} // namespace MatrixImp
/** \addtogroup ISTL_SPMV
\{
*/
/** \brief A generic dynamic dense matrix
*/
template<class T, class A=std::allocator<T> >
class Matrix
{
public:
/** \brief Export the type representing the underlying field */
typedef typename T::field_type field_type;
/** \brief Export the type representing the components */
typedef T block_type;
/** \brief Export the allocator */
typedef A allocator_type;
/** \brief The type implementing a matrix row */
typedef typename MatrixImp::DenseMatrixBase<T,A>::window_type row_type;
/** \brief Type for indices and sizes */
typedef typename A::size_type size_type;
/** \brief Iterator over the matrix rows */
typedef typename MatrixImp::DenseMatrixBase<T,A>::Iterator RowIterator;
/** \brief Iterator for the entries of each row */
typedef typename row_type::iterator ColIterator;
/** \brief Const iterator over the matrix rows */
typedef typename MatrixImp::DenseMatrixBase<T,A>::ConstIterator ConstRowIterator;
/** \brief Const iterator for the entries of each row */
typedef typename row_type::const_iterator ConstColIterator;
enum {
//! The number of nesting levels the matrix contains.
blocklevel = T::blocklevel+1
};
/** \brief Create empty matrix */
Matrix() : data_(0,0), cols_(0)
{}
/** \brief Create uninitialized matrix of size rows x cols
*/
Matrix(size_type rows, size_type cols) : data_(rows,cols), cols_(cols)
{}
/** \brief Change the matrix size
*
* The way the data is handled is unpredictable.
*/
void setSize(size_type rows, size_type cols) {
data_.resize(rows,cols);
cols_ = cols;
}
/** \brief Get iterator to first row */
RowIterator begin()
{
return data_.begin();
}
/** \brief Get iterator to one beyond last row */
RowIterator end()
{
return data_.end();
}
//! @returns an iterator that is positioned before
//! the end iterator of the rows, i.e. at the last row.
RowIterator beforeEnd ()
{
return data_.beforeEnd();
}
//! @returns an iterator that is positioned before
//! the first row of the matrix.
RowIterator beforeBegin ()
{
return data_.beforeBegin();
}
/** \brief Get const iterator to first row */
ConstRowIterator begin() const
{
return data_.begin();
}
/** \brief Get const iterator to one beyond last row */
ConstRowIterator end() const
{
return data_.end();
}
//! @returns an iterator that is positioned before
//! the end iterator of the rows. i.e. at the last row.
ConstRowIterator beforeEnd() const
{
return data_.beforeEnd();
}
//! @returns an iterator that is positioned before
//! the first row if the matrix.
ConstRowIterator beforeBegin () const
{
return data_.beforeBegin();
}
/** \brief Assignment from scalar */
Matrix& operator= (const field_type& t)
{
data_ = t;
return *this;
}
/** \brief The index operator */
row_type operator[](size_type row) {
#ifdef DUNE_ISTL_WITH_CHECKING
if (row<0)
DUNE_THROW(ISTLError, "Can't access negative rows!");
if (row>=N())
DUNE_THROW(ISTLError, "Row index out of range!");
#endif
return data_[row];
}
/** \brief The const index operator */
const row_type operator[](size_type row) const {
#ifdef DUNE_ISTL_WITH_CHECKING
if (row<0)
DUNE_THROW(ISTLError, "Can't access negative rows!");
if (row>=N())
DUNE_THROW(ISTLError, "Row index out of range!");
#endif
return data_[row];
}
/** \brief Return the number of rows */
size_type N() const {
return data_.N();
}
/** \brief Return the number of columns */
size_type M() const {
return cols_;
}
/** \brief Multiplication with a scalar */
Matrix<T>& operator*=(const field_type& scalar) {
data_ *= scalar;
return (*this);
}
/** \brief Division by a scalar */
Matrix<T>& operator/=(const field_type& scalar) {
data_ /= scalar;
return (*this);
}
/*! \brief Add the entries of another matrix to this one.
*
* \param b The matrix to add to this one. Its size has to
* be the same as the size of this matrix.
*/
Matrix& operator+= (const Matrix& b) {
#ifdef DUNE_ISTL_WITH_CHECKING
if(N()!=b.N() || M() != b.M())
DUNE_THROW(RangeError, "Matrix sizes do not match!");
#endif
data_ += b.data_;
return (*this);
}
/*! \brief Subtract the entries of another matrix from this one.
*
* \param b The matrix to subtract from this one. Its size has to
* be the same as the size of this matrix.
*/
Matrix& operator-= (const Matrix& b) {
#ifdef DUNE_ISTL_WITH_CHECKING
if(N()!=b.N() || M() != b.M())
DUNE_THROW(RangeError, "Matrix sizes do not match!");
#endif
data_ -= b.data_;
return (*this);
}
/** \brief Return the transpose of the matrix */
Matrix transpose() const {
Matrix out(M(), N());
for (size_type i=0; i<N(); i++)
for (size_type j=0; j<M(); j++)
out[j][i] = (*this)[i][j];
return out;
}
/// Generic matrix multiplication.
friend Matrix<T> operator*(const Matrix<T>& m1, const Matrix<T>& m2) {
Matrix<T> out(m1.N(), m2.M());
out = 0;
for (size_type i=0; i<out.N(); i++ ) {
for ( size_type j=0; j<out.M(); j++ )
for (size_type k=0; k<m1.M(); k++)
out[i][j] += m1[i][k]*m2[k][j];
}
return out;
}
/// Generic matrix-vector multiplication.
template <class X, class Y>
friend Y operator*(const Matrix<T>& m, const X& vec) {
#ifdef DUNE_ISTL_WITH_CHECKING
if (m.M()!=vec.size())
DUNE_THROW(ISTLError, "Vector size doesn't match the number of matrix columns!");
#endif
Y out(m.N());
out = 0;
for (size_type i=0; i<out.size(); i++ ) {
for ( size_type j=0; j<vec.size(); j++ )
out[i] += m[i][j]*vec[j];
}
return out;
}
//! y = A x
template <class X, class Y>
void mv(const X& x, Y& y) const
{
#ifdef DUNE_ISTL_WITH_CHECKING
if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
#endif
for (size_type i=0; i<data_.N(); i++) {
y[i]=0;
for (size_type j=0; j<cols_; j++)
(*this)[i][j].umv(x[j], y[i]);
}
}
//! y = A^T x
template<class X, class Y>
void mtv (const X& x, Y& y) const
{
#ifdef DUNE_ISTL_WITH_CHECKING
if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
#endif
for(size_type i=0; i<y.N(); ++i)
y[i]=0;
umtv(x,y);
}
//! y += A x
template <class X, class Y>
void umv(const X& x, Y& y) const
{
#ifdef DUNE_ISTL_WITH_CHECKING
if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
#endif
for (size_type i=0; i<data_.N(); i++) {
for (size_type j=0; j<cols_; j++)
(*this)[i][j].umv(x[j], y[i]);
}
}
//! y -= A x
template<class X, class Y>
void mmv (const X& x, Y& y) const
{
#ifdef DUNE_ISTL_WITH_CHECKING
if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
#endif
ConstRowIterator endi=end();
for (ConstRowIterator i=begin(); i!=endi; ++i)
{
ConstColIterator endj = (*i).end();
for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
(*j).mmv(x[j.index()],y[i.index()]);
}
}
/** \brief \f$ y += \alpha A x \f$ */
template <class X, class Y>
void usmv(const field_type& alpha, const X& x, Y& y) const
{
#ifdef DUNE_ISTL_WITH_CHECKING
if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
#endif
for (size_type i=0; i<data_.N(); i++) {
for (size_type j=0; j<cols_; j++)
(*this)[i][j].usmv(alpha, x[j], y[i]);
}
}
//! y += A^T x
template<class X, class Y>
void umtv (const X& x, Y& y) const
{
#ifdef DUNE_ISTL_WITH_CHECKING
if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
#endif
ConstRowIterator endi=end();
for (ConstRowIterator i=begin(); i!=endi; ++i)
{
ConstColIterator endj = (*i).end();
for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
(*j).umtv(x[i.index()],y[j.index()]);
}
}
//! y -= A^T x
template<class X, class Y>
void mmtv (const X& x, Y& y) const
{
#ifdef DUNE_ISTL_WITH_CHECKING
if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
#endif
ConstRowIterator endi=end();
for (ConstRowIterator i=begin(); i!=endi; ++i)
{
ConstColIterator endj = (*i).end();
for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
(*j).mmtv(x[i.index()],y[j.index()]);
}
}
//! y += alpha A^T x
template<class X, class Y>
void usmtv (const field_type& alpha, const X& x, Y& y) const
{
#ifdef DUNE_ISTL_WITH_CHECKING
if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
#endif
ConstRowIterator endi=end();
for (ConstRowIterator i=begin(); i!=endi; ++i)
{
ConstColIterator endj = (*i).end();
for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
(*j).usmtv(alpha,x[i.index()],y[j.index()]);
}
}
//! y += A^H x
template<class X, class Y>
void umhv (const X& x, Y& y) const
{
#ifdef DUNE_ISTL_WITH_CHECKING
if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
#endif
ConstRowIterator endi=end();
for (ConstRowIterator i=begin(); i!=endi; ++i)
{
ConstColIterator endj = (*i).end();
for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
(*j).umhv(x[i.index()],y[j.index()]);
}
}
//! y -= A^H x
template<class X, class Y>
void mmhv (const X& x, Y& y) const
{
#ifdef DUNE_ISTL_WITH_CHECKING
if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
#endif
ConstRowIterator endi=end();
for (ConstRowIterator i=begin(); i!=endi; ++i)
{
ConstColIterator endj = (*i).end();
for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
(*j).mmhv(x[i.index()],y[j.index()]);
}
}
//! y += alpha A^H x
template<class X, class Y>
void usmhv (const field_type& alpha, const X& x, Y& y) const
{
#ifdef DUNE_ISTL_WITH_CHECKING
if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
#endif
ConstRowIterator endi=end();
for (ConstRowIterator i=begin(); i!=endi; ++i)
{
ConstColIterator endj = (*i).end();
for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
(*j).usmhv(alpha,x[i.index()],y[j.index()]);
}
}
//===== norms
//! frobenius norm: sqrt(sum over squared values of entries)
typename FieldTraits<field_type>::real_type frobenius_norm () const
{
return std::sqrt(frobenius_norm2());
}
//! square of frobenius norm, need for block recursion
typename FieldTraits<field_type>::real_type frobenius_norm2 () const
{
double sum=0;
for (size_type i=0; i<N(); ++i)
for (size_type j=0; j<M(); ++j)
sum += data_[i][j].frobenius_norm2();
return sum;
}
//! infinity norm (row sum norm, how to generalize for blocks?)
template <typename ft = field_type,
typename std::enable_if<!has_nan<ft>::value, int>::type = 0>
typename FieldTraits<ft>::real_type infinity_norm() const {
using real_type = typename FieldTraits<ft>::real_type;
using std::max;
real_type norm = 0;
for (auto const &x : *this) {
real_type sum = 0;
for (auto const &y : x)
sum += y.infinity_norm();
norm = max(sum, norm);
}
return norm;
}
//! simplified infinity norm (uses Manhattan norm for complex values)
template <typename ft = field_type,
typename std::enable_if<!has_nan<ft>::value, int>::type = 0>
typename FieldTraits<ft>::real_type infinity_norm_real() const {
using real_type = typename FieldTraits<ft>::real_type;
using std::max;
real_type norm = 0;
for (auto const &x : *this) {
real_type sum = 0;
for (auto const &y : x)
sum += y.infinity_norm_real();
norm = max(sum, norm);
}
return norm;
}
//! infinity norm (row sum norm, how to generalize for blocks?)
template <typename ft = field_type,
typename std::enable_if<has_nan<ft>::value, int>::type = 0>
typename FieldTraits<ft>::real_type infinity_norm() const {
using real_type = typename FieldTraits<ft>::real_type;
using std::max;
real_type norm = 0;
real_type isNaN = 1;
for (auto const &x : *this) {
real_type sum = 0;
for (auto const &y : x)
sum += y.infinity_norm();
norm = max(sum, norm);
isNaN += sum;
}
isNaN /= isNaN;
return norm * isNaN;
}
//! simplified infinity norm (uses Manhattan norm for complex values)
template <typename ft = field_type,
typename std::enable_if<has_nan<ft>::value, int>::type = 0>
typename FieldTraits<ft>::real_type infinity_norm_real() const {
using real_type = typename FieldTraits<ft>::real_type;
using std::max;
real_type norm = 0;
real_type isNaN = 1;
for (auto const &x : *this) {
real_type sum = 0;
for (auto const &y : x)
sum += y.infinity_norm_real();
norm = max(sum, norm);
isNaN += sum;
}
isNaN /= isNaN;
return norm * isNaN;
}
//===== query
//! return true if (i,j) is in pattern
bool exists (size_type i, size_type j) const
{
#ifdef DUNE_ISTL_WITH_CHECKING
if (i<0 || i>=N()) DUNE_THROW(ISTLError,"row index out of range");
if (j<0 || i>=M()) DUNE_THROW(ISTLError,"column index out of range");
#else
DUNE_UNUSED_PARAMETER(i); DUNE_UNUSED_PARAMETER(j);
#endif
return true;
}
protected:
/** \brief Abuse DenseMatrixBase as an engine for a 2d array ISTL-style
*/
MatrixImp::DenseMatrixBase<T,A> data_;
/** \brief Number of columns of the matrix
In general you can extract the same information from the data_ member. However if you
want to be able to properly handle 0xn matrices then you need a separate member.
*/
size_type cols_;
};
/** \} */
} // end namespace Dune
#endif
|