CppNoddy  0.92
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Public Attributes | List of all members
CppNoddy::BandedMatrix< _Type > Class Template Reference

A matrix class that constructs a BANDED matrix. More...

#include <BandedMatrix.h>

Inheritance diagram for CppNoddy::BandedMatrix< _Type >:
CppNoddy::Sequential_Matrix_base< _Type >

Public Types

typedef DenseVector< _Type >::elt_iter elt_iter
 

Public Member Functions

 BandedMatrix ()
 Empty constructor .. should we stop this? More...
 
void blank ()
 
 BandedMatrix (const std::size_t &rows, const std::size_t &offdiag, const _Type &fill)
 Noddy Banded Matrix constructor. More...
 
 BandedMatrix (const BandedMatrix &source)
 Copy constructor. More...
 
BandedMatrixoperator= (const BandedMatrix &source)
 Assignment operator. More...
 
const _Type & operator() (const std::size_t &row, const std::size_t &col) const
 Access operator. More...
 
_Type & operator() (const std::size_t &row, const std::size_t &col)
 Access operator. More...
 
const _Type & get (const std::size_t &row, const std::size_t &col) const
 Access operator. More...
 
_Type & set (const std::size_t &row, const std::size_t &col)
 Access operator. More...
 
std::size_t nrows () const
 Get the number of rows. More...
 
std::size_t ncols () const
 Get the number of columns. More...
 
std::size_t nelts () const
 Get the number of elements. More...
 
void scale (const _Type &mult)
 Scale all entries in the matrix by a scalar. More...
 
void transpose ()
 Transpose the matrix. More...
 
double one_norm () const
 
double two_norm () const
 
double inf_norm () const
 
double frob_norm () const
 
DenseVector< _Type > multiply (const DenseVector< _Type > &X) const
 Right multiply the matrix by a DENSE vector. More...
 
void dump () const
 Output the matrix contents to std::cout. More...
 
void assign (_Type elt)
 Assign a value the matrix so that it has the same geometry, but zero entries in all locations including those reserved for pivotting. More...
 
std::size_t noffdiag () const
 Get the number of off-diagonal elements where the total INPUT band width is 2*noffdiag+1 since the band structure is symmetric. More...
 
void row_swap (const std::size_t &row1, const std::size_t &row2)
 Exchange rows in the matrix. More...
 
double * base ()
 Allow direct access to the vector m_storage. More...
 
elt_iter get_elt_iter (std::size_t row, std::size_t col)
 
double * base ()
 
double * base ()
 
- Public Member Functions inherited from CppNoddy::Sequential_Matrix_base< _Type >
 Sequential_Matrix_base ()
 An empty constructor. More...
 
virtual ~Sequential_Matrix_base ()
 
virtual const _Type & operator() (const std::size_t &row, const std::size_t &col) const =0
 
virtual _Type & operator() (const std::size_t &row, const std::size_t &col)=0
 
virtual const _Type & get (const std::size_t &row, const std::size_t &col) const =0
 
virtual _Type & set (const std::size_t &row, const std::size_t &col)=0
 
virtual std::size_t nrows () const =0
 
virtual std::size_t ncols () const =0
 
virtual std::size_t nelts () const =0
 
virtual void scale (const _Type &mult)=0
 
virtual void dump () const =0
 

Public Attributes

DenseVector< _Type > m_storage
 A contiguous vector. More...
 
std::size_t m_N
 The number of rows/cols in the matrix. More...
 
std::size_t m_L
 Max number of (INPUT) bands above OR below the main diagonal. More...
 

Detailed Description

template<typename _Type>
class CppNoddy::BandedMatrix< _Type >

A matrix class that constructs a BANDED matrix.

Definition at line 16 of file BandedMatrix.h.

Member Typedef Documentation

◆ elt_iter

template<typename _Type >
typedef DenseVector<_Type>::elt_iter CppNoddy::BandedMatrix< _Type >::elt_iter

Definition at line 20 of file BandedMatrix.h.

Constructor & Destructor Documentation

◆ BandedMatrix() [1/3]

template<typename _Type >
CppNoddy::BandedMatrix< _Type >::BandedMatrix ( )
inline

Empty constructor .. should we stop this?

Definition at line 23 of file BandedMatrix.h.

24 {}

◆ BandedMatrix() [2/3]

template<typename _Type >
CppNoddy::BandedMatrix< _Type >::BandedMatrix ( const std::size_t &  rows,
const std::size_t &  offdiag,
const _Type &  fill 
)

