oomph::DenseDoubleMatrix Class Reference

#include <matrices.h>

+ Inheritance diagram for oomph::DenseDoubleMatrix:

Public Member Functions

 DenseDoubleMatrix ()
 Constructor, set the default linear solver. More...
 
 DenseDoubleMatrix (const unsigned long &n)
 Constructor to build a square n by n matrix. More...
 
 DenseDoubleMatrix (const unsigned long &n, const unsigned long &m)
 Constructor to build a matrix with n rows and m columns. More...
 
 DenseDoubleMatrix (const unsigned long &n, const unsigned long &m, const double &initial_val)
 
 DenseDoubleMatrix (const DenseDoubleMatrix &matrix)=delete
 Broken copy constructor. More...
 
void operator= (const DenseDoubleMatrix &)=delete
 Broken assignment operator. More...
 
unsigned long nrow () const
 Return the number of rows of the matrix. More...
 
unsigned long ncol () const
 Return the number of columns of the matrix. More...
 
double operator() (const unsigned long &i, const unsigned long &j) const
 
doubleoperator() (const unsigned long &i, const unsigned long &j)
 
virtual ~DenseDoubleMatrix ()
 Destructor. More...
 
virtual void ludecompose ()
 LU decomposition using DenseLU (default linea solver) More...
 
virtual void lubksub (DoubleVector &rhs)
 LU backsubstitution. More...
 
virtual void lubksub (Vector< double > &rhs)
 LU backsubstitution. More...
 
void eigenvalues_by_jacobi (Vector< double > &eigen_val, DenseMatrix< double > &eigen_vect) const
 
void multiply (const DoubleVector &x, DoubleVector &soln) const
 Multiply the matrix by the vector x: soln=Ax. More...
 
void multiply_transpose (const DoubleVector &x, DoubleVector &soln) const
 Multiply the transposed matrix by the vector x: soln=A^T x. More...
 
void matrix_reduction (const double &alpha, DenseDoubleMatrix &reduced_matrix)
 
void multiply (const DenseDoubleMatrix &matrix_in, DenseDoubleMatrix &result)
 Function to multiply this matrix by a DenseDoubleMatrix matrix_in. More...
 
- Public Member Functions inherited from oomph::DoubleMatrixBase
 DoubleMatrixBase ()
 (Empty) constructor. More...
 
 DoubleMatrixBase (const DoubleMatrixBase &matrix)=delete
 Broken copy constructor. More...
 
void operator= (const DoubleMatrixBase &)=delete
 Broken assignment operator. More...
 
virtual ~DoubleMatrixBase ()
 virtual (empty) destructor More...
 
LinearSolver *& linear_solver_pt ()
 Return a pointer to the linear solver object. More...
 
LinearSolver *const & linear_solver_pt () const
 Return a pointer to the linear solver object (const version) More...
 
void solve (DoubleVector &rhs)
 
void solve (const DoubleVector &rhs, DoubleVector &soln)
 
void solve (Vector< double > &rhs)
 
void solve (const Vector< double > &rhs, Vector< double > &soln)
 
virtual void residual (const DoubleVector &x, const DoubleVector &b, DoubleVector &residual_)
 Find the residual, i.e. r=b-Ax the residual. More...
 
virtual double max_residual (const DoubleVector &x, const DoubleVector &rhs)
 
- Public Member Functions inherited from oomph::DenseMatrix< double >
 DenseMatrix ()
 Empty constructor, simply assign the lengths N and M to 0. More...
 
 DenseMatrix (const DenseMatrix &source_matrix)
 Copy constructor: Deep copy! More...
 
 DenseMatrix (const unsigned long &n)
 Constructor to build a square n by n matrix. More...
 
 DenseMatrix (const unsigned long &n, const unsigned long &m)
 Constructor to build a matrix with n rows and m columns. More...
 
 DenseMatrix (const unsigned long &n, const unsigned long &m, const double &initial_val)
 
