ʻOhana
Population structure, admixture history, and selection using learning methods.
Public Types | Public Member Functions | Static Public Member Functions
jade::basic_matrix< TValue > Class Template Reference

A template class for a matrix. More...

#include <jade.matrix.hpp>

+ Inheritance diagram for jade::basic_matrix< TValue >:
+ Collaboration diagram for jade::basic_matrix< TValue >:

Public Types

typedef TValue value_type
 The value type. More...
 
typedef basic_blas< value_typeblas_type
 The BLAS type. More...
 
typedef basic_lapack< value_typelapack_type
 The LAPACK type. More...
 
typedef std::initializer_list< std::initializer_list< value_type > > initializer_list_type
 The initializer list type. More...
 

Public Member Functions

 basic_matrix ()
 Initializes a new instance of the class with no size. More...
 
 basic_matrix (const size_t cy, const size_t cx)
 Initializes a new instance of the class based on the specified width and height. More...
 
 basic_matrix (const char *const path)
 Initializes a new instance of the class based on values from the specified file. More...
 
 basic_matrix (const std::string &path)
 Initializes a new instance of the class based on values from the specified file. More...
 
 basic_matrix (std::istream &in)
 Initializes a new instance of the class based on values from the specified input stream. More...
 
 basic_matrix (const initializer_list_type &values)
 Initializes a new instance of the class based on the specified values. More...
 
template<typename TPredicate >
bool all_of (const TPredicate predicate) const
 
template<typename TPredicate >
bool any_of (const TPredicate predicate) const
 
void clamp (const value_type min, const value_type max)
 Clamps all values in the matrix to the specified range. More...
 
void clamp_column (const size_t column, const value_type min, const value_type max)
 Clamps all values in a column to the specified range. More...
 
void clamp_row (const size_t row, const value_type min, const value_type max)
 Clamps all values in a row to the specified range. More...
 
bool contains_inf () const
 
bool contains_nan () const
 
basic_matrix copy_column (const size_t column) const
 
void copy_column (const size_t column, basic_matrix &out) const
 Copies a column into the specified vector. More...
 
void copy_lower_to_upper ()
 Copies the elements in the lower triangle of the square matrix to the corresponding upper triangle elements. More...
 
basic_matrix copy_row (const size_t row) const
 
void copy_row (const size_t row, basic_matrix &out) const
 Copies a row into the specified matrix. More...
 
void copy_upper_to_lower ()
 Copies the elements in the upper triangle of the square matrix to the corresponding lower triangle elements. More...
 
basic_matrix create_transpose () const
 
void create_transpose (basic_matrix &out) const
 Stores the transpose of the matrix. More...
 
bool gesv ()
 Computes the solution to the system of linear equations with a square coefficient matrix A and multiple right-hand sides. Specifically, this method solves for X the system of linear equations A*X = B, where A is an n-by-n matrix, the columns of matrix B are individual right-hand sides, and the columns of X are the corresponding solutions. More...
 
value_type get_column_sum (const size_t column) const
 
const value_typeget_data () const
 
value_typeget_data ()
 
const value_typeget_data (const size_t index) const
 
value_typeget_data (const size_t index)
 
const value_typeget_data (const size_t row, const size_t column) const
 
value_typeget_data (const size_t row, const size_t column)
 
size_t get_height () const
 
size_t get_index (const size_t row, const size_t column) const
 
size_t get_length () const
 
value_type get_max_value () const
 Returns the maximum element value. This method should not be called if the matrix is empty. More...
 
value_type get_min_value () const
 Returns the minimum element value. This method should not be called if the matrix is empty. More...
 
bool get_min_max (value_type &min, value_type &max) const
 Returns the minimum and maximum elements in the matrix. If the matrix is empty, this method returns false. Otherwise, this method stores the results in the specified arguments and returns true. More...
 
bool get_min_max_column (const size_t column, value_type &min, value_type &max) const
 Returns the minimum and maximum elements in a column. If the matrix is empty, this method returns false. Otherwise, this method stores the results in the specified arguments and returns true. More...
 
value_type get_row_sum (const size_t row) const
 
std::string get_size_str () const
 
value_type get_sum () const
 
size_t get_width () const
 
value_type get_value (const size_t index) const
 
value_type get_value (const size_t row, const size_t column) const
 
bool invert (value_type &log_det)
 Computes and stores the inverse of this matrix using the Cholesky square root method. Only the lower-left triangle is used to perform the calculation. If this method fails, indicating this matrix is not positive semidefinite, then the contents of this instance and the parameter are left undefined, and the method returns false. Otherwise, the lower triangle is copied to the upper triangle, the parameter is assigned the log of the determinant of the original matrix values, and the method returns true. More...
 
bool invert ()
 Computes and stores the inverse of this matrix using the Cholesky square root method. Only the lower-left triangle is used to perform the calculation. If this method fails, indicating this matrix is not positive semidefinite, then the contents of this instance is left undefined, and the method returns false. Otherwise, the lower triangle is copied to the upper triangle, and the method returns true. More...
 
bool is_column_vector () const
 
bool is_empty () const
 
bool is_row_vector () const
 
bool is_length (const basic_matrix &other) const
 
bool is_length (const size_t length) const
 
bool is_size (const basic_matrix &other) const
 
bool is_size (const size_t height, const size_t width) const
 
bool is_square () const
 
bool is_vector () const
 
void multiply_column (const size_t column, const value_type value)
 Multiplies a column by a specified value. More...
 
value_type multiply_column (const size_t column, const basic_matrix &vector)
 
void multiply_row (const size_t row, const value_type value)
 Multiplies a row by a specified value. More...
 
value_type multiply_row (const size_t row, const basic_matrix &vector) const
 
bool potrf_lower ()
 Forms the Cholesky factorization of a symmetric positive-definite matrix: [A = L * L^T] where L is the lower triangular portion of the matrix. More...
 
bool potri_lower ()
 Computes the inverse inv(A) of a symmetric positive definite matrix A. Before calling this routine, call dpotrf to factorize A. More...
 
void read (const std::string &path)
 Reads the matrix values from the specified file. More...
 
void read (const char *const path)
 Reads the matrix values from the specified file. More...
 
void read (std::istream &in)
 Reads the matrix values from the specified stream. More...
 
void resize (const size_t height, const size_t width)
 Resizes the matrix to the specified dimensions. The values of the matrix are not reset. More...
 
void set_value (const size_t index, const value_type &value)
 Sets a single value of the vector. More...
 
void set_value (const size_t row, const size_t column, const value_type &value)
 Sets a single value of the matrix. More...
 
void set_values (const value_type value)
 Sets all values of the matrix to the specified value. More...
 
std::string str () const
 
void swap (basic_matrix &other)
 Swaps this matrix and another matrix. More...
 
void transpose ()
 Transposes this matrix. More...
 
void write (const std::string &path) const
 Writes this matrix to the specified file. More...
 
void write (const char *const path) const
 Writes this matrix to the specified file. More...
 
void write (std::ostream &out) const
 Writes this matrix to the specified output stream. More...
 
const value_typeoperator[] (const size_t rhs) const
 
value_typeoperator[] (const size_t rhs)
 
const value_typeoperator() (const size_t rhs_1, const size_t rhs_2) const
 
value_typeoperator() (const size_t rhs_1, const size_t rhs_2)
 
basic_matrixoperator+= (const basic_matrix &rhs)
 Adds the specified matrix to this instance. More...
 
basic_matrixoperator+= (const value_type rhs)
 Adds the specified scalar to this instance. More...
 
basic_matrixoperator-= (const basic_matrix &rhs)
 Subtracts the specified matrix from this instance. More...
 
basic_matrixoperator-= (const value_type rhs)
 Subtracts the specified scalar from this instance. More...
 
basic_matrixoperator*= (const basic_matrix &rhs)
 Multiplies this matrix and the specified matrix and assigns the product to this instance. More...
 
basic_matrixoperator*= (const value_type rhs)
 Multiplies this matrix and the specified scalar and assigns the product to this instance. More...
 
basic_matrixoperator/= (const value_type rhs)
 Divides this matrix by the specified scalar and assigns the quotient to this instance. More...
 
