#include <AutoDiffMatrix.hpp>

Public Types

typedef std::vector< double > DiagRep
 
typedef Eigen::SparseMatrix
< double > 
SparseRep
 

Public Member Functions

 AutoDiffMatrix ()
 
 AutoDiffMatrix (const int num_rows, const int num_cols)
 
 AutoDiffMatrix (const Eigen::DiagonalMatrix< double, Eigen::Dynamic > &d)
 
 AutoDiffMatrix (const Eigen::SparseMatrix< double > &s)
 
 AutoDiffMatrix (const AutoDiffMatrix &other)=default
 
AutoDiffMatrixoperator= (const AutoDiffMatrix &other)=default
 
 AutoDiffMatrix (AutoDiffMatrix &&other)
 
AutoDiffMatrixoperator= (AutoDiffMatrix &&other)
 
void swap (AutoDiffMatrix &other)
 
AutoDiffMatrix operator+ (const AutoDiffMatrix &rhs) const
 
AutoDiffMatrix operator* (const AutoDiffMatrix &rhs) const
 
AutoDiffMatrixoperator+= (const AutoDiffMatrix &rhs)
 
AutoDiffMatrixoperator-= (const AutoDiffMatrix &rhs)
 
AutoDiffMatrix operator* (const double rhs) const
 
AutoDiffMatrix operator/ (const double rhs) const
 
Eigen::VectorXd operator* (const Eigen::VectorXd &rhs) const
 
template<class Scalar , int Options, class Index >
void toSparse (Eigen::SparseMatrix< Scalar, Options, Index > &s) const
 
int rows () const
 
int cols () const
 
int nonZeros () const
 
double coeff (const int row, const int col) const
 
const SparseRepgetSparse () const
 

Static Public Member Functions

static AutoDiffMatrix createIdentity (const int num_rows_cols)
 