Noddy Banded Matrix constructor.

Parameters
rowsThe number of rows in the matrix.
offdiagThe maximum number of bands above OR below the diagonal. The total number stored will be 3 * offdiag + 1; i.e. assumes the number of offdiagonal elts is the same both above and below, and we must an extra 'offdiag' because of fill-in during pivotting.
fillThe initial value to be placed in each element of the banded matrix.

Definition at line 16 of file BandedMatrix.cpp.

18 :
19 Sequential_Matrix_base<_Type>(),
20 m_N(rows),
21 m_L(upper_offdiag_bands) {
22 // we'll assume that the upper & lower offdiags are of equal size.
23 // logical m_storage is the main diagonal + upper offdiagonal + lower offdiagonal
24 // = 2 * m_L + 1.
25 // However, allowing for row interchanging during pivotting leads to additional
26 // padding of m_L, to make a total of 3*m_L+1.
27 m_storage = DenseVector<_Type>(m_N * (3 * m_L + 1), 0.0);
28 }
std::size_t m_L
Max number of (INPUT) bands above OR below the main diagonal.
Definition: BandedMatrix.h:138
std::size_t m_N
The number of rows/cols in the matrix.
Definition: BandedMatrix.h:136
DenseVector< _Type > m_storage
A contiguous vector.
Definition: BandedMatrix.h:134

References CppNoddy::BandedMatrix< _Type >::m_L, CppNoddy::BandedMatrix< _Type >::m_N, and CppNoddy::BandedMatrix< _Type >::m_storage.

◆ BandedMatrix() [3/3]

template<typename _Type >
CppNoddy::BandedMatrix< _Type >::BandedMatrix ( const BandedMatrix< _Type > &  source)

Copy constructor.

Parameters
sourceThe source object to be copied

Definition at line 31 of file BandedMatrix.cpp.

31 :
32 Sequential_Matrix_base<_Type>() {
33 *this = source;
34 }
double source(const double &x, const double &y, const double &t)
Definition: IBVPLinear.cpp:26

Member Function Documentation

◆ assign()

template<typename _Type >
void CppNoddy::BandedMatrix< _Type >::assign ( _Type  elt)
inline

Assign a value the matrix so that it has the same geometry, but zero entries in all locations including those reserved for pivotting.

Parameters
eltThe value to be assigned to all entries

Definition at line 106 of file BandedMatrix.h.

106 {
107 m_storage.assign(m_storage.size(), elt);
108 }

References CppNoddy::BandedMatrix< _Type >::m_storage.

Referenced by CppNoddy::PDE_IBVP< _Type >::assemble_matrix_problem().

◆ base() [1/3]

template<typename _Type >
double * CppNoddy::BandedMatrix< _Type >::base ( )

Allow direct access to the vector m_storage.

Dangerous, but used for passing to LAPACK.

◆ base() [2/3]

double * CppNoddy::BandedMatrix< double >::base ( )

Definition at line 160 of file BandedMatrix.cpp.

160 {
161 return &(m_storage[0]);
162 }

◆ base() [3/3]

double * CppNoddy::BandedMatrix< std::complex< double > >::base ( )

Definition at line 165 of file BandedMatrix.cpp.

165 {
166 //return &( m_storage[0].real() );
167 return &reinterpret_cast<double(&)[2]>(m_storage[0])[0];
168 }

◆ blank()

template<typename _Type >
void CppNoddy::BandedMatrix< _Type >::blank ( )
inline

Definition at line 26 of file BandedMatrix.h.

26 {
27 m_storage = DenseVector<_Type>(m_N * (3 * m_L + 1), 0.0);
28 }

References CppNoddy::BandedMatrix< _Type >::m_L, CppNoddy::BandedMatrix< _Type >::m_N, and CppNoddy::BandedMatrix< _Type >::m_storage.

◆ dump()

template<typename _Type >
void CppNoddy::BandedMatrix< _Type >::dump
virtual

Output the matrix contents to std::cout.

Implements CppNoddy::Sequential_Matrix_base< _Type >.

Definition at line 137 of file BandedMatrix.cpp.