basic_matrix operator+ (const basic_matrix &rhs) const
 Adds a matrix. More...
 
basic_matrix operator+ (const value_type rhs) const
 Adds a value. More...
 
basic_matrix operator- (const basic_matrix &rhs) const
 Subtracts a matrix. More...
 
basic_matrix operator- (const value_type rhs) const
 Subtracts a value. More...
 
basic_matrix operator- () const
 
basic_matrix operator* (const basic_matrix &rhs) const
 Multiplies a matrix. More...
 
basic_matrix operator* (const value_type rhs) const
 Multiplies a scalar. More...
 
basic_matrix operator/ (const value_type rhs) const
 Divides a scalar. More...
 

Static Public Member Functions

static void gemm (const basic_matrix &lhs, const basic_matrix &rhs, basic_matrix &dst, const value_type alpha=value_type(1), const value_type beta=value_type(0))
 Multiplies a left matrix by a right matrix and stores the result into a destination matrix. If the beta parameter is provided, the initial values of the destination matrix are scaled by beta, and the product of the left and right matrices is added to the destination matrix. If the alpha parameter is provided, the product of the left and right matrices is scaled by alpha. More...
 
static void set_high_precision (std::ostream &out)
 Sets the specified stream to scientific notation with a high precision. More...
 

Detailed Description

template<typename TValue>
class jade::basic_matrix< TValue >

A template class for a matrix.

Definition at line 20 of file jade.matrix.hpp.

Member Typedef Documentation

◆ blas_type

template<typename TValue >
typedef basic_blas<value_type> jade::basic_matrix< TValue >::blas_type

The BLAS type.

Definition at line 27 of file jade.matrix.hpp.

◆ initializer_list_type

template<typename TValue >
typedef std::initializer_list< std::initializer_list<value_type> > jade::basic_matrix< TValue >::initializer_list_type

The initializer list type.

Definition at line 35 of file jade.matrix.hpp.

◆ lapack_type

template<typename TValue >
typedef basic_lapack<value_type> jade::basic_matrix< TValue >::lapack_type

The LAPACK type.

Definition at line 30 of file jade.matrix.hpp.

◆ value_type

template<typename TValue >
typedef TValue jade::basic_matrix< TValue >::value_type

The value type.

Definition at line 24 of file jade.matrix.hpp.

Constructor & Destructor Documentation

◆ basic_matrix() [1/6]

template<typename TValue >
jade::basic_matrix< TValue >::basic_matrix ( )
inline

Initializes a new instance of the class with no size.

Definition at line 40 of file jade.matrix.hpp.

41  : basic_matrix(0, 0)
42  {
43  }
+ Here is the caller graph for this function:

◆ basic_matrix() [2/6]

template<typename TValue >
jade::basic_matrix< TValue >::basic_matrix ( const size_t  cy,
const size_t  cx 
)
inline

Initializes a new instance of the class based on the specified width and height.

Parameters
cyThe height.
cxThe width.

Definition at line 49 of file jade.matrix.hpp.

52  : _cy (cy)
53  , _cx (cx)
54  , _m ()
55  {
56  assert((cx == 0 && cy == 0) || (cx != 0 && cy != 0));
57  _m.resize(cx * cy);
58  }

◆ basic_matrix() [3/6]

template<typename TValue >
jade::basic_matrix< TValue >::basic_matrix ( const char *const  path)
inlineexplicit

Initializes a new instance of the class based on values from the specified file.

Parameters
pathThe path to the file.

Definition at line 64 of file jade.matrix.hpp.

66  : basic_matrix()
67  {
68  read(path);
69  }

◆ basic_matrix() [4/6]

template<typename TValue >
jade::basic_matrix< TValue >::basic_matrix ( const std::string &  path)
inlineexplicit

Initializes a new instance of the class based on values from the specified file.

Parameters
pathThe path to the file.

Definition at line 75 of file jade.matrix.hpp.

77  : basic_matrix(path.c_str())
78  {
79  }

◆ basic_matrix() [5/6]

template<typename TValue >
jade::basic_matrix< TValue >::basic_matrix ( std::istream &  in)
inlineexplicit

Initializes a new instance of the class based on values from the specified input stream.

Parameters
inThe input stream.

Definition at line 85 of file jade.matrix.hpp.

87  : basic_matrix()
88  {
89  read(in);
90  }

◆ basic_matrix() [6/6]

template<typename TValue >
jade::basic_matrix< TValue >::basic_matrix ( const initializer_list_type values)
inlineexplicit

Initializes a new instance of the class based on the specified values.

Parameters
valuesThe two-dimensional values.

Definition at line 98 of file jade.matrix.hpp.

99  : _cy (values.size())
100  , _cx (_cy == 0 ? 0 : values.begin()->size())
101  , _m ()
102  {
103  assert(std::all_of(values.begin(), values.end(), [this](
104  const std::initializer_list<value_type> & row_values) -> bool {
105  return row_values.size() == _cx; }));
106 
107  _m.reserve(_cy * _cx);
108  for (const auto & row : values)
109  for (const auto n : row)
110  _m.push_back(n);
111  }

Member Function Documentation

◆ all_of()

template<typename TValue >
template<typename TPredicate >
bool jade::basic_matrix< TValue >::all_of ( const TPredicate  predicate) const
inline
Returns
True if the specified predicate is true for all values.
Parameters
predicateThe predicate function.

Definition at line 117 of file jade.matrix.hpp.

120  {
121  return std::all_of(_m.begin(), _m.end(), predicate);
122  }
+ Here is the caller graph for this function:

◆ any_of()

template<typename TValue >
template<typename TPredicate >
bool jade::basic_matrix< TValue >::any_of ( const TPredicate  predicate) const
inline
Returns
True if the specified predicate is true for any value.
Parameters
predicateThe predicate function.

Definition at line 128 of file jade.matrix.hpp.

131  {
132  return std::any_of(_m.begin(), _m.end(), predicate);
133  }

◆ clamp()

template<typename TValue >
void jade::basic_matrix< TValue >::clamp ( const value_type  min,
const value_type  max 
)
inline

Clamps all values in the matrix to the specified range.

Parameters
minThe minimum value.
maxThe maximum value.

Definition at line 138 of file jade.matrix.hpp.

141  {
142  assert(min <= max);
143 
144  auto dst_ptr = get_data();
145  const auto dst_end = dst_ptr + get_length();
146 
147  while (dst_ptr != dst_end)
148  {
149  *dst_ptr = std::min(std::max(min, *dst_ptr), max);
150  dst_ptr++;
151  }
152  }

◆ clamp_column()

template<typename TValue >
void jade::basic_matrix< TValue >::clamp_column ( const size_t  column,
const value_type  min,
const value_type  max 
)
inline

Clamps all values in a column to the specified range.

Parameters
columnThe column to clamp.
minThe minimum value.
maxThe maximum value.

Definition at line 157 of file jade.matrix.hpp.

161  {
162  assert(column < get_width());
163  assert(min <= max);
164 
165  auto dst_ptr = get_data(0, column);
166  const auto dst_end = dst_ptr + get_length();
167 
168  while (dst_ptr != dst_end)
169  {
170  *dst_ptr = std::min(std::max(min, *dst_ptr), max);
171  dst_ptr += get_width();
172  }
173  }

◆ clamp_row()

template<typename TValue >
void jade::basic_matrix< TValue >::clamp_row ( const size_t  row,
const value_type  min,
const value_type  max 
)
inline

Clamps all values in a row to the specified range.

Parameters
rowThe row to clamp.
minThe minimum value.
maxThe maximum value.

Definition at line 178 of file jade.matrix.hpp.

182  {
183  assert(row < get_height());
184  assert(min <= max);
185 
186  auto dst_ptr = get_data(row, 0);
187  const auto dst_end = dst_ptr + get_width();
188 
189  while (dst_ptr != dst_end)
190  {
191  *dst_ptr = std::min(std::max(min, *dst_ptr), max);
192  dst_ptr++;
193  }
194  }
+ Here is the caller graph for this function:

◆ contains_inf()

template<typename TValue >
bool jade::basic_matrix< TValue >::contains_inf ( ) const
inline
Returns
True if any element in the matrix is infinite.