DenseMatrixoperator= (const DenseMatrix &source_matrix)
 Copy assignment. More...
 
doubleentry (const unsigned long &i, const unsigned long &j)
 
double get_entry (const unsigned long &i, const unsigned long &j) const
 
virtual ~DenseMatrix ()
 Destructor, clean up the matrix data. More...
 
unsigned long nrow () const
 Return the number of rows of the matrix. More...
 
unsigned long ncol () const
 Return the number of columns of the matrix. More...
 
void resize (const unsigned long &n)
 
void resize (const unsigned long &n, const unsigned long &m)
 
void resize (const unsigned long &n, const unsigned long &m, const double &initial_value)
 
void initialise (const double &val)
 Initialize all values in the matrix to val. More...
 
void output (std::ostream &outfile) const
 Output function to print a matrix row-by-row to the stream outfile. More...
 
void output (std::string filename) const
 Output function to print a matrix row-by-row to a file. Specify filename. More...
 
void indexed_output (std::ostream &outfile) const
 Indexed output as i,j,a(i,j) More...
 
void indexed_output (std::string filename) const
 
void output_bottom_right_zero_helper (std::ostream &outfile) const
 
void sparse_indexed_output_helper (std::ostream &outfile) const
 Sparse indexed output as i,j,a(i,j) for a(i,j)!=0 only. More...
 
- Public Member Functions inherited from oomph::Matrix< T, MATRIX_TYPE >
 Matrix ()
 (Empty) constructor More...
 
 Matrix (const Matrix &matrix)=delete
 Broken copy constructor. More...
 
void operator= (const Matrix &)=delete
 Broken assignment operator. More...
 
virtual ~Matrix ()
 Virtual (empty) destructor. More...
 
T operator() (const unsigned long &i, const unsigned long &j) const
 
Toperator() (const unsigned long &i, const unsigned long &j)
 
void sparse_indexed_output (std::ostream &outfile, const unsigned &precision=0, const bool &output_bottom_right_zero=false) const
 
void sparse_indexed_output (std::string filename, const unsigned &precision=0, const bool &output_bottom_right_zero=false) const
 

Additional Inherited Members

- Protected Member Functions inherited from oomph::Matrix< T, MATRIX_TYPE >
void range_check (const unsigned long &i, const unsigned long &j) const
 
- Protected Attributes inherited from oomph::DoubleMatrixBase
LinearSolverLinear_solver_pt
 
LinearSolverDefault_linear_solver_pt
 
- Protected Attributes inherited from oomph::DenseMatrix< double >
doubleMatrixdata
 Internal representation of matrix as a pointer to data. More...
 
unsigned long N
 Number of rows. More...
 
unsigned long M
 Number of columns. More...
 

Detailed Description

Class of matrices containing doubles, and stored as a DenseMatrix<double>, but with solving functionality inherited from the abstract DoubleMatrix class.

Constructor & Destructor Documentation

◆ DenseDoubleMatrix() [1/5]

oomph::DenseDoubleMatrix::DenseDoubleMatrix ( )

Constructor, set the default linear solver.

///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// Constructor, set the default linear solver to be the DenseLU solver

140  {
142  }
LinearSolver * Default_linear_solver_pt
Definition: matrices.h:267
LinearSolver * Linear_solver_pt
Definition: matrices.h:264

References oomph::DoubleMatrixBase::Default_linear_solver_pt, and oomph::DoubleMatrixBase::Linear_solver_pt.

◆ DenseDoubleMatrix() [2/5]

oomph::DenseDoubleMatrix::DenseDoubleMatrix ( const unsigned long &  n)

Constructor to build a square n by n matrix.

Constructor to build a square n by n matrix. Set the default linear solver to be DenseLU