137 {
138 std::cout << "BANDED mtx size = " << nrows() << " x " << ncols()
139 << " with INPUT no. of bands = " << 2 * m_L + 1
140 << " and total bands STORED (for pivotting) = " << 3 * m_L + 1 << "\n";
141 std::cout.precision(4);
142 std::cout << "- start matrix \n";
143 for(std::size_t i = 0; i < nrows(); ++i) {
144 std::cout << " row " << i << " = ";
145 for(std::size_t j = 0; j < nrows(); ++j) {
146 if(i == j)
147 std::cout << "*";
148 if(((int) j < (int) i - (int) m_L) || ((int) j > (int) i + (int)(2*m_L))) {
149 std::cout << 0 << ", ";
150 } else {
151 std::cout << operator()(i, j) << ", ";
152 }
153 }
154 std::cout << "\n";
155 }
156 std::cout << "- end matrix \n";
157 }
const _Type & operator()(const std::size_t &row, const std::size_t &col) const
Access operator.
Definition: BandedMatrix.h:147
std::size_t nrows() const
Get the number of rows.
Definition: BandedMatrix.h:217
std::size_t ncols() const
Get the number of columns.
Definition: BandedMatrix.h:222

◆ frob_norm()

template<typename _Type >
double CppNoddy::BandedMatrix< _Type >::frob_norm
Returns
The sum of the two_norm of all rows

Definition at line 126 of file BandedMatrix.cpp.

126 {
127 double sum(0.0);
128 for(unsigned row = 0; row < m_N; ++row) {
129 // copy the row into a temp vector
130 DenseVector<_Type> temp(3 * m_L + 1, &m_storage[ row * (3 * m_L + 1) ]);
131 sum += temp.two_norm();
132 }
133 return sum;
134 }

References CppNoddy::DenseVector< _Type >::two_norm().

◆ get()

template<typename _Type >
const _Type & CppNoddy::BandedMatrix< _Type >::get ( const std::size_t &  row,
const std::size_t &  col 
) const
inlinevirtual

Access operator.

Implements CppNoddy::Sequential_Matrix_base< _Type >.

Definition at line 205 of file BandedMatrix.h.

206 {
207 return operator()(row, col);
208 }

◆ get_elt_iter()

template<typename _Type >
elt_iter CppNoddy::BandedMatrix< _Type >::get_elt_iter ( std::size_t  row,
std::size_t  col 
)
inline

Definition at line 127 of file BandedMatrix.h.

127 {
128 return m_storage.begin() + m_L * (3 * col + 2) + row;
129 }

References CppNoddy::BandedMatrix< _Type >::m_L, and CppNoddy::BandedMatrix< _Type >::m_storage.

Referenced by CppNoddy::PDE_IBVP< _Type >::assemble_matrix_problem().

◆ inf_norm()

template<typename _Type >
double CppNoddy::BandedMatrix< _Type >::inf_norm
Returns
The maximum inf_norm of all rows

Definition at line 114 of file BandedMatrix.cpp.

114 {
115 double max(0.0);
116 for(unsigned row = 0; row < m_N; ++row) {
117 // copy the row into a temp vector
118 DenseVector<_Type> temp(3 * m_L + 1, &m_storage[ row * (3 * m_L + 1) ]);
119 max = std::max(max, temp.inf_norm());
120 }
121 return max;
122 }

References CppNoddy::DenseVector< _Type >::inf_norm().

◆ multiply()

template<typename _Type >
DenseVector< _Type > CppNoddy::BandedMatrix< _Type >::multiply ( const DenseVector< _Type > &  X) const

Right multiply the matrix by a DENSE vector.

Parameters
XThe DENSE vector to multiply by
Returns
A DENSE vector of length ncols() produced from the multiplication

Definition at line 58 of file BandedMatrix.cpp.

58 {
59 std::string problem;
60 problem = " The multiply method of the BandedMatrix class has \n";
61 problem += " not been implemented. \n";
62 throw ExceptionRuntime(problem);
63 return m_storage; // dummy return
64 }

◆ ncols()

template<typename _Type >
std::size_t CppNoddy::BandedMatrix< _Type >::ncols
inlinevirtual

Get the number of columns.

Returns
The number of columns in the matrix

Implements CppNoddy::Sequential_Matrix_base< _Type >.

Definition at line 222 of file BandedMatrix.h.

222 {
223 return m_N;
224 }

◆ nelts()

template<typename _Type >
std::size_t CppNoddy::BandedMatrix< _Type >::nelts
virtual

Get the number of elements.

Returns
The number of elements

Implements CppNoddy::Sequential_Matrix_base< _Type >.

Definition at line 47 of file BandedMatrix.cpp.

47 {
48 return m_storage.size();
49 }

◆ noffdiag()

template<typename _Type >
std::size_t CppNoddy::BandedMatrix< _Type >::noffdiag
inline