Definition at line 199 of file jade.matrix.hpp.

200  {
201  auto src_ptr = get_data();
202  const auto src_end = src_ptr + get_length();
203 
204  while (src_ptr != src_end)
205  if (std::isinf(*src_ptr++))
206  return true;
207 
208  return false;
209  }

◆ contains_nan()

template<typename TValue >
bool jade::basic_matrix< TValue >::contains_nan ( ) const
inline
Returns
True if any element in the matrix is not a number.

Definition at line 214 of file jade.matrix.hpp.

215  {
216  auto src_ptr = get_data();
217  const auto src_end = src_ptr + get_length();
218 
219  while (src_ptr != src_end)
220  if (std::isnan(*src_ptr++))
221  return true;
222 
223  return false;
224  }

◆ copy_column() [1/2]

template<typename TValue >
basic_matrix jade::basic_matrix< TValue >::copy_column ( const size_t  column) const
inline
Returns
A new column vector based on the values from the specified column.
Parameters
columnThe column to copy.

Definition at line 230 of file jade.matrix.hpp.

233  {
234  basic_matrix out;
235  copy_column(column, out);
236  return out;
237  }
+ Here is the caller graph for this function:

◆ copy_column() [2/2]

template<typename TValue >
void jade::basic_matrix< TValue >::copy_column ( const size_t  column,
basic_matrix< TValue > &  out 
) const
inline

Copies a column into the specified vector.

Parameters
columnThe column to copy.
outThe output vector.

Definition at line 242 of file jade.matrix.hpp.

246  {
247  assert(column < get_width());
248 
249  out.resize(get_height(), 1);
250 
251  auto dst_ptr = out.get_data();
252  auto src_ptr = get_data(0, column);
253  const auto src_end = src_ptr + get_length();
254 
255  while (src_ptr != src_end)
256  {
257  *dst_ptr++ = *src_ptr;
258  src_ptr += get_width();
259  }
260  }

◆ copy_lower_to_upper()

template<typename TValue >
void jade::basic_matrix< TValue >::copy_lower_to_upper ( )
inline

Copies the elements in the lower triangle of the square matrix to the corresponding upper triangle elements.

Definition at line 266 of file jade.matrix.hpp.

267  {
268  assert(is_square());
269 
270  const auto n = get_height();
271  if (n < 2)
272  return;
273 
274  //
275  // Prepare pointers to the data.
276  //
277  const auto dst_end = get_data() + get_length();
278  auto dst_ptr = get_data() + 1;
279  auto src_end = get_data() + get_length();
280  auto src_ptr = get_data() + n;
281 
282  //
283  // Loop until the entire triangle is copied.
284  //
285  while (dst_ptr != dst_end)
286  {
287  const auto dst_anchor = dst_ptr;
288  const auto src_anchor = src_ptr;
289 
290  //
291  // Advance down the row of the lower triangle and across the
292  // column of the upper triangle.
293  //
294  while (src_ptr != src_end)
295  {
296  *dst_ptr++ = *src_ptr;
297  src_ptr += n;
298  }
299 
300  //
301  // Adjust the pointers for the next row and column.
302  //
303  dst_ptr = dst_anchor + n + 1;
304  src_ptr = src_anchor + n + 1;
305  src_end++;
306  }
307  }
+ Here is the caller graph for this function:

◆ copy_row() [1/2]

template<typename TValue >
basic_matrix jade::basic_matrix< TValue >::copy_row ( const size_t  row) const
inline
Returns
A new row vector based on the values from the specified row.
Parameters
rowThe row to copy.

Definition at line 312 of file jade.matrix.hpp.

315  {
316  basic_matrix out;
317  copy_row(row, out);
318  return out;
319  }
+ Here is the caller graph for this function:

◆ copy_row() [2/2]

template<typename TValue >
void jade::basic_matrix< TValue >::copy_row ( const size_t  row,
basic_matrix< TValue > &  out 
) const
inline

Copies a row into the specified matrix.

Parameters
rowThe row to copy.
outThe output vector.

Definition at line 324 of file jade.matrix.hpp.

328  {
329  assert(row < get_height());
330 
331  out.resize(1, get_width());
332 
333  auto dst_ptr = out.get_data();
334  auto src_ptr = get_data(row, 0);
335  const auto src_end = src_ptr + get_width();
336 
337  while (src_ptr != src_end)
338  *dst_ptr++ = *src_ptr++;
339  }

◆ copy_upper_to_lower()

template<typename TValue >
void jade::basic_matrix< TValue >::copy_upper_to_lower ( )
inline

Copies the elements in the upper triangle of the square matrix to the corresponding lower triangle elements.

Definition at line 345 of file jade.matrix.hpp.

346  {
347  assert(is_square());
348 
349  const auto n = get_height();
350  if (n < 2)
351  return;
352 
353  //
354  // Prepare pointers to the data.
355  //
356  const auto src_end = get_data() + get_length();
357  auto src_ptr = get_data() + 1;
358  auto dst_end = get_data() + get_length();
359  auto dst_ptr = get_data() + n;
360 
361  //
362  // Loop until the entire triangle is copied.
363  //
364  while (src_ptr != src_end)
365  {
366  const auto dst_anchor = dst_ptr;
367  const auto src_anchor = src_ptr;
368 
369  //
370  // Advance down the row of the lower triangle and across the
371  // column of the upper triangle.
372  //
373  while (dst_ptr != dst_end)
374  {
375  *dst_ptr = *src_ptr++;
376  dst_ptr += n;
377  }
378 
379  //
380  // Adjust the pointers for the next row and column.
381  //
382  dst_ptr = dst_anchor + n + 1;
383  src_ptr = src_anchor + n + 1;
384  dst_end++;
385  }
386  }

◆ create_transpose() [1/2]

template<typename TValue >
basic_matrix jade::basic_matrix< TValue >::create_transpose ( ) const
inline
Returns
The transpose of the matrix.

Definition at line 391 of file jade.matrix.hpp.

392  {
393  basic_matrix out;
394  create_transpose(out);
395  return out;
396  }
+ Here is the caller graph for this function:

◆ create_transpose() [2/2]

template<typename TValue >
void jade::basic_matrix< TValue >::create_transpose ( basic_matrix< TValue > &  out) const
inline

Stores the transpose of the matrix.

Parameters
outThe matrix receiving the transpose.

Definition at line 401 of file jade.matrix.hpp.

404  {
405  if (this == &out)
406  {
407  if (is_empty())
408  return;
409 
410  if (is_vector())
411  {
412  std::swap(out._cx, out._cy);
413  return;
414  }
415 
416  if (is_square())
417  {
418  for (size_t i = 1; i < _cy; i++)
419  for (size_t j = 0; j < i; j++)
420  std::swap(out(i, j), out(j, i));
421  return;
422  }
423 
424  basic_matrix tmp;
425  create_transpose(tmp);
426  out.swap(tmp);
427  return;
428  }
429 
430  if (is_empty() || is_vector() || is_square())
431  {
432  out = *this;
433  out.create_transpose(out);
434  return;
435  }
436 
437  out.resize(_cx, _cy);
438 
439  auto & src = *this;
440  for (size_t i = 0; i < _cy; i++)
441  for (size_t j = 0; j < _cx; j++)
442  out(j, i) = src(i, j);
443  }

◆ gemm()

template<typename TValue >
static void jade::basic_matrix< TValue >::gemm ( const basic_matrix< TValue > &  lhs,
const basic_matrix< TValue > &  rhs,
basic_matrix< TValue > &  dst,
const value_type  alpha = value_type(1),
const value_type  beta = value_type(0) 
)
inlinestatic

Multiplies a left matrix by a right matrix and stores the result into a destination matrix. If the beta parameter is provided, the initial values of the destination matrix are scaled by beta, and the product of the left and right matrices is added to the destination matrix. If the alpha parameter is provided, the product of the left and right matrices is scaled by alpha.

Parameters
lhsThe matrix left of the operator.
rhsThe matrix right of the operator.
dstThe matrix receiving the product.
alphaThe scalar applied to the left matrix.
betaThe scalar applied to the right matrix.

Definition at line 459 of file jade.matrix.hpp.