static AutoDiffMatrix addII (const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
 
static AutoDiffMatrix addDI (const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
 
static AutoDiffMatrix addDD (const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
 
static AutoDiffMatrix addSI (const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
 
static AutoDiffMatrix addSD (const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
 
static AutoDiffMatrix addSS (const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
 
static AutoDiffMatrix mulDD (const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
 
static AutoDiffMatrix mulDS (const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
 
static AutoDiffMatrix mulSD (const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
 
static AutoDiffMatrix mulSS (const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
 

Detailed Description

AutoDiffMatrix is a wrapper class that optimizes matrix operations. Internally, an AutoDiffMatrix can be either Zero, Identity, Diagonal, or Sparse, and we utilize this to perform faster matrix operations.

Member Typedef Documentation

typedef std::vector<double> Opm::AutoDiffMatrix::DiagRep
typedef Eigen::SparseMatrix<double> Opm::AutoDiffMatrix::SparseRep

Constructor & Destructor Documentation

Opm::AutoDiffMatrix::AutoDiffMatrix ( )
inline

Creates an empty zero matrix

Referenced by createIdentity(), and operator*().

Opm::AutoDiffMatrix::AutoDiffMatrix ( const int  num_rows,
const int  num_cols 
)
inline

Creates a zero matrix with num_rows x num_cols entries

Opm::AutoDiffMatrix::AutoDiffMatrix ( const Eigen::DiagonalMatrix< double, Eigen::Dynamic > &  d)
inlineexplicit

Creates a diagonal matrix from an Eigen diagonal matrix

Opm::AutoDiffMatrix::AutoDiffMatrix ( const Eigen::SparseMatrix< double > &  s)
inlineexplicit

Creates a sparse matrix from an Eigen sparse matrix

Opm::AutoDiffMatrix::AutoDiffMatrix ( const AutoDiffMatrix other)
default
Opm::AutoDiffMatrix::AutoDiffMatrix ( AutoDiffMatrix &&  other)
inline

References swap().

Member Function Documentation

static AutoDiffMatrix Opm::AutoDiffMatrix::addDD ( const AutoDiffMatrix lhs,
const AutoDiffMatrix rhs 
)
inlinestatic

Referenced by operator+().

static AutoDiffMatrix Opm::AutoDiffMatrix::addDI ( const AutoDiffMatrix lhs,
const AutoDiffMatrix rhs 
)
inlinestatic

Referenced by operator+().

static AutoDiffMatrix Opm::AutoDiffMatrix::addII ( const AutoDiffMatrix lhs,
const AutoDiffMatrix rhs 
)
inlinestatic

Referenced by operator+().

static AutoDiffMatrix Opm::AutoDiffMatrix::addSD ( const AutoDiffMatrix lhs,
const AutoDiffMatrix rhs 
)
inlinestatic

Referenced by operator+().

static AutoDiffMatrix Opm::AutoDiffMatrix::addSI ( const AutoDiffMatrix lhs,
const AutoDiffMatrix rhs 
)
inlinestatic

Referenced by operator+().

static AutoDiffMatrix Opm::AutoDiffMatrix::addSS ( const AutoDiffMatrix lhs,
const AutoDiffMatrix rhs 
)
inlinestatic

Referenced by operator+().

double Opm::AutoDiffMatrix::coeff ( const int  row,
const int  col 
) const
inline

Returns element (row, col) in the matrix

int Opm::AutoDiffMatrix::cols ( ) const
inline

Returns number of columns in the matrix

Referenced by Opm::operator*(), and Opm::vertcatCollapseJacs().

static AutoDiffMatrix Opm::AutoDiffMatrix::createIdentity ( const int  num_rows_cols)
inlinestatic

Creates an identity matrix with num_rows_cols x num_rows_cols entries

References AutoDiffMatrix().

Referenced by Opm::AutoDiffBlock< double >::variable().

const SparseRep& Opm::AutoDiffMatrix::getSparse ( ) const
inline

Returns the sparse representation of this matrix. Note that this might be an expensive operation to perform if the internal structure is not sparse.

If we are not a sparse matrix, our internal variable sparse_ is undefined, and hence changing it so that it happens to be a sparse representation of our true data does not change our true data, and hence justifies that we do not really violate the const qualifier.

References toSparse().

static AutoDiffMatrix Opm::AutoDiffMatrix::mulDD ( const AutoDiffMatrix lhs,
const AutoDiffMatrix rhs 
)
inlinestatic

Referenced by operator*().

static AutoDiffMatrix Opm::AutoDiffMatrix::mulDS ( const AutoDiffMatrix lhs,
const AutoDiffMatrix rhs 
)
inlinestatic

References Opm::fastDiagSparseProduct().

Referenced by operator*().

static AutoDiffMatrix Opm::AutoDiffMatrix::mulSD ( const AutoDiffMatrix lhs,
const AutoDiffMatrix rhs 
)
inlinestatic

References Opm::fastSparseDiagProduct().

Referenced by operator*().

static AutoDiffMatrix Opm::AutoDiffMatrix::mulSS ( const AutoDiffMatrix lhs,
const AutoDiffMatrix rhs 
)
inlinestatic

References Opm::fastSparseProduct().

Referenced by operator*().

int Opm::AutoDiffMatrix::nonZeros ( ) const
inline

Returns number of non-zero elements in the matrix. Optimizes internally by exploiting that e.g., an n*n identity matrix has n non-zeros. Note that an n*n diagonal matrix is defined to have n non-zeros, even though several diagonal elements might be 0.0.

AutoDiffMatrix Opm::AutoDiffMatrix::operator* ( const AutoDiffMatrix rhs) const
inline

Multiplies two AutoDiffMatrices. Internally, this function optimizes the multiplication operation based on the structure of the matrix, e.g., multiplying M with a zero matrix will obviously yield a zero matrix.

References AutoDiffMatrix(), mulDD(), mulDS(), mulSD(), and mulSS().

AutoDiffMatrix Opm::AutoDiffMatrix::operator* ( const double  rhs) const
inline

Multiplies an AutoDiffMatrix with a scalar. Optimizes internally by exploiting that e.g., an identity matrix multiplied by a scalar x yields a diagonal matrix with x the diagonal.

Eigen::VectorXd Opm::AutoDiffMatrix::operator* ( const Eigen::VectorXd &  rhs) const
inline

Multiplies an AutoDiffMatrix with a vector. Optimizes internally by exploiting that e.g., an identity matrix multiplied by a vector yields the vector itself.

AutoDiffMatrix Opm::AutoDiffMatrix::operator+ ( const AutoDiffMatrix rhs) const
inline

Adds two AutoDiffMatrices. Internally, this function optimizes the addition operation based on the structure of the matrix, e.g., adding two zero matrices will obviously yield a zero matrix, and so on.

References addDD(), addDI(), addII(), addSD(), addSI(), and addSS().

AutoDiffMatrix& Opm::AutoDiffMatrix::operator+= ( const AutoDiffMatrix rhs)
inline
AutoDiffMatrix& Opm::AutoDiffMatrix::operator-= ( const AutoDiffMatrix rhs)
inline
AutoDiffMatrix Opm::AutoDiffMatrix::operator/ ( const double  rhs) const
inline

Divides an AutoDiffMatrix by a scalar. Optimizes internally by exploiting that e.g., an identity matrix divided by a scalar x yields a diagonal matrix with 1/x on the diagonal.

AutoDiffMatrix& Opm::AutoDiffMatrix::operator= ( const AutoDiffMatrix other)
default
AutoDiffMatrix& Opm::AutoDiffMatrix::operator= ( AutoDiffMatrix &&  other)
inline

References swap().

int Opm::AutoDiffMatrix::rows ( ) const
inline

Returns number of rows in the matrix

void Opm::AutoDiffMatrix::swap ( AutoDiffMatrix other)
inline

Referenced by AutoDiffMatrix(), and operator=().

template<class Scalar , int Options, class Index >
void Opm::AutoDiffMatrix::toSparse ( Eigen::SparseMatrix< Scalar, Options, Index > &  s) const
inline

Converts the AutoDiffMatrix to an Eigen SparseMatrix.This might be an expensive operation to perform for e.g., an identity matrix or a diagonal matrix.

Referenced by Opm::collapseJacs(), and getSparse().


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