Get the number of off-diagonal elements where the total INPUT band width is 2*noffdiag+1 since the band structure is symmetric.

Returns
The number of upper or lower offdiagonal bands; NOTE: the number of offdiagonals STORED will typically be three times this to allow for pivotting

Definition at line 227 of file BandedMatrix.h.

227 {
228 return m_L;
229 }

Referenced by CppNoddy::PDE_IBVP< _Type >::assemble_matrix_problem(), and main().

◆ nrows()

template<typename _Type >
std::size_t CppNoddy::BandedMatrix< _Type >::nrows
inlinevirtual

Get the number of rows.

Returns
The number of rows in the matrix

Implements CppNoddy::Sequential_Matrix_base< _Type >.

Definition at line 217 of file BandedMatrix.h.

217 {
218 return m_N;
219 }

◆ one_norm()

template<typename _Type >
double CppNoddy::BandedMatrix< _Type >::one_norm
Returns
The maximum one_norm of all rows

Definition at line 90 of file BandedMatrix.cpp.

90 {
91 double max(0.0);
92 for(unsigned row = 0; row < m_N; ++row) {
93 // copy the row into a temp vector
94 DenseVector<_Type> temp(3 * m_L + 1, &m_storage[ row * (3 * m_L + 1) ]);
95 max = std::max(max, temp.one_norm());
96 }
97 return max;
98 }

References CppNoddy::DenseVector< _Type >::one_norm().

◆ operator()() [1/2]

template<typename _Type >
_Type & CppNoddy::BandedMatrix< _Type >::operator() ( const std::size_t &  row,
const std::size_t &  col 
)
inlinevirtual

Access operator.

Implements CppNoddy::Sequential_Matrix_base< _Type >.

Definition at line 175 of file BandedMatrix.h.

175 {
176#ifdef PARANOID
177 // if outside the m_Nxm_N matrix
178 if((row >= m_N) || (row < 0) || (col >= m_N) || (col < 0)) {
179 std::string problem;
180 problem = " The operator() of BandedMatrix has been called \n";
181 problem += " with a (row, col) index that is outside \n";
182 problem += " the square m_Nxm_N matrix.\n";
183 throw ExceptionRange(problem, m_N, m_N, row, col);
184 }
185 // check if the subscripts are out of the band
186 if(!((col + m_L >= row) && (col <= 2*m_L + row))) {
187 std::string problem;
188 problem = " The operator() of BandedMatrix has been called \n";
189 problem += " with a (row, col) index that is outside \n";
190 problem += " the band structure.\n";
191 std::cout << " m_L = " << m_L << "\n";
192 std::cout << " row = " << row << "\n";
193 std::cout << " col = " << col << "\n";
194 throw ExceptionRange(problem, m_N, m_L, row, col);
195 }
196#endif
197 // MOSTLY WE PUSH BANDED MATRICES TO LAPACK - so we'll keep column major
198 // internal storage ... a la Fortran
199 // return m_storage[ col * ( 3 * m_L + 1 ) + ( row - col ) + 2 * m_L ];
200 // or equiv.
201 return m_storage[ m_L * (3 * col + 2) + row ];
202 }

◆ operator()() [2/2]

template<typename _Type >
const _Type & CppNoddy::BandedMatrix< _Type >::operator() ( const std::size_t &  row,
const std::size_t &  col 
) const
inlinevirtual

Access operator.

Implements CppNoddy::Sequential_Matrix_base< _Type >.

Definition at line 147 of file BandedMatrix.h.

147 {
148#ifdef PARANOID
149 // if outside the m_Nxm_N matrix
150 if((row >= m_N) || (row < 0) || (col >= m_N) || (col < 0)) {
151 std::string problem;
152 problem = " The const operator() of BandedMatrix has been called \n";
153 problem += " with a (row, col) index that is outside \n";
154 problem += " the square m_Nxm_N matrix.\n";
155 throw ExceptionGeom(problem, m_N, m_N, row, col);
156 }
157 // check if the subscripts are out of the band
158 if(!((col + m_L >= row) && (col <= 2*m_L + row))) {
159 std::string problem;
160 problem = " The const operator() of BandedMatrix has been called \n";
161 problem += " with a (row, col) index that is outside \n";
162 problem += " the band structure. Bandwidth and offset from\n";
163 problem += " the diagonal information follows.\n";
164 throw ExceptionGeom(problem, m_L, col - row);
165 }
166#endif
167 // MOSTLY WE PUSH BANDED MATRICES TO LAPACK - so we'll keep column major
168 // internal storage ... a la Fortran
169 // return m_storage[ col * ( 3 * m_L + 1 ) + ( row - col ) + 2 * m_L ];
170 // or equiv.
171 return m_storage[ m_L * (3 * col + 2) + row ];
172 }