465  {
466  assert(lhs.get_width() == rhs.get_height());
467  assert(dst.is_size(lhs.get_height(), rhs.get_width()));
468 
470  CblasRowMajor, // Layout
471  CblasNoTrans, // transa
472  CblasNoTrans, // transb
473  int(lhs.get_height()), // m
474  int(rhs.get_width()), // n
475  int(lhs.get_width()), // k
476  alpha, // alpha
477  lhs.get_data(), // a
478  int(lhs.get_width()), // lda
479  rhs.get_data(), // b
480  int(rhs.get_width()), // ldb
481  beta, // beta
482  dst.get_data(), // c
483  int(dst.get_width())); // ldc
484  }
+ Here is the caller graph for this function:

◆ gesv()

template<typename TValue >
bool jade::basic_matrix< TValue >::gesv ( )
inline

Computes the solution to the system of linear equations with a square coefficient matrix A and multiple right-hand sides. Specifically, this method solves for X the system of linear equations A*X = B, where A is an n-by-n matrix, the columns of matrix B are individual right-hand sides, and the columns of X are the corresponding solutions.

Returns
True if successful; otherwise, false.

Definition at line 496 of file jade.matrix.hpp.

497  {
498  const auto cx = get_width();
499  const auto cy = get_height();
500 
501  assert(cx > cy);
502 
503  std::vector<int> ipiv;
504  ipiv.resize(cy);
505 
506  return 0 == lapack_type::gesv(
507  lapack_type::row_major, // layout
508  int(cy), // n
509  int(cx - cy), // nrhs
510  get_data(), // a
511  int(cx), // lda
512  ipiv.data(), // ipiv
513  get_data() + cy, // b
514  int(cx)); // ldb
515  }

◆ get_column_sum()

template<typename TValue >
value_type jade::basic_matrix< TValue >::get_column_sum ( const size_t  column) const
inline
Returns
The sum of a column.
Parameters
columnThe column to sum.

Definition at line 520 of file jade.matrix.hpp.

523  {
524  assert(column < get_width());
525 
526  auto sum = value_type(0);
527  auto ptr = get_data(0, column);
528  const auto end = ptr + get_length();
529 
530  while (ptr != end)
531  {
532  sum += *ptr;
533  ptr += get_width();
534  }
535 
536  return sum;
537  }

◆ get_data() [1/6]

template<typename TValue >
value_type* jade::basic_matrix< TValue >::get_data ( )
inline
Returns
A pointer to the first element of the data.

Definition at line 551 of file jade.matrix.hpp.

552  {
553  assert(!is_empty());
554  return _m.data();
555  }

◆ get_data() [2/6]

template<typename TValue >
const value_type* jade::basic_matrix< TValue >::get_data ( ) const
inline
Returns
A pointer to the first element of the data.

Definition at line 542 of file jade.matrix.hpp.

543  {
544  assert(!is_empty());
545  return _m.data();
546  }
+ Here is the caller graph for this function:

◆ get_data() [3/6]

template<typename TValue >
value_type* jade::basic_matrix< TValue >::get_data ( const size_t  index)
inline
Returns
A pointer to the specified element of the vector data.
Parameters
indexThe index of the element.

Definition at line 572 of file jade.matrix.hpp.

574  {
575  assert(index < get_length());
576  return _m.data() + index;
577  }

◆ get_data() [4/6]

template<typename TValue >
const value_type* jade::basic_matrix< TValue >::get_data ( const size_t  index) const
inline
Returns
A pointer to the specified element of the vector data.
Parameters
indexThe index of the element.

Definition at line 560 of file jade.matrix.hpp.

563  {
564  assert(is_vector());
565  assert(index < get_length());
566  return _m.data() + index;
567  }

◆ get_data() [5/6]

template<typename TValue >
value_type* jade::basic_matrix< TValue >::get_data ( const size_t  row,
const size_t  column 
)
inline
Returns
A pointer to the specified element of the vector data.
Parameters
rowThe row of the element.
columnThe column of the element.

Definition at line 593 of file jade.matrix.hpp.

596  {
597  return _m.data() + get_index(row, column);
598  }

◆ get_data() [6/6]

template<typename TValue >
const value_type* jade::basic_matrix< TValue >::get_data ( const size_t  row,
const size_t  column 
) const
inline
Returns
AReturns a pointer to the specified element of the vector data.
Parameters
rowThe row.
columnThe column.

Definition at line 582 of file jade.matrix.hpp.

586  {
587  return _m.data() + get_index(row, column);
588  }

◆ get_height()

template<typename TValue >
size_t jade::basic_matrix< TValue >::get_height ( ) const
inline
Returns
The height of the matrix.

Definition at line 603 of file jade.matrix.hpp.

604  {
605  return _cy;
606  }
+ Here is the caller graph for this function:

◆ get_index()

template<typename TValue >
size_t jade::basic_matrix< TValue >::get_index ( const size_t  row,
const size_t  column 
) const
inline
Returns
The element index for the specified row and column.
Parameters
rowThe row.
columnThe column.

Definition at line 611 of file jade.matrix.hpp.

615  {
616  assert(row < _cy);
617  assert(column < _cx);
618  return row * _cx + column;
619  }
+ Here is the caller graph for this function:

◆ get_length()

template<typename TValue >
size_t jade::basic_matrix< TValue >::get_length ( ) const
inline
Returns
The length of the matrix.

Definition at line 624 of file jade.matrix.hpp.

625  {
626  return _m.size();
627  }
+ Here is the caller graph for this function:

◆ get_max_value()

template<typename TValue >
value_type jade::basic_matrix< TValue >::get_max_value ( ) const
inline

Returns the maximum element value. This method should not be called if the matrix is empty.

Returns
The maximum element value.

Definition at line 635 of file jade.matrix.hpp.

636  {
637  assert(!is_empty());
638  return *std::max_element(_m.begin(), _m.end());
639  }

◆ get_min_max()

template<typename TValue >
bool jade::basic_matrix< TValue >::get_min_max ( value_type min,
value_type max 
) const
inline

Returns the minimum and maximum elements in the matrix. If the matrix is empty, this method returns false. Otherwise, this method stores the results in the specified arguments and returns true.

Parameters
minThe minimum value (output).
maxThe maximum value (output).
Returns
True if successful; otherwise, false.

Definition at line 663 of file jade.matrix.hpp.

664  {
665  if (_m.empty())
666  return false;
667 
668  const auto p = *std::minmax_element(_m.begin(), _m.end());
669  min = p.first;
670  max = p.second;
671  return true;
672  }

◆ get_min_max_column()

template<typename TValue >
bool jade::basic_matrix< TValue >::get_min_max_column ( const size_t  column,
value_type min,
value_type max 
) const
inline

Returns the minimum and maximum elements in a column. If the matrix is empty, this method returns false. Otherwise, this method stores the results in the specified arguments and returns true.

Parameters
columnThe column index.
minThe minimum value (output).
maxThe maximum value (output).
Returns
True if successful; otherwise, false.
Parameters
columnThe column index.
minThe minimum column element.
maxThe maximum column element.

Definition at line 686 of file jade.matrix.hpp.

691  {
692  assert(column < get_width());
693 
694  if (_m.empty())
695  return false;
696 
697  auto ptr = get_data(0, column);
698  auto const end = ptr + get_length();
699 
700  min = max = *ptr;
701 
702  for (ptr += get_width(); ptr < end; ptr += get_width())
703  {
704  const auto v = *ptr;
705  min = std::min(min, v);
706  max = std::max(max, v);
707  }
708 
709  return true;
710  }
+ Here is the caller graph for this function:

◆ get_min_value()

template<typename TValue >
value_type jade::basic_matrix< TValue >::get_min_value ( ) const
inline

Returns the minimum element value. This method should not be called if the matrix is empty.

Returns
The minimum element value.

Definition at line 647 of file jade.matrix.hpp.

648  {
649  assert(!is_empty());
650  return *std::min_element(_m.begin(), _m.end());
651  }

◆ get_row_sum()

template<typename TValue >
value_type jade::basic_matrix< TValue >::get_row_sum ( const size_t  row) const
inline
Returns
The sum of a row.
Parameters
rowThe row to sum.

Definition at line 716 of file jade.matrix.hpp.

