Cantera  3.2.0
Loading...
Searching...
No Matches
GeneralMatrix.h
Go to the documentation of this file.
1/**
2 * @file GeneralMatrix.h
3 * Declarations for the class GeneralMatrix which is a virtual base class for matrices handled by solvers
4 * (see class @ref matrices and @link Cantera::GeneralMatrix GeneralMatrix@endlink).
5 */
6
7// This file is part of Cantera. See License.txt in the top-level directory or
8// at https://cantera.org/license.txt for license and copyright information.
9
10#ifndef CT_GENERALMATRIX_H
11#define CT_GENERALMATRIX_H
12
15#include "cantera/base/global.h"
16
17namespace Cantera
18{
19
20//! Generic matrix
21//! @ingroup matrices
23{
24public:
25 //! Base Constructor
26 GeneralMatrix() = default;
27
28 virtual ~GeneralMatrix() = default;
29
30 //! Zero the matrix elements
31 virtual void zero() = 0;
32
33 //! Multiply A*b and write result to prod.
34 /*!
35 * @param b Vector to do the rh multiplication
36 * @param prod OUTPUT vector to receive the result
37 */
38 virtual void mult(const double* b, double* prod) const = 0;
39
40 //! Multiply b*A and write result to prod.
41 /*!
42 * @param b Vector to do the lh multiplication
43 * @param prod OUTPUT vector to receive the result
44 */
45 virtual void leftMult(const double* const b, double* const prod) const = 0;
46
47 //! Factors the A matrix, overwriting A.
48 /*!
49 * We flip m_factored boolean to indicate that the matrix is now A-1.
50 */
51 virtual int factor() = 0;
52
53 //! Factors the A matrix using the QR algorithm, overwriting A
54 /*!
55 * we set m_factored to 2 to indicate the matrix is now QR factored
56 *
57 * @returns the info variable from LAPACK
58 */
59 virtual int factorQR() {
60 throw NotImplementedError("GeneralMatrix::factorQR");
61 }
62
63 //! Returns an estimate of the inverse of the condition number for the matrix
64 /*!
65 * The matrix must have been previously factored using the QR algorithm
66 *
67 * @returns the inverse of the condition number
68 */
69 virtual double rcondQR() {
70 throw NotImplementedError("GeneralMatrix::rcondQR");
71 }
72
73 //! Returns an estimate of the inverse of the condition number for the matrix
74 /*!
75 * The matrix must have been previously factored using the LU algorithm
76 *
77 * @param a1norm Norm of the matrix
78 * @returns the inverse of the condition number
79 */
80 virtual double rcond(double a1norm) = 0;
81
82 //! Change the way the matrix is factored
83 /*!
84 * @param fAlgorithm integer
85 * 0 LU factorization
86 * 1 QR factorization
87 */
88 virtual void useFactorAlgorithm(int fAlgorithm) {
89 throw NotImplementedError("GeneralMatrix::useFactorAlgorithm");
90 };
91
92 //! Return the factor algorithm used
93 virtual int factorAlgorithm() const = 0;
94
95 //! Calculate the one norm of the matrix
96 virtual double oneNorm() const = 0;
97
98 //! Return the number of rows in the matrix
99 virtual size_t nRows() const = 0;
100
101 //! clear the factored flag
102 virtual void clearFactorFlag() {
103 m_factored = 0;
104 };
105
106 //! Solves the Ax = b system returning x in the b spot.
107 /*!
108 * @param b Vector for the RHS of the equation system
109 * @param nrhs Number of right-hand sides to solve, default 1
110 * @param ldb Leading dimension of the right-hand side array. Defaults to
111 * nRows()
112 * @deprecated After %Cantera 3.2, the return type will be `void`
113 */
114 virtual int solve(double* b, size_t nrhs=1, size_t ldb=0) = 0;
115
116 //! true if the current factorization is up to date with the matrix
117 virtual bool factored() const {
118 return (m_factored != 0);
119 }
120
121 //! Return a pointer to the top of column j, columns are assumed to be
122 //! contiguous in memory
123 /*!
124 * @param j Value of the column
125 * @returns a pointer to the top of the column
126 */
127 virtual double* ptrColumn(size_t j) = 0;
128
129 //! Index into the (i,j) element
130 /*!
131 * @param i row
132 * @param j column
133 * @returns a changeable reference to the matrix entry
134 */
135 virtual double& operator()(size_t i, size_t j) = 0;
136
137 //! Constant Index into the (i,j) element
138 /*!
139 * @param i row
140 * @param j column
141 * @returns an unchangeable reference to the matrix entry
142 */
143 virtual double operator()(size_t i, size_t j) const = 0;
144
145 //! Return a vector of const pointers to the columns
146 /*!
147 * Note the value of the pointers are protected by their being const.
148 * However, the value of the matrix is open to being changed.
149 *
150 * @returns a vector of pointers to the top of the columns of the matrices.
151 */
152 virtual double* const* colPts() = 0;
153
154 //! Check to see if we have any zero rows in the Jacobian
155 /*!
156 * This utility routine checks to see if any rows are zero. The smallest row
157 * is returned along with the largest coefficient in that row
158 *
159 * @param valueSmall OUTPUT value of the largest coefficient in the smallest row
160 * @return index of the row that is most nearly zero
161 */
162 virtual size_t checkRows(double& valueSmall) const = 0;
163
164 //! Check to see if we have any zero columns in the Jacobian
165 /*!
166 * This utility routine checks to see if any columns are zero. The smallest
167 * column is returned along with the largest coefficient in that column
168 *
169 * @param valueSmall OUTPUT value of the largest coefficient in the smallest column
170 * @return index of the column that is most nearly zero
171 */
172 virtual size_t checkColumns(double& valueSmall) const = 0;
173
174protected:
175 //! Indicates whether the matrix is factored. 0 for unfactored; Non-zero
176 //! values indicate a particular factorization (LU=1, QR=2).
177 int m_factored = false;
178};
179
180}
181#endif
virtual double rcondQR()
Returns an estimate of the inverse of the condition number for the matrix.
GeneralMatrix()=default
Base Constructor.
virtual int factor()=0
Factors the A matrix, overwriting A.
virtual size_t checkColumns(double &valueSmall) const =0
Check to see if we have any zero columns in the Jacobian.
virtual void clearFactorFlag()
clear the factored flag
virtual size_t checkRows(double &valueSmall) const =0
Check to see if we have any zero rows in the Jacobian.
virtual double rcond(double a1norm)=0
Returns an estimate of the inverse of the condition number for the matrix.
virtual int factorQR()
Factors the A matrix using the QR algorithm, overwriting A.
virtual int factorAlgorithm() const =0
Return the factor algorithm used.
virtual int solve(double *b, size_t nrhs=1, size_t ldb=0)=0
Solves the Ax = b system returning x in the b spot.
virtual double oneNorm() const =0
Calculate the one norm of the matrix.
virtual size_t nRows() const =0
Return the number of rows in the matrix.
virtual double & operator()(size_t i, size_t j)=0
Index into the (i,j) element.
virtual double operator()(size_t i, size_t j) const =0
Constant Index into the (i,j) element.
virtual double *const * colPts()=0
Return a vector of const pointers to the columns.
virtual void zero()=0
Zero the matrix elements.
int m_factored
Indicates whether the matrix is factored.
virtual void useFactorAlgorithm(int fAlgorithm)
Change the way the matrix is factored.
virtual double * ptrColumn(size_t j)=0
Return a pointer to the top of column j, columns are assumed to be contiguous in memory.
virtual void leftMult(const double *const b, double *const prod) const =0
Multiply b*A and write result to prod.
virtual void mult(const double *b, double *prod) const =0
Multiply A*b and write result to prod.
virtual bool factored() const
true if the current factorization is up to date with the matrix
An error indicating that an unimplemented function has been called.
This file contains definitions of constants, types and terms that are used in internal routines and a...
Definitions for the classes that are thrown when Cantera experiences an error condition (also contain...
This file contains definitions for utility functions and text for modules, inputfiles and logging,...
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595