◆ operator=()

template<typename _Type >
BandedMatrix< _Type > & CppNoddy::BandedMatrix< _Type >::operator= ( const BandedMatrix< _Type > &  source)
inline

Assignment operator.

Parameters
sourceThe source object for the assignment
Returns
The newly assigned object

Definition at line 37 of file BandedMatrix.cpp.

37 {
38 if(this == &source)
39 return * this;
40 m_storage = source.m_storage;
41 m_N = source.m_N;
42 m_L = source.m_L;
43 return *this;
44 }

◆ row_swap()

template<typename _Type >
void CppNoddy::BandedMatrix< _Type >::row_swap ( const std::size_t &  row1,
const std::size_t &  row2 
)

Exchange rows in the matrix.

Parameters
row1First row to be swapped
row2Second row to be swapped
Todo:
MAKE ME PRIVATE OR ELSE!

Definition at line 79 of file BandedMatrix.cpp.

79 {
80 // WE MUST HAVE row2 > row1 & row1 must have no elements left of the diagonal!
81 // only of any use to native Gaussian elimination solver
82 /// \todo MAKE ME PRIVATE OR ELSE!
83 for(std::size_t col = row1; col <= std::min(row1 + 2 * m_L, m_N - 1); ++col) {
84 std::swap(operator()(row1, col), operator()(row2, col));
85 }
86 }

◆ scale()

template<typename _Type >
void CppNoddy::BandedMatrix< _Type >::scale ( const _Type &  mult)
virtual

Scale all entries in the matrix by a scalar.

Parameters
multThe scalar multiplier

Implements CppNoddy::Sequential_Matrix_base< _Type >.

Definition at line 52 of file BandedMatrix.cpp.

52 {
53 m_storage *= mult;
54 }

◆ set()

template<typename _Type >
_Type & CppNoddy::BandedMatrix< _Type >::set ( const std::size_t &  row,
const std::size_t &  col 
)
inlinevirtual

Access operator.

Implements CppNoddy::Sequential_Matrix_base< _Type >.

Definition at line 211 of file BandedMatrix.h.

212 {
213 return operator()(row, col);
214 }

◆ transpose()

template<typename _Type >
void CppNoddy::BandedMatrix< _Type >::transpose

Transpose the matrix.

Definition at line 68 of file BandedMatrix.cpp.

68 {
69 for(std::size_t row = 0; row < m_N; ++row) {
70 for(std::size_t col = std::max(0, int(row) - int(m_L)); col < row; ++col) {
71 // swap elements
72 std::swap(operator()(row, col), operator()(col, row));
73 }
74 }
75 }

◆ two_norm()

template<typename _Type >
double CppNoddy::BandedMatrix< _Type >::two_norm
Returns
The maximum two_norm of all rows

Definition at line 102 of file BandedMatrix.cpp.

102 {
103 double max(0.0);
104 for(unsigned row = 0; row < m_N; ++row) {
105 // copy the row into a temp vector
106 DenseVector<_Type> temp(3 * m_L + 1, &m_storage[ row * (3 * m_L + 1) ]);
107 max = std::max(max, temp.two_norm());
108 }
109 return max;
110 }

References CppNoddy::DenseVector< _Type >::two_norm().

Member Data Documentation

◆ m_L

template<typename _Type >
std::size_t CppNoddy::BandedMatrix< _Type >::m_L

Max number of (INPUT) bands above OR below the main diagonal.

Definition at line 138 of file BandedMatrix.h.

Referenced by CppNoddy::BandedMatrix< _Type >::BandedMatrix(), CppNoddy::BandedMatrix< _Type >::blank(), and CppNoddy::BandedMatrix< _Type >::get_elt_iter().

◆ m_N

template<typename _Type >
std::size_t CppNoddy::BandedMatrix< _Type >::m_N

The number of rows/cols in the matrix.

Definition at line 136 of file BandedMatrix.h.

Referenced by CppNoddy::BandedMatrix< _Type >::BandedMatrix(), and CppNoddy::BandedMatrix< _Type >::blank().

◆ m_storage

template<typename _Type >
DenseVector<_Type> CppNoddy::BandedMatrix< _Type >::m_storage

The documentation for this class was generated from the following files:

© 2012

R.E. Hewitt