150  {
152  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11

References oomph::DoubleMatrixBase::Default_linear_solver_pt, and oomph::DoubleMatrixBase::Linear_solver_pt.

◆ DenseDoubleMatrix() [3/5]

oomph::DenseDoubleMatrix::DenseDoubleMatrix ( const unsigned long &  n,
const unsigned long &  m 
)

Constructor to build a matrix with n rows and m columns.

Constructor to build a matrix with n rows and m columns. Set the default linear solver to be DenseLU

162  {
164  }
int * m
Definition: level2_cplx_impl.h:294

References oomph::DoubleMatrixBase::Default_linear_solver_pt, and oomph::DoubleMatrixBase::Linear_solver_pt.

◆ DenseDoubleMatrix() [4/5]

oomph::DenseDoubleMatrix::DenseDoubleMatrix ( const unsigned long &  n,
const unsigned long &  m,
const double initial_val 
)

Constructor to build a matrix with n rows and m columns, with initial value initial_val

Constructor to build a matrix with n rows and m columns, with initial value initial_val Set the default linear solver to be DenseLU

174  : DenseMatrix<double>(n, m, initial_val)
175  {
177  }

References oomph::DoubleMatrixBase::Default_linear_solver_pt, and oomph::DoubleMatrixBase::Linear_solver_pt.

◆ DenseDoubleMatrix() [5/5]

oomph::DenseDoubleMatrix::DenseDoubleMatrix ( const DenseDoubleMatrix matrix)
delete

Broken copy constructor.

◆ ~DenseDoubleMatrix()

oomph::DenseDoubleMatrix::~DenseDoubleMatrix ( )
virtual

Destructor.

Destructor delete the default linear solver.

183  {
184  // Delete the default linear solver
186  }

References oomph::DoubleMatrixBase::Default_linear_solver_pt.

Member Function Documentation

◆ eigenvalues_by_jacobi()

void oomph::DenseDoubleMatrix::eigenvalues_by_jacobi ( Vector< double > &  eigen_vals,
DenseMatrix< double > &  eigen_vect 
) const

Determine eigenvalues and eigenvectors, using Jacobi rotations. Only for symmetric matrices. Nothing gets overwritten!

  • eigen_vect(i,j) = j-th component of i-th eigenvector.
  • eigen_val(i) is the i-th eigenvalue; same ordering as in eigenvectors

Determine eigenvalues and eigenvectors, using Jacobi rotations. Only for symmetric matrices. Nothing gets overwritten!

  • eigen_vect(i,j) = j-th component of i-th eigenvector.
  • eigen_val[i] is the i-th eigenvalue; same ordering as in eigenvectors
226  {
227 #ifdef PARANOID
228  // Check Matrix is square
229  if (N != M)
230  {
231  throw OomphLibError(
232  "This matrix is not square, the matrix MUST be square!",
235  }
236 #endif
237  // Make a copy of the matrix & check that it's symmetric
238 
239  // Check that the sizes of eigen_vals and eigen_vect are correct. If not
240  // correct them.
241  if (eigen_vals.size() != N)
242  {
243  eigen_vals.resize(N);
244  }
245  if (eigen_vect.ncol() != N || eigen_vect.nrow() != N)
246  {
247  eigen_vect.resize(N);
248  }
249 
250  DenseDoubleMatrix working_matrix(N);
251  for (unsigned long i = 0; i < N; i++)
252  {
253  for (unsigned long j = 0; j < M; j++)
254  {
255 #ifdef PARANOID
256  if (Matrixdata[M * i + j] != Matrixdata[M * j + i])
257  {
258  throw OomphLibError(
259  "Matrix needs to be symmetric for eigenvalues_by_jacobi()",
262  }
263 #endif
264  working_matrix(i, j) = (*this)(i, j);
265  }
266  }
267 
268  DenseDoubleMatrix aux_eigen_vect(N);
269 
270  throw OomphLibError("Sorry JacobiEigenSolver::jacobi() removed because of "
271  "licencing problems.",
274 
275  // // Call eigensolver
276  // unsigned long nrot;
277  // JacobiEigenSolver::jacobi(working_matrix, eigen_vals, aux_eigen_vect,
278  // nrot);
279 
280  // Copy across (and transpose)
281  for (unsigned long i = 0; i < N; i++)
282  {
283  for (unsigned long j = 0; j < M; j++)
284  {
285  eigen_vect(i, j) = aux_eigen_vect(j, i);
286  }
287  }
288  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
DenseDoubleMatrix()
Constructor, set the default linear solver.
Definition: matrices.cc:139
unsigned long nrow() const
Return the number of rows of the matrix.
Definition: matrices.h:485
unsigned long N
Number of rows.
Definition: matrices.h:392
double * Matrixdata
Internal representation of matrix as a pointer to data.
Definition: matrices.h:389
unsigned long ncol() const
Return the number of columns of the matrix.
Definition: matrices.h:491
unsigned long M
Number of columns.
Definition: matrices.h:395
void resize(const unsigned long &n)
Definition: matrices.h:498
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References i, j, oomph::DenseMatrix< double >::M, oomph::DenseMatrix< double >::Matrixdata, oomph::DenseMatrix< double >::N, oomph::DenseMatrix< T >::ncol(), oomph::DenseMatrix< T >::nrow(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::DenseMatrix< T >::resize().

◆ lubksub() [1/2]

void oomph::DenseDoubleMatrix::lubksub ( DoubleVector rhs)
virtual

LU backsubstitution.

Back substitute an LU decomposed matrix.

203  {
204  // Use the default (DenseLU) solver to perform the backsubstitution
205  static_cast<DenseLU*>(Default_linear_solver_pt)->backsub(rhs, rhs);
206  }

References oomph::DoubleMatrixBase::Default_linear_solver_pt.

Referenced by oomph::DGElement::get_inverse_mass_matrix_times_residuals(), oomph::Z2ErrorEstimator::get_recovered_flux_in_patch(), and oomph::VorticitySmoother< ELEMENT >::get_recovered_vorticity_in_patch().

◆ lubksub() [2/2]

void oomph::DenseDoubleMatrix::lubksub ( Vector< double > &  rhs)
virtual

LU backsubstitution.

Back substitute an LU decomposed matrix.

212  {
213  // Use the default (DenseLU) solver to perform the backsubstitution
214  static_cast<DenseLU*>(Default_linear_solver_pt)->backsub(rhs, rhs);
215  }

References oomph::DoubleMatrixBase::Default_linear_solver_pt.

◆ ludecompose()

void oomph::DenseDoubleMatrix::ludecompose ( )
virtual

LU decomposition using DenseLU (default linea solver)

LU decompose a matrix, by using the default linear solver (DenseLU)

193  {
194  // Use the default (DenseLU) solver to ludecompose the matrix
195  static_cast<DenseLU*>(Default_linear_solver_pt)->factorise(this);
196  }

References oomph::DoubleMatrixBase::Default_linear_solver_pt.

Referenced by oomph::DGElement::get_inverse_mass_matrix_times_residuals(), oomph::Z2ErrorEstimator::get_recovered_flux_in_patch(), oomph::VorticitySmoother< ELEMENT >::get_recovered_vorticity_in_patch(), and oomph::DGElement::pre_compute_mass_matrix().

◆ matrix_reduction()

void oomph::DenseDoubleMatrix::matrix_reduction ( const double alpha,
DenseDoubleMatrix reduced_matrix 
)

For every row, find the maximum absolute value of the entries in this row. Set all values that are less than alpha times this maximum to zero and return the resulting matrix in reduced_matrix. Note: Diagonal entries are retained regardless of their size.

483  {
484  reduced_matrix.resize(N, M, 0.0);
485  // maximum value in a row
486  double max_row;
487 
488  // Loop over rows
489  for (unsigned i = 0; i < N; i++)
490  {
491  // Initialise max value in row
492  max_row = 0.0;
493 
494  // Loop over entries in columns
495  for (unsigned long j = 0; j < M; j++)
496  {
497  // Find max. value in row
498  if (std::fabs(Matrixdata[M * i + j]) > max_row)
499  {
500  max_row = std::fabs(Matrixdata[M * i + j]);
501  }
502  }
503 
504  // Decide if we need to retain the entries in the row
505  for (unsigned long j = 0; j < M; j++)
506  {
507  // If we're on the diagonal or the value is sufficiently large: retain
508  // i.e. copy across.
509  if (i == j || std::fabs(Matrixdata[M * i + j]) > alpha * max_row)
510  {
511  reduced_matrix(i, j) = Matrixdata[M * i + j];
512  }
513  }
514  }
515  }
RealScalar alpha
Definition: level1_cplx_impl.h:151
Real fabs(const Real &a)
Definition: boostmultiprec.cpp:117

References alpha, boost::multiprecision::fabs(), i, j, oomph::DenseMatrix< double >::M, oomph::DenseMatrix< double >::Matrixdata, oomph::DenseMatrix< double >::N, and oomph::DenseMatrix< T >::resize().

◆ multiply() [1/2]

void oomph::DenseDoubleMatrix::multiply ( const DenseDoubleMatrix matrix_in,
DenseDoubleMatrix result 
)

Function to multiply this matrix by a DenseDoubleMatrix matrix_in.

Function to multiply this matrix by the DenseDoubleMatrix matrix_in.

523  {
524 #ifdef PARANOID
525  // check matrix dimensions are compatable
526  if (this->ncol() != matrix_in.nrow())
527  {
528  std::ostringstream error_message;
529  error_message
530  << "Matrix dimensions incompatable for matrix-matrix multiplication"
531  << "ncol() for first matrix:" << this->ncol()
532  << "nrow() for second matrix: " << matrix_in.nrow();
533 
534  throw OomphLibError(
535  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
536  }
537 #endif
538 
539  // NB N is number of rows!
540  unsigned long n_row = this->nrow();
541  unsigned long m_col = matrix_in.ncol();
542 
543  // resize and intialize result
544  result.resize(n_row, m_col, 0.0);
545 
546  // clock_t clock1 = clock();
547 
548  // do calculation
549  unsigned long n_col = this->ncol();
550  for (unsigned long k = 0; k < n_col; k++)
551  {
552  for (unsigned long i = 0; i < n_row; i++)
553  {
554  for (unsigned long j = 0; j < m_col; j++)
555  {
556  result(i, j) += Matrixdata[m_col * i + k] * matrix_in(k, j);
557  }
558  }
559  }
560  }
unsigned long nrow() const
Return the number of rows of the matrix.
Definition: matrices.h:1295
unsigned long ncol() const
Return the number of columns of the matrix.
Definition: matrices.h:1301
char char char int int * k
Definition: level2_impl.h:374

References i, j, k, oomph::DenseMatrix< double >::Matrixdata, ncol(), nrow(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::DenseMatrix< T >::resize().

◆ multiply() [2/2]

void oomph::DenseDoubleMatrix::multiply ( const DoubleVector x,
DoubleVector soln 
) const
virtual

Multiply the matrix by the vector x: soln=Ax.

Implements oomph::DoubleMatrixBase.

296  {
297 #ifdef PARANOID
298  // Check to see if x.size() = ncol().
299  if (x.nrow() != this->ncol())
300  {
301  std::ostringstream error_message_stream;
302  error_message_stream << "The x vector is not the right size. It is "
303  << x.nrow() << ", it should be " << this->ncol()
304  << std::endl;
305  throw OomphLibError(error_message_stream.str(),
308  }
309  // check that x is not distributed
310  if (x.distributed())
311  {
312  std::ostringstream error_message_stream;
313  error_message_stream
314  << "The x vector cannot be distributed for DenseDoubleMatrix "
315  << "matrix-vector multiply" << std::endl;
316  throw OomphLibError(error_message_stream.str(),
319  }
320  // if soln is setup...
321  if (soln.built())
322  {
323  // check that soln is not distributed
324  if (soln.distributed())
325  {
326  std::ostringstream error_message_stream;
327  error_message_stream
328  << "The x vector cannot be distributed for DenseDoubleMatrix "
329  << "matrix-vector multiply" << std::endl;
330  throw OomphLibError(error_message_stream.str(),
333  }
334  if (soln.nrow() != this->nrow())
335  {
336  std::ostringstream error_message_stream;
337  error_message_stream
338  << "The soln vector is setup and therefore must have the same "
339  << "number of rows as the matrix";
340  throw OomphLibError(error_message_stream.str(),
343  }
344  if (*x.distribution_pt()->communicator_pt() !=
345  *soln.distribution_pt()->communicator_pt())
346  {
347  std::ostringstream error_message_stream;
348  error_message_stream
349  << "The soln vector and the x vector must have the same communicator"
350  << std::endl;
351  throw OomphLibError(error_message_stream.str(),
354  }
355  }
356 #endif
357 
358  // if soln is not setup then setup the distribution
359  if (!soln.built())
360  {
361  LinearAlgebraDistribution dist(
362  x.distribution_pt()->communicator_pt(), this->nrow(), false);
363  soln.build(&dist, 0.0);
364  }
365 
366  // Initialise the solution
367  soln.initialise(0.0);
368 
369  // Multiply the matrix A, by the vector x
370  const double* x_pt = x.values_pt();
371  double* soln_pt = soln.values_pt();
372  for (unsigned long i = 0; i < N; i++)
373  {
374  for (unsigned long j = 0; j < M; j++)
375  {
376  soln_pt[i] += Matrixdata[M * i + j] * x_pt[j];
377  }
378  }
379  }
list x
Definition: plotDoE.py:28

References oomph::DoubleVector::build(), oomph::DoubleVector::built(), oomph::LinearAlgebraDistribution::communicator_pt(), oomph::DistributableLinearAlgebraObject::distributed(), oomph::DistributableLinearAlgebraObject::distribution_pt(), i, oomph::DoubleVector::initialise(), j, oomph::DenseMatrix< double >::M, oomph::DenseMatrix< double >::Matrixdata, oomph::DenseMatrix< double >::N, ncol(), oomph::DistributableLinearAlgebraObject::nrow(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::DoubleVector::values_pt(), and plotDoE::x.

Referenced by main().

◆ multiply_transpose()

void oomph::DenseDoubleMatrix::multiply_transpose ( const DoubleVector x,
DoubleVector soln 
) const
virtual

Multiply the transposed matrix by the vector x: soln=A^T x.

Implements oomph::DoubleMatrixBase.

387  {
388 #ifdef PARANOID
389  // Check to see if x.size() = ncol().
390  if (x.nrow() != this->nrow())
391  {
392  std::ostringstream error_message_stream;
393  error_message_stream << "The x vector is not the right size. It is "
394  << x.nrow() << ", it should be " << this->nrow()
395  << std::endl;
396  throw OomphLibError(error_message_stream.str(),
399  }
400  // check that x is not distributed
401  if (x.distributed())
402  {
403  std::ostringstream error_message_stream;
404  error_message_stream
405  << "The x vector cannot be distributed for DenseDoubleMatrix "
406  << "matrix-vector multiply" << std::endl;
407  throw OomphLibError(error_message_stream.str(),
410  }
411  // if soln is setup...
412  if (soln.built())
413  {
414  // check that soln is not distributed
415  if (soln.distributed())
416  {
417  std::ostringstream error_message_stream;
418  error_message_stream
419  << "The x vector cannot be distributed for DenseDoubleMatrix "
420  << "matrix-vector multiply" << std::endl;
421  throw OomphLibError(error_message_stream.str(),
424  }
425  if (soln.nrow() != this->ncol())
426  {
427  std::ostringstream error_message_stream;
428  error_message_stream
429  << "The soln vector is setup and therefore must have the same "
430  << "number of columns as the matrix";
431  throw OomphLibError(error_message_stream.str(),
434  }
435  if (*soln.distribution_pt()->communicator_pt() !=
436  *x.distribution_pt()->communicator_pt())
437  {
438  std::ostringstream error_message_stream;
439  error_message_stream
440  << "The soln vector and the x vector must have the same communicator"
441  << std::endl;
442  throw OomphLibError(error_message_stream.str(),
445  }
446  }
447 #endif
448 
449  // if soln is not setup then setup the distribution
450  if (!soln.built())
451  {
452  LinearAlgebraDistribution* dist_pt = new LinearAlgebraDistribution(
453  x.distribution_pt()->communicator_pt(), this->ncol(), false);
454  soln.build(dist_pt, 0.0);
455  delete dist_pt;
456  }
457 
458  // Initialise the solution
459  soln.initialise(0.0);
460 
461  // Matrix vector product
462  double* soln_pt = soln.values_pt();
463  const double* x_pt = x.values_pt();
464  for (unsigned long i = 0; i < N; i++)
465  {
466  for (unsigned long j = 0; j < M; j++)
467  {
468  soln_pt[j] += Matrixdata[N * i + j] * x_pt[i];
469  }
470  }
471  }

References oomph::DoubleVector::build(), oomph::DoubleVector::built(), oomph::LinearAlgebraDistribution::communicator_pt(), oomph::DistributableLinearAlgebraObject::distributed(), oomph::DistributableLinearAlgebraObject::distribution_pt(), i, oomph::DoubleVector::initialise(), j, oomph::DenseMatrix< double >::M, oomph::DenseMatrix< double >::Matrixdata, oomph::DenseMatrix< double >::N, oomph::DistributableLinearAlgebraObject::nrow(), nrow(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::DoubleVector::values_pt(), and plotDoE::x.

◆ ncol()

unsigned long oomph::DenseDoubleMatrix::ncol ( ) const
inlinevirtual

Return the number of columns of the matrix.

Implements oomph::DoubleMatrixBase.

1302  {
1303  return DenseMatrix<double>::ncol();
1304  }

References oomph::DenseMatrix< T >::ncol().

Referenced by main(), and multiply().

◆ nrow()

unsigned long oomph::DenseDoubleMatrix::nrow ( ) const
inlinevirtual

Return the number of rows of the matrix.

Implements oomph::DoubleMatrixBase.

1296  {
1297  return DenseMatrix<double>::nrow();
1298  }

References oomph::DenseMatrix< T >::nrow().

Referenced by main(), multiply(), and multiply_transpose().

◆ operator()() [1/2]

double& oomph::DenseDoubleMatrix::operator() ( const unsigned long &  i,
const unsigned long &  j 
)
inline

Overload the non-const version of the round-bracket access operator for read-write access

1317  {
1318  return DenseMatrix<double>::entry(i, j);
1319  }
T & entry(const unsigned long &i, const unsigned long &j)
Definition: matrices.h:447

References oomph::DenseMatrix< T >::entry(), i, and j.

◆ operator()() [2/2]

double oomph::DenseDoubleMatrix::operator() ( const unsigned long &  i,
const unsigned long &  j 
) const
inlinevirtual

Overload the const version of the round-bracket access operator for read-only access.

Implements oomph::DoubleMatrixBase.

1310  {
1312  }
T get_entry(const unsigned long &i, const unsigned long &j) const
Definition: matrices.h:457

References oomph::DenseMatrix< T >::get_entry(), i, and j.

◆ operator=()

void oomph::DenseDoubleMatrix::operator= ( const DenseDoubleMatrix )
delete

Broken assignment operator.


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