Cantera  3.2.0a1
Loading...
Searching...
No Matches
MultiJac.cpp
Go to the documentation of this file.
1//! @file MultiJac.cpp Implementation file for class MultiJac
2
3// This file is part of Cantera. See License.txt in the top-level directory or
4// at https://cantera.org/license.txt for license and copyright information.
5
7#include <ctime>
8
9namespace Cantera
10{
11
12MultiJac::MultiJac(OneDim& r)
13{
14 warn_deprecated("MultiJac::MultiJac(OneDim&)",
15 "Non-default constructor to be removed after Cantera 3.2.");
16 m_resid = &r;
19}
20
22{
23 m_mat.bfill(0.0);
24 m_age = 10000;
25}
26
27void MultiJac::initialize(size_t nVars)
28{
29 m_dim = nVars;
31 m_ssdiag.resize(m_dim);
32 m_mask.resize(m_dim);
33}
34
36{
37 m_mat.resize(m_dim, bw, bw);
38}
39
40void MultiJac::setValue(size_t row, size_t col, double value)
41{
42 m_mat(row, col) = value;
43 if (row == col) {
44 m_ssdiag[row] = value;
45 }
46}
47
48void MultiJac::updateTransient(double rdt, integer* mask)
49{
50 for (size_t n = 0; n < m_dim; n++) {
51 m_mat.value(n,n) = m_ssdiag[n] - mask[n]*rdt;
52 }
53 factorize();
54}
55
56void MultiJac::eval(double* x0, double* resid0, double rdt)
57{
58 warn_deprecated("MultiJac::eval", "To be removed after Cantera 3.2. "
59 "Jacobian evaluation moved to OneDim::evalJacobian().");
60 if (!m_resid) {
61 throw CanteraError("MultiJac::eval", "Can only be used in combination with "
62 "(deprecated) constructor that takes a OneDim object.");
63 }
65}
66
67} // namespace
size_t nSubDiagonals() const
Number of subdiagonals.
void resize(size_t n, size_t kl, size_t ku, double v=0.0)
Resize the matrix problem.
void bfill(double v=0.0)
Fill or zero the matrix.
size_t nSuperDiagonals() const
Number of superdiagonals.
double & value(size_t i, size_t j)
Return a changeable reference to element (i,j).
Base class for exceptions thrown by Cantera classes.
BandMatrix m_mat
Underlying matrix storage.
Definition MultiJac.h:100
void setBandwidth(size_t bw) override
Used to provide system bandwidth for implementations that use banded matrix storage.
Definition MultiJac.cpp:35
void setValue(size_t row, size_t col, double value) override
Set a value at the specified row and column of the jacobian triplet vector.
Definition MultiJac.cpp:40
void factorize() override
Factorize the system matrix.
Definition MultiJac.h:59
vector< int > m_mask
Transient mask for transient terms, 1 if transient, 0 if steady-state.
Definition MultiJac.h:105
void eval(double *x0, double *resid0, double rdt)
Evaluates the Jacobian at x0 using finite differences.
Definition MultiJac.cpp:56
vector< double > m_ssdiag
Diagonal of the steady-state Jacobian.
Definition MultiJac.h:101
OneDim * m_resid
Residual evaluator for this Jacobian.
Definition MultiJac.h:98
void reset() override
Reset parameters as needed.
Definition MultiJac.cpp:21
void updateTransient(double rdt, integer *mask) override
Update the diagonal terms in the Jacobian by using the transient mask .
Definition MultiJac.cpp:48
void initialize(size_t nVars) override
Called during setup for any processes that need to be completed prior to setup functions used in sund...
Definition MultiJac.cpp:27
Container class for multiple-domain 1D problems.
Definition OneDim.h:28
size_t bandwidth() const
Jacobian bandwidth.
Definition OneDim.h:147
void evalJacobian(double *x0)
Evaluates the Jacobian at x0 using finite differences.
Definition OneDim.cpp:293
size_t m_dim
Dimension of the system.
int m_age
Age of the Jacobian (times incrementAge() has been called)
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
void warn_deprecated(const string &source, const AnyBase &node, const string &message)
A deprecation warning for syntax in an input file.
Definition AnyMap.cpp:1997