Cantera  2.5.1
MultiJac.h
Go to the documentation of this file.
1 //! @file MultiJac.h
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 
6 #ifndef CT_MULTIJAC_H
7 #define CT_MULTIJAC_H
8 
10 #include "OneDim.h"
11 
12 namespace Cantera
13 {
14 
15 /**
16  * Class MultiJac evaluates the Jacobian of a system of equations defined by a
17  * residual function supplied by an instance of class OneDim. The residual
18  * function may consist of several linked 1D domains, with different variables
19  * in each domain.
20  * @ingroup onedim
21  */
22 class MultiJac : public BandMatrix
23 {
24 public:
25  MultiJac(OneDim& r);
26 
27  /**
28  * Evaluate the Jacobian at x0. The unperturbed residual function is resid0,
29  * which must be supplied on input. The third parameter 'rdt' is the
30  * reciprocal of the time step. If zero, the steady-state Jacobian is
31  * evaluated.
32  */
33  void eval(doublereal* x0, doublereal* resid0, double rdt);
34 
35  //! Elapsed CPU time spent computing the Jacobian.
36  doublereal elapsedTime() const {
37  return m_elapsed;
38  }
39 
40  //! Number of Jacobian evaluations.
41  int nEvals() const {
42  return m_nevals;
43  }
44 
45  //! Number of times 'incrementAge' has been called since the last evaluation
46  int age() const {
47  return m_age;
48  }
49 
50  //! Increment the Jacobian age.
51  void incrementAge() {
52  m_age++;
53  }
54 
55  void updateTransient(doublereal rdt, integer* mask);
56 
57  //! Set the Jacobian age.
58  void setAge(int age) {
59  m_age = age;
60  }
61 
62  vector_int& transientMask() {
63  return m_mask;
64  }
65 
66  void incrementDiagonal(int j, doublereal d);
67 
68 protected:
69  //! Residual evaluator for this Jacobian
70  /*!
71  * This is a pointer to the residual evaluator. This object isn't owned by
72  * this Jacobian object.
73  */
75 
76  vector_fp m_r1;
77  doublereal m_rtol, m_atol;
78  doublereal m_elapsed;
79  vector_fp m_ssdiag;
80  vector_int m_mask;
81  int m_nevals;
82  int m_age;
83  size_t m_size;
84  size_t m_points;
85 };
86 }
87 
88 #endif
Declarations for the class BandMatrix which is a child class of GeneralMatrix for banded matrices han...
A class for banded matrices, involving matrix inversion processes.
Definition: BandMatrix.h:35
Class MultiJac evaluates the Jacobian of a system of equations defined by a residual function supplie...
Definition: MultiJac.h:23
int nEvals() const
Number of Jacobian evaluations.
Definition: MultiJac.h:41
void setAge(int age)
Set the Jacobian age.
Definition: MultiJac.h:58
void eval(doublereal *x0, doublereal *resid0, double rdt)
Evaluate the Jacobian at x0.
Definition: MultiJac.cpp:43
void incrementAge()
Increment the Jacobian age.
Definition: MultiJac.h:51
doublereal elapsedTime() const
Elapsed CPU time spent computing the Jacobian.
Definition: MultiJac.h:36
int age() const
Number of times 'incrementAge' has been called since the last evaluation.
Definition: MultiJac.h:46
OneDim * m_resid
Residual evaluator for this Jacobian.
Definition: MultiJac.h:74
Container class for multiple-domain 1D problems.
Definition: OneDim.h:26
std::vector< int > vector_int
Vector of ints.
Definition: ct_defs.h:182
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
Definition: ct_defs.h:180
Namespace for the Cantera kernel.
Definition: AnyMap.cpp:264