719  {
720  assert(row < get_height());
721 
722  auto sum = value_type(0);
723  auto ptr = get_data(row, 0);
724  const auto end = ptr + get_width();
725 
726  while (ptr < end)
727  sum += *ptr++;
728 
729  return sum;
730  }
+ Here is the caller graph for this function:

◆ get_size_str()

template<typename TValue >
std::string jade::basic_matrix< TValue >::get_size_str ( ) const
inline
Returns
A string representation of the size of the matrix.

Definition at line 735 of file jade.matrix.hpp.

736  {
737  std::ostringstream out;
738  out << '[' << _cy << 'x' << _cx << ']';
739  return out.str();
740  }
+ Here is the caller graph for this function:

◆ get_sum()

template<typename TValue >
value_type jade::basic_matrix< TValue >::get_sum ( ) const
inline
Returns
The sum of all elements in the matrix.

Definition at line 745 of file jade.matrix.hpp.

746  {
747  return std::accumulate(
748  _m.begin(),
749  _m.end(),
750  value_type(0),
751  std::plus<value_type>());
752  }

◆ get_value() [1/2]

template<typename TValue >
value_type jade::basic_matrix< TValue >::get_value ( const size_t  index) const
inline
Returns
Thevalue at the specified index in the vector data.
Parameters
indexThe index in the vector data.

Definition at line 765 of file jade.matrix.hpp.

768  {
769  assert(index < get_length());
770  return _m[index];
771  }

◆ get_value() [2/2]

template<typename TValue >
value_type jade::basic_matrix< TValue >::get_value ( const size_t  row,
const size_t  column 
) const
inline
Returns
The value at the specified row and column.
Parameters
rowThe row index.
columnThe column index.

Definition at line 776 of file jade.matrix.hpp.

780  {
781  const auto index = get_index(row, column);
782  return _m[index];
783  }

◆ get_width()

template<typename TValue >
size_t jade::basic_matrix< TValue >::get_width ( ) const
inline
Returns
The height of the matrix.

Definition at line 757 of file jade.matrix.hpp.

758  {
759  return _cx;
760  }
+ Here is the caller graph for this function:

◆ invert() [1/2]

template<typename TValue >
bool jade::basic_matrix< TValue >::invert ( )
inline

Computes and stores the inverse of this matrix using the Cholesky square root method. Only the lower-left triangle is used to perform the calculation. If this method fails, indicating this matrix is not positive semidefinite, then the contents of this instance is left undefined, and the method returns false. Otherwise, the lower triangle is copied to the upper triangle, and the method returns true.

Returns
True if successful; otherwise, false.

Definition at line 843 of file jade.matrix.hpp.

844  {
845  value_type log_det;
846  return invert(log_det);
847  }
+ Here is the caller graph for this function:

◆ invert() [2/2]

template<typename TValue >
bool jade::basic_matrix< TValue >::invert ( value_type log_det)
inline

Computes and stores the inverse of this matrix using the Cholesky square root method. Only the lower-left triangle is used to perform the calculation. If this method fails, indicating this matrix is not positive semidefinite, then the contents of this instance and the parameter are left undefined, and the method returns false. Otherwise, the lower triangle is copied to the upper triangle, the parameter is assigned the log of the determinant of the original matrix values, and the method returns true.

Returns
True if successful; otherwise, false.
Parameters
log_detThe log of the determinant.

Definition at line 796 of file jade.matrix.hpp.

798  {
799  assert(is_square());
800 
801  //
802  // Compute the Cholesky square root. If it fails, return false to
803  // indicate this is an unacceptable set of parameters.
804  //
805  if (!potrf_lower())
806  return false;
807 
808  //
809  // Calculate the log of the determinant by summing twice the log of
810  // the diagonal entries.
811  //
812  log_det = value_type(0.0);
813  {
814  const auto end = get_data() + get_length() + _cx;
815  for (auto ptr = get_data(); ptr != end; ptr += _cx + 1)
816  log_det += value_type(2.0) * std::log(*ptr);
817  }
818 
819  //
820  // Compute the inverse. If it fails, return false to indicate
821  // this is an unacceptable set of parameters.
822  //
823  if (!potri_lower())
824  return false;
825 
826  //
827  // Mirror the values from the lower triangle to the upper triangle.
828  //
830  return true;
831  }
+ Here is the caller graph for this function:

◆ is_column_vector()

template<typename TValue >
bool jade::basic_matrix< TValue >::is_column_vector ( ) const
inline
Returns
True if this is a column vector.

Definition at line 852 of file jade.matrix.hpp.

853  {
854  return _cx == 1;
855  }
+ Here is the caller graph for this function:

◆ is_empty()

template<typename TValue >
bool jade::basic_matrix< TValue >::is_empty ( ) const
inline
Returns
True if this is empty.

Definition at line 860 of file jade.matrix.hpp.

861  {
862  return _m.empty();
863  }
+ Here is the caller graph for this function:

◆ is_length() [1/2]

template<typename TValue >
bool jade::basic_matrix< TValue >::is_length ( const basic_matrix< TValue > &  other) const
inline
Returns
True if this has the same length as a specified vector.
Parameters
otherThe other matrix.

Definition at line 876 of file jade.matrix.hpp.

879  {
880  return get_length() == other.get_length();
881  }

◆ is_length() [2/2]

template<typename TValue >
bool jade::basic_matrix< TValue >::is_length ( const size_t  length) const
inline
Returns
True if this has the specified length.
Parameters
lengthThe length to compare.

Definition at line 886 of file jade.matrix.hpp.

889  {
890  return get_length() == length;
891  }

◆ is_row_vector()

template<typename TValue >
bool jade::basic_matrix< TValue >::is_row_vector ( ) const
inline
Returns
True if this is a row vector.

Definition at line 868 of file jade.matrix.hpp.

869  {
870  return _cy == 1;
871  }
+ Here is the caller graph for this function:

◆ is_size() [1/2]

template<typename TValue >
bool jade::basic_matrix< TValue >::is_size ( const basic_matrix< TValue > &  other) const
inline
Returns
True if this has the same size as a specified matrix.
Parameters
otherThe other matrix.

Definition at line 896 of file jade.matrix.hpp.

899  {
900  return (_cx == other._cx) && (_cy == other._cy);
901  }
+ Here is the caller graph for this function:

◆ is_size() [2/2]

template<typename TValue >
bool jade::basic_matrix< TValue >::is_size ( const size_t  height,
const size_t  width 
) const
inline
Returns
True if this has the same size as the specified dimensions.
Parameters
heightThe height of the matrix.
widthThe width of the matrix.

Definition at line 906 of file jade.matrix.hpp.

910  {
911  return (_cy == height) && (_cx == width);
912  }

◆ is_square()

template<typename TValue >
bool jade::basic_matrix< TValue >::is_square ( ) const
inline
Returns
True if this matrix is square.

Definition at line 917 of file jade.matrix.hpp.

918  {
919  return _cx == _cy;
920  }
+ Here is the caller graph for this function:

◆ is_vector()

template<typename TValue >
bool jade::basic_matrix< TValue >::is_vector ( ) const
inline
Returns
True if this is a vector.

Definition at line 925 of file jade.matrix.hpp.

926  {
927  return is_column_vector() || is_row_vector();
928  }
+ Here is the caller graph for this function:

◆ multiply_column() [1/2]

template<typename TValue >
value_type jade::basic_matrix< TValue >::multiply_column ( const size_t  column,
const basic_matrix< TValue > &  vector 
)
inline
Returns
The product of the column and the vector.
Parameters
columnThe column to multiply.
vectorThe vector to multiply.

Definition at line 952 of file jade.matrix.hpp.

955  {
956  assert(column < get_height());
957  assert(vector.is_vector());
958 
959  auto a_ptr = get_data(0, column);
960  const auto a_end = a_ptr + get_length();
961  auto b_ptr = vector.get_data();
962  auto sum = value_type(0);
963 
964  while (a_ptr != a_end)
965  {
966  sum += *a_ptr * *b_ptr++;
967  a_ptr += get_width();
968  }
969 
970  return sum;
971  }

◆ multiply_column() [2/2]

template<typename TValue >
void jade::basic_matrix< TValue >::multiply_column ( const size_t  column,
const value_type  value 
)
inline

Multiplies a column by a specified value.

Parameters
columnThe column to multiply.
valueThe value to multiply.

Definition at line 933 of file jade.matrix.hpp.

936  {
937  assert(column < get_width());
938 
939  auto dst_ptr = get_data(0, column);
940  const auto dst_end = dst_ptr + get_length();
941 
942  while (dst_ptr != dst_end)
943  {
944  *dst_ptr *= value;
945  dst_ptr += get_width();
946  }
947  }

◆ multiply_row() [1/2]

template<typename TValue >
value_type jade::basic_matrix< TValue >::multiply_row ( const size_t  row,
const basic_matrix< TValue > &  vector 
) const
inline
Returns
The product of the row and the vector.
Parameters
rowThe row to multiply.
vectorThe vector to multiply.

Definition at line 992 of file jade.matrix.hpp.

996  {
997  assert(row < get_height());
998  assert(vector.is_vector());
999 
1000  auto a_ptr = get_data(row, 0);
1001  const auto a_end = a_ptr + get_width();
1002  auto b_ptr = vector.get_data();
1003  auto sum = value_type(0);
1004 
1005  while (a_ptr != a_end)
1006  sum += *a_ptr++ * *b_ptr++;
1007 
1008  return sum;
1009  }

◆ multiply_row() [2/2]

template<typename TValue >
void jade::basic_matrix< TValue >::multiply_row ( const size_t  row,
const value_type  value 
)
inline

Multiplies a row by a specified value.

Parameters
rowThe row to multiply.
valueThe value to multiply.

Definition at line 976 of file jade.matrix.hpp.

979  {
980  assert(row < get_height());
981 
982  auto dst_ptr = get_data(row, 0);
983  const auto dst_end = dst_ptr + get_width();
984 
985  while (dst_ptr != dst_end)
986  *dst_ptr++ *= value;
987  }
+ Here is the caller graph for this function:

◆ operator()() [1/2]

template<typename TValue >
value_type& jade::basic_matrix< TValue >::operator() ( const size_t  rhs_1,
const size_t  rhs_2 
)
inline
Returns
The element at the specified row and column.
Parameters
rhs_1The row.
rhs_2The column.

Definition at line 1317 of file jade.matrix.hpp.

1320  {
1321  const auto index = get_index(rhs_1, rhs_2);
1322  return _m[index];
1323  }

◆ operator()() [2/2]

template<typename TValue >
const value_type& jade::basic_matrix< TValue >::operator() ( const size_t  rhs_1,
const size_t  rhs_2 
) const
inline
Returns
The element at the specified row and column.
Parameters
rhs_1The row.
rhs_2The column.

Definition at line 1305 of file jade.matrix.hpp.

1309  {
1310  const auto index = get_index(rhs_1, rhs_2);
1311  return _m[index];
1312  }

◆ operator*() [1/2]

template<typename TValue >
basic_matrix jade::basic_matrix< TValue >::operator* ( const basic_matrix< TValue > &  rhs) const
inline

Multiplies a matrix.

Parameters
rhsThe matrix to multiply.
Returns
A new matrix.

Definition at line 1519 of file jade.matrix.hpp.

1520  {
1521  basic_matrix dst (get_height(), rhs.get_width());
1522  gemm(*this, rhs, dst);
1523  return dst;
1524  }

◆ operator*() [2/2]

template<typename TValue >
basic_matrix jade::basic_matrix< TValue >::operator* ( const value_type  rhs) const
inline

Multiplies a scalar.

Parameters
rhsThe scalar to multiply.
Returns
A new matrix.

Definition at line 1532 of file jade.matrix.hpp.

1533  {
1534  return basic_matrix(*this) *= rhs;
1535  }

◆ operator*=() [1/2]

template<typename TValue >
basic_matrix& jade::basic_matrix< TValue >::operator*= ( const basic_matrix< TValue > &  rhs)
inline

Multiplies this matrix and the specified matrix and assigns the product to this instance.

Parameters
rhsThe matrix to multiply.
Returns
This instance.

Definition at line 1408 of file jade.matrix.hpp.

1409  {
1410  basic_matrix dst (get_height(), rhs.get_width());
1411  gemm(*this, rhs, dst);
1412  swap(dst);
1413  return *this;
1414  }

◆ operator*=() [2/2]

template<typename TValue >
basic_matrix& jade::basic_matrix< TValue >::operator*= ( const value_type  rhs)
inline

Multiplies this matrix and the specified scalar and assigns the product to this instance.

Parameters
rhsThe scalar to multiply.
Returns
This instance.

Definition at line 1423 of file jade.matrix.hpp.

1424  {
1425  auto dst_ptr = get_data();
1426  const auto dst_end = dst_ptr + get_length();
1427 
1428  while (dst_ptr != dst_end)
1429  *dst_ptr++ *= rhs;
1430 
1431  return *this;
1432  }

◆ operator+() [1/2]

template<typename TValue >
basic_matrix jade::basic_matrix< TValue >::operator+ ( const basic_matrix< TValue > &  rhs) const
inline

Adds a matrix.

Parameters
rhsThe matrix to add.
Returns
A new matrix.

Definition at line 1458 of file jade.matrix.hpp.

1459  {
1460  return basic_matrix(*this) += rhs;
1461  }

◆ operator+() [2/2]

template<typename TValue >
basic_matrix jade::basic_matrix< TValue >::operator+ ( const value_type  rhs) const
inline

Adds a value.

Parameters
rhsThe value to add.
Returns
A new matrix.

Definition at line 1469 of file jade.matrix.hpp.

1470  {
1471  return basic_matrix(*this) += rhs;
1472  }

◆ operator+=() [1/2]

template<typename TValue >
basic_matrix& jade::basic_matrix< TValue >::operator+= ( const basic_matrix< TValue > &  rhs)
inline

Adds the specified matrix to this instance.

Parameters
rhsThe matrix to add.
Returns
This instance.

Definition at line 1331 of file jade.matrix.hpp.

1332  {
1333  assert(get_width() == rhs.get_width());
1334  assert(get_height() == rhs.get_height());
1335 
1336  auto dst_ptr = get_data();
1337  auto src_ptr = rhs.get_data();
1338  const auto dst_end = dst_ptr + get_length();
1339 
1340  while (dst_ptr != dst_end)
1341  *dst_ptr++ += *src_ptr++;
1342 
1343  return *this;
1344  }

◆ operator+=() [2/2]

template<typename TValue >
basic_matrix& jade::basic_matrix< TValue >::operator+= ( const value_type  rhs)
inline

Adds the specified scalar to this instance.

Parameters
rhsThe scalar to add.
Returns
This instance.

Definition at line 1352 of file jade.matrix.hpp.

1353  {
1354  auto dst_ptr = get_data();
1355  const auto dst_end = dst_ptr + get_length();
1356 
1357  while (dst_ptr != dst_end)
1358  *dst_ptr++ += rhs;
1359 
1360  return *this;
1361  }

◆ operator-() [1/3]

template<typename TValue >
basic_matrix jade::basic_matrix< TValue >::operator- ( ) const
inline
Returns
The negation of this instance.

Definition at line 1499 of file jade.matrix.hpp.

1500  {
1501  basic_matrix out (get_height(), get_width());
1502 
1503  auto dst = out.get_data();
1504  auto src = get_data();
1505  const auto end = get_data() + get_length();
1506 
1507  while (src != end)
1508  *dst++ = -(*src++);
1509 
1510  return out;
1511  }

◆ operator-() [2/3]

template<typename TValue >
basic_matrix jade::basic_matrix< TValue >::operator- ( const basic_matrix< TValue > &  rhs) const
inline

Subtracts a matrix.

Parameters
rhsThe matrix to subtract.
Returns
A new matrix.

Definition at line 1480 of file jade.matrix.hpp.

1481  {
1482  return basic_matrix(*this) -= rhs;
1483  }

◆ operator-() [3/3]

template<typename TValue >
basic_matrix jade::basic_matrix< TValue >::operator- ( const value_type  rhs) const
inline

Subtracts a value.

Parameters
rhsThe value to subtract.
Returns
A new matrix.

Definition at line 1491 of file jade.matrix.hpp.

1492  {
1493  return basic_matrix(*this) -= rhs;
1494  }

◆ operator-=() [1/2]

template<typename TValue >
basic_matrix& jade::basic_matrix< TValue >::operator-= ( const basic_matrix< TValue > &  rhs)
inline

Subtracts the specified matrix from this instance.

Parameters
rhsThe matrix to subtract.
Returns
This instance.

Definition at line 1369 of file jade.matrix.hpp.

1370  {
1371  assert(get_width() == rhs.get_width());
1372  assert(get_height() == rhs.get_height());
1373 
1374  auto dst_ptr = get_data();
1375  auto src_ptr = rhs.get_data();
1376  const auto dst_end = dst_ptr + get_length();
1377 
1378  while (dst_ptr != dst_end)
1379  *dst_ptr++ -= *src_ptr++;
1380 
1381  return *this;
1382  }

◆ operator-=() [2/2]

template<typename TValue >
basic_matrix& jade::basic_matrix< TValue >::operator-= ( const value_type  rhs)
inline

Subtracts the specified scalar from this instance.

Parameters
rhsThe scalar to subtract.
Returns
This instance.

Definition at line 1390 of file jade.matrix.hpp.

1391  {
1392  auto dst_ptr = get_data();
1393  const auto dst_end = dst_ptr + get_length();
1394 
1395  while (dst_ptr != dst_end)
1396  *dst_ptr++ -= rhs;
1397 
1398  return *this;
1399  }

◆ operator/()

template<typename TValue >
basic_matrix jade::basic_matrix< TValue >::operator/ ( const value_type  rhs) const
inline

Divides a scalar.

Parameters
rhsThe scalar to divides.
Returns
A new matrix.

Definition at line 1543 of file jade.matrix.hpp.

1544  {
1545  return basic_matrix(*this) /= rhs;
1546  }

◆ operator/=()

template<typename TValue >
basic_matrix& jade::basic_matrix< TValue >::operator/= ( const value_type  rhs)
inline

Divides this matrix by the specified scalar and assigns the quotient to this instance.

Parameters
rhsThe scalar to multiply.
Returns
This instance.

Definition at line 1441 of file jade.matrix.hpp.

1442  {
1443  auto dst_ptr = get_data();
1444  const auto dst_end = dst_ptr + get_length();
1445 
1446  while (dst_ptr != dst_end)
1447  *dst_ptr++ /= rhs;
1448 
1449  return *this;
1450  }

◆ operator[]() [1/2]

template<typename TValue >
value_type& jade::basic_matrix< TValue >::operator[] ( const size_t  rhs)
inline
Returns
The element at the specified index.
Parameters
rhsThe element index.

Definition at line 1295 of file jade.matrix.hpp.

1297  {
1298  assert(rhs < get_length());
1299  return _m[rhs];
1300  }

◆ operator[]() [2/2]

template<typename TValue >
const value_type& jade::basic_matrix< TValue >::operator[] ( const size_t  rhs) const
inline
Returns
The element at the specified index.
Parameters
rhsThe element index.

Definition at line 1284 of file jade.matrix.hpp.

1287  {
1288  assert(rhs < get_length());
1289  return _m[rhs];
1290  }

◆ potrf_lower()

template<typename TValue >
bool jade::basic_matrix< TValue >::potrf_lower ( )
inline

Forms the Cholesky factorization of a symmetric positive-definite matrix: [A = L * L^T] where L is the lower triangular portion of the matrix.

Returns
True if the method is successful.

Definition at line 1018 of file jade.matrix.hpp.

1019  {
1020  assert(is_square());
1021  assert(!is_empty());
1022 
1023  const auto n = int(get_width());
1024 
1025  return 0 == lapack_type::potrf(
1026  lapack_type::row_major, // layout
1027  'L', // uplo
1028  n, // n
1029  get_data(), // a
1030  n); // lda
1031  }
+ Here is the caller graph for this function:

◆ potri_lower()

template<typename TValue >
bool jade::basic_matrix< TValue >::potri_lower ( )
inline

Computes the inverse inv(A) of a symmetric positive definite matrix A. Before calling this routine, call dpotrf to factorize A.

Returns
True if the method is successful.

Definition at line 1039 of file jade.matrix.hpp.

1040  {
1041  assert(is_square());
1042  assert(!is_empty());
1043 
1044  const auto n = int(get_width());
1045 
1046  return 0 == lapack_type::potri(
1047  lapack_type::row_major, // layout
1048  'L', // uplo
1049  n, // n
1050  get_data(), // a
1051  n); // lda
1052  }
+ Here is the caller graph for this function:

◆ read() [1/3]

template<typename TValue >
void jade::basic_matrix< TValue >::read ( const char *const  path)
inline

Reads the matrix values from the specified file.

Exceptions
Anexception if there is an error reading the file.
Parameters
pathThe path to the file.

Definition at line 1070 of file jade.matrix.hpp.

1072  {
1073  assert(path != nullptr);
1074 
1075  std::ifstream in (path);
1076  if (!in.good())
1077  throw error() << "failed to open matrix '" << path << "'";
1078 
1079  try
1080  {
1081  read(in);
1082  }
1083  catch (const std::exception & e)
1084  {
1085  throw error()
1086  << "failed to read matrix '"
1087  << path << "': " << e.what();
1088  }
1089  }

◆ read() [2/3]

template<typename TValue >
void jade::basic_matrix< TValue >::read ( const std::string &  path)
inline

Reads the matrix values from the specified file.

Exceptions
Anexception if there is an error reading the file.
Parameters
pathThe path to the file.

Definition at line 1059 of file jade.matrix.hpp.

1061  {
1062  read(path.c_str());
1063  }
+ Here is the caller graph for this function:

◆ read() [3/3]

template<typename TValue >
void jade::basic_matrix< TValue >::read ( std::istream &  in)
inline

Reads the matrix values from the specified stream.

Exceptions
Anexception if there is an error reading the stream.
Parameters
inThe input stream.

Definition at line 1096 of file jade.matrix.hpp.

1098  {
1099  size_t cx, cy;
1100  if (!(in >> cy >> cx))
1101  throw error()
1102  << "failed to parse matrix size";
1103 
1104  if ((cx == 0 && cy != 0) || (cx != 0 && cy == 0))
1105  throw error()
1106  << "invalid matrix size ["
1107  << cy << "x" << cx << "]";
1108 
1109  vector_type m;
1110  m.reserve(cx * cy);
1111 
1112  for (size_t y = 0; y < cy; y++)
1113  {
1114  for (size_t x = 0; x < cx; x++)
1115  {
1116  value_type value;
1117  if (!(in >> value))
1118  throw error()
1119  << "failed to parse matrix value at cell ["
1120  << y+1 << "," << x+1 << "]";
1121 
1122  m.push_back(value);
1123  }
1124  }
1125 
1126  _cx = cx;
1127  _cy = cy;
1128  _m.swap(m);
1129  }

◆ resize()

template<typename TValue >
void jade::basic_matrix< TValue >::resize ( const size_t  height,
const size_t  width 
)
inline

Resizes the matrix to the specified dimensions. The values of the matrix are not reset.

Parameters
heightThe new height.
widthThe new width.

Definition at line 1135 of file jade.matrix.hpp.

1138  {
1139  assert((height == 0 && width == 0)
1140  || (height != 0 && width != 0));
1141 
1142  if (height == _cy && width == _cx)
1143  return;
1144 
1145  _cx = width;
1146  _cy = height;
1147  _m.resize(width * height);
1148  }
+ Here is the caller graph for this function:

◆ set_high_precision()

template<typename TValue >
static void jade::basic_matrix< TValue >::set_high_precision ( std::ostream &  out)
inlinestatic

Sets the specified stream to scientific notation with a high precision.

Parameters
outThe output stream.

Definition at line 1154 of file jade.matrix.hpp.

1156  {
1157  static const auto precision =
1158  1 + std::numeric_limits<value_type>::digits10;
1159 
1160  out << std::scientific << std::setprecision(precision);
1161  }

◆ set_value() [1/2]

template<typename TValue >
void jade::basic_matrix< TValue >::set_value ( const size_t  index,
const value_type value 
)
inline

Sets a single value of the vector.

Parameters
indexThe index of the value.
valueThe new value.

Definition at line 1166 of file jade.matrix.hpp.

1169  {
1170  assert(index < _m.size());
1171  _m[index] = value;
1172  }

◆ set_value() [2/2]

template<typename TValue >
void jade::basic_matrix< TValue >::set_value ( const size_t  row,
const size_t  column,
const value_type value 
)
inline

Sets a single value of the matrix.

Parameters
rowThe row to assign.
columnThe colum to assign.
valueThe value to assign.

Definition at line 1177 of file jade.matrix.hpp.

1181  {
1182  const auto index = get_index(row, column);
1183  _m[index] = value;
1184  }

◆ set_values()

template<typename TValue >
void jade::basic_matrix< TValue >::set_values ( const value_type  value)
inline

Sets all values of the matrix to the specified value.

Parameters
valueThe value to assign.

Definition at line 1189 of file jade.matrix.hpp.

1191  {
1192  for (auto & m : _m)
1193  m = value;
1194  }
+ Here is the caller graph for this function:

◆ str()

template<typename TValue >
std::string jade::basic_matrix< TValue >::str ( ) const
inline
Returns
A string representation of this instance.

Definition at line 1199 of file jade.matrix.hpp.

1200  {
1201  std::ostringstream out;
1202  write(out);
1203  return out.str();
1204  }
+ Here is the caller graph for this function:

◆ swap()

template<typename TValue >
void jade::basic_matrix< TValue >::swap ( basic_matrix< TValue > &  other)
inline

Swaps this matrix and another matrix.

Parameters
otherThe matrix to swap.

Definition at line 1209 of file jade.matrix.hpp.

1211  {
1212  std::swap(_cx, other._cx);
1213  std::swap(_cy, other._cy);
1214  _m.swap(other._m);
1215  }
+ Here is the caller graph for this function:

◆ transpose()

template<typename TValue >
void jade::basic_matrix< TValue >::transpose ( )
inline

Transposes this matrix.

Definition at line 1220 of file jade.matrix.hpp.

1221  {
1222  create_transpose(*this);
1223  }

◆ write() [1/3]

template<typename TValue >
void jade::basic_matrix< TValue >::write ( const char *const  path) const
inline

Writes this matrix to the specified file.

Exceptions
Anexception if there is an error writing to the file.
Parameters
pathThe path to the file.

Definition at line 1242 of file jade.matrix.hpp.

1245  {
1246  assert(path != nullptr);
1247 
1248  std::ofstream out (path);
1249 
1250  if (!out.good())
1251  throw error() << "failed to create matrix '" << path << "'";
1252 
1253  write(out);
1254  }

◆ write() [2/3]

template<typename TValue >
void jade::basic_matrix< TValue >::write ( const std::string &  path) const
inline

Writes this matrix to the specified file.

Exceptions
Anexception if there is an error writing to the file.
Parameters
pathThe path to the file.

Definition at line 1230 of file jade.matrix.hpp.

1233  {
1234  write(path.c_str());
1235  }
+ Here is the caller graph for this function:

◆ write() [3/3]

template<typename TValue >
void jade::basic_matrix< TValue >::write ( std::ostream &  out) const
inline

Writes this matrix to the specified output stream.

Parameters
outThe output stream.

Definition at line 1259 of file jade.matrix.hpp.

1262  {
1263  value_type const * ptr = _m.data();
1264 
1265  out << _cy << ' ' << _cx << '\n';
1266 
1267  if (_cx == 0)
1268  return;
1269 
1270  for (size_t y = 0; y < _cy; y++)
1271  {
1272  out << *ptr++;
1273 
1274  for (size_t x = 1; x < _cx; x++)
1275  out << '\t' << *ptr++;
1276 
1277  out << '\n';
1278  }
1279  }

The documentation for this class was generated from the following file:
jade::basic_matrix::is_square
bool is_square() const
Definition: jade.matrix.hpp:917
jade::basic_matrix::invert
bool invert()
Computes and stores the inverse of this matrix using the Cholesky square root method....
Definition: jade.matrix.hpp:843
jade::basic_matrix::get_data
const value_type * get_data() const
Definition: jade.matrix.hpp:542
jade::basic_matrix::get_length
size_t get_length() const
Definition: jade.matrix.hpp:624
jade::basic_matrix::get_index
size_t get_index(const size_t row, const size_t column) const
Definition: jade.matrix.hpp:611
jade::basic_matrix::get_width
size_t get_width() const
Definition: jade.matrix.hpp:757
jade::basic_lapack::potrf
static int potrf(layout_type layout, char uplo, int n, value_type *a, int lda)
Computes the Cholesky factorization of a symmetric (Hermitian) positive-definite matrix.
jade::basic_matrix::is_vector
bool is_vector() const
Definition: jade.matrix.hpp:925
jade::basic_matrix::swap
void swap(basic_matrix &other)
Swaps this matrix and another matrix.
Definition: jade.matrix.hpp:1209
jade::basic_matrix::write
void write(const std::string &path) const
Writes this matrix to the specified file.
Definition: jade.matrix.hpp:1230
jade::basic_matrix::value_type
TValue value_type
The value type.
Definition: jade.matrix.hpp:24
jade::basic_matrix::read
void read(const std::string &path)
Reads the matrix values from the specified file.
Definition: jade.matrix.hpp:1059
jade::basic_matrix::potrf_lower
bool potrf_lower()
Forms the Cholesky factorization of a symmetric positive-definite matrix: [A = L * L^T] where L is th...
Definition: jade.matrix.hpp:1018
jade::basic_lapack::potri
static int potri(layout_type layout, char uplo, int n, value_type *a, int lda)
Computes the inverse of a symmetric (Hermitian) positive-definite matrix using the Cholesky factoriza...
jade::basic_matrix::is_empty
bool is_empty() const
Definition: jade.matrix.hpp:860
jade::basic_matrix::basic_matrix
basic_matrix()
Initializes a new instance of the class with no size.
Definition: jade.matrix.hpp:40
jade::basic_matrix::get_height
size_t get_height() const
Definition: jade.matrix.hpp:603
jade::basic_lapack::gesv
static int gesv(layout_type layout, int n, int nrhs, value_type *a, int lda, int *ipiv, value_type *b, int ldb)
Computes the solution to the system of linear equations with a square coefficient matrix A and multip...
jade::basic_matrix::is_row_vector
bool is_row_vector() const
Definition: jade.matrix.hpp:868
jade::basic_matrix::is_column_vector
bool is_column_vector() const
Definition: jade.matrix.hpp:852
jade::basic_matrix::copy_lower_to_upper
void copy_lower_to_upper()
Copies the elements in the lower triangle of the square matrix to the corresponding upper triangle el...
Definition: jade.matrix.hpp:266
jade::basic_blas::gemm
static void gemm(const layout_type Order, const transpose_type TransA, const transpose_type TransB, const int M, const int N, const int K, const value_type alpha, const value_type *A, const int lda, const value_type *B, const int ldb, const value_type beta, value_type *C, const int ldc)
Computes a matrix-matrix product with general matrices.
jade::basic_matrix::copy_column
basic_matrix copy_column(const size_t column) const
Definition: jade.matrix.hpp:230
jade::basic_matrix::copy_row
basic_matrix copy_row(const size_t row) const
Definition: jade.matrix.hpp:312
jade::basic_matrix::potri_lower
bool potri_lower()
Computes the inverse inv(A) of a symmetric positive definite matrix A. Before calling this routine,...
Definition: jade.matrix.hpp:1039
jade::basic_matrix::gemm
static void gemm(const basic_matrix &lhs, const basic_matrix &rhs, basic_matrix &dst, const value_type alpha=value_type(1), const value_type beta=value_type(0))
Multiplies a left matrix by a right matrix and stores the result into a destination matrix....
Definition: jade.matrix.hpp:459
jade::basic_matrix::create_transpose
basic_matrix create_transpose() const
Definition: jade.matrix.hpp:391
jade::basic_lapack::row_major
@ row_major
A row-major order.
Definition: jade.lapack.hpp:37