Cantera
3.1.0a1
|
Methods for calculating analytical and/or numerical derivatives.
Classes | |
class | MultiJac |
Class MultiJac evaluates the Jacobian of a system of equations defined by a residual function supplied by an instance of class OneDim. More... | |
Functions | |
Eigen::SparseMatrix< double > | jacobian () |
Return semi-analytical Jacobian of an AdaptivePreconditioner object. More... | |
MultiJac & | jacobian () |
Return a reference to the Jacobian evaluator of an OneDim object. More... | |
virtual Eigen::SparseMatrix< double > | jacobian () |
Calculate the Jacobian of a specific Reactor specialization. More... | |
Kinetics derivatives are calculated with respect to temperature, pressure, molar concentrations and species mole fractions for forward/reverse/net rates of progress as well as creation/destruction and net production of species. The following suffixes are used to indicate derivatives:
Source term derivatives are based on a generic rate-of-progress expression for the \( i \)-th reaction \( R_i \), which is a function of temperature \( T \), pressure \( P \) and molar concentrations \( C_j \): \[ R_i = k_{f,i} C_M^{\nu_{M,i}} \prod_j C_j^{\nu_{ji}^\prime} - k_{r,i} C_M^{\nu_{M,i}} \prod_j C_j^{\nu_{ji}^{\prime\prime}} \] Forward/reverse rate expressions \( k_{f,i} \) and \( k_{r,i} \) are implemented by ReactionRate specializations; forward/reverse stoichiometric coefficients are \( \nu_{ji}^\prime \) and \( \nu_{ji}^{\prime\prime} \). Unless the reaction involves third-body colliders, \( \nu_{M,i} = 0 \). For three-body reactions, effective ThirdBody collider concentrations \( C_M \) are considered with \( \nu_{M,i} = 1 \). For more detailed information on relevant theory, see, for example, Perini, et al. [31] or Niemeyer, et al. [29], although specifics of Cantera's implementation may differ. Partial derivatives are obtained from the product rule, where resulting terms consider reaction rate derivatives, derivatives of the concentration product term, and, if applicable, third-body term derivatives. ReactionRate specializations may implement exact derivatives (example: ArrheniusRate::ddTScaledFromStruct) or approximate them numerically (examples: ReactionData::perturbTemperature, PlogData::perturbPressure, FalloffData::perturbThirdBodies). Derivatives of concentration and third-body terms are based on analytic expressions. Species creation and destruction rates are obtained by multiplying rate-of-progress vectors by stoichiometric coefficient matrices. As this is a linear operation, it is possible to calculate derivatives the same way. All derivatives are calculated for source terms while holding other properties constant, independent of whether equation of state or \( \sum X_k = 1 \) constraints are satisfied. Thus, derivatives deviate from Jacobians and numerical derivatives that implicitly enforce these constraints. Depending on application and equation of state, derivatives can nevertheless be used to obtain Jacobians, for example:
While some applications require exact derivatives, others can tolerate approximate derivatives that neglect terms to increase computational speed and/or improve Jacobian sparsity (example: AdaptivePreconditioner). Derivative evaluations settings are accessible by keyword/value pairs using the methods getDerivativeSettings() and setDerivativeSettings(). For BulkKinetics, the following keyword/value pairs are supported:
For InterfaceKinetics, the following keyword/value pairs are supported:
| |
virtual void | getDerivativeSettings (AnyMap &settings) const |
Retrieve derivative settings. More... | |
virtual void | setDerivativeSettings (const AnyMap &settings) |
Set/modify derivative settings. More... | |
virtual void | getFwdRateConstants_ddT (double *dkfwd) |
Calculate derivatives for forward rate constants with respect to temperature at constant pressure, molar concentration and mole fractions. More... | |
virtual void | getFwdRateConstants_ddP (double *dkfwd) |
Calculate derivatives for forward rate constants with respect to pressure at constant temperature, molar concentration and mole fractions. More... | |
virtual void | getFwdRateConstants_ddC (double *dkfwd) |
Calculate derivatives for forward rate constants with respect to molar concentration at constant temperature, pressure and mole fractions. More... | |
virtual void | getFwdRatesOfProgress_ddT (double *drop) |
Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions. More... | |
virtual void | getFwdRatesOfProgress_ddP (double *drop) |
Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions. More... | |
virtual void | getFwdRatesOfProgress_ddC (double *drop) |
Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions. More... | |
virtual Eigen::SparseMatrix< double > | fwdRatesOfProgress_ddX () |
Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration. More... | |
virtual Eigen::SparseMatrix< double > | fwdRatesOfProgress_ddCi () |
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant temperature, pressure and remaining species concentrations. More... | |
virtual void | getRevRatesOfProgress_ddT (double *drop) |
Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions. More... | |
virtual void | getRevRatesOfProgress_ddP (double *drop) |
Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions. More... | |
virtual void | getRevRatesOfProgress_ddC (double *drop) |
Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions. More... | |
virtual Eigen::SparseMatrix< double > | revRatesOfProgress_ddX () |
Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration. More... | |
virtual Eigen::SparseMatrix< double > | revRatesOfProgress_ddCi () |
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant temperature, pressure and remaining species concentrations. More... | |
virtual void | getNetRatesOfProgress_ddT (double *drop) |
Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions. More... | |
virtual void | getNetRatesOfProgress_ddP (double *drop) |
Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions. More... | |
virtual void | getNetRatesOfProgress_ddC (double *drop) |
Calculate derivatives for net rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions. More... | |
virtual Eigen::SparseMatrix< double > | netRatesOfProgress_ddX () |
Calculate derivatives for net rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration. More... | |
virtual Eigen::SparseMatrix< double > | netRatesOfProgress_ddCi () |
Calculate derivatives for net rates-of-progress with respect to species concentration at constant temperature, pressure, and remaining species concentrations. More... | |
void | getCreationRates_ddT (double *dwdot) |
Calculate derivatives for species creation rates with respect to temperature at constant pressure, molar concentration and mole fractions. More... | |
void | getCreationRates_ddP (double *dwdot) |
Calculate derivatives for species creation rates with respect to pressure at constant temperature, molar concentration and mole fractions. More... | |
void | getCreationRates_ddC (double *dwdot) |
Calculate derivatives for species creation rates with respect to molar concentration at constant temperature, pressure and mole fractions. More... | |
Eigen::SparseMatrix< double > | creationRates_ddX () |
Calculate derivatives for species creation rates with respect to species mole fractions at constant temperature, pressure and molar concentration. More... | |
Eigen::SparseMatrix< double > | creationRates_ddCi () |
Calculate derivatives for species creation rates with respect to species concentration at constant temperature, pressure, and concentration of all other species. More... | |
void | getDestructionRates_ddT (double *dwdot) |
Calculate derivatives for species destruction rates with respect to temperature at constant pressure, molar concentration and mole fractions. More... | |
void | getDestructionRates_ddP (double *dwdot) |
Calculate derivatives for species destruction rates with respect to pressure at constant temperature, molar concentration and mole fractions. More... | |
void | getDestructionRates_ddC (double *dwdot) |
Calculate derivatives for species destruction rates with respect to molar concentration at constant temperature, pressure and mole fractions. More... | |
Eigen::SparseMatrix< double > | destructionRates_ddX () |
Calculate derivatives for species destruction rates with respect to species mole fractions at constant temperature, pressure and molar concentration. More... | |
Eigen::SparseMatrix< double > | destructionRates_ddCi () |
Calculate derivatives for species destruction rates with respect to species concentration at constant temperature, pressure, and concentration of all other species. More... | |
void | getNetProductionRates_ddT (double *dwdot) |
Calculate derivatives for species net production rates with respect to temperature at constant pressure, molar concentration and mole fractions. More... | |
void | getNetProductionRates_ddP (double *dwdot) |
Calculate derivatives for species net production rates with respect to pressure at constant temperature, molar concentration and mole fractions. More... | |
void | getNetProductionRates_ddC (double *dwdot) |
Calculate derivatives for species net production rates with respect to molar concentration at constant temperature, pressure and mole fractions. More... | |
Eigen::SparseMatrix< double > | netProductionRates_ddX () |
Calculate derivatives for species net production rates with respect to species mole fractions at constant temperature, pressure and molar concentration. More... | |
Eigen::SparseMatrix< double > | netProductionRates_ddCi () |
Calculate derivatives for species net production rates with respect to species concentration at constant temperature, pressure, and concentration of all other species. More... | |
|
inlinevirtual |
Retrieve derivative settings.
settings | AnyMap containing settings determining derivative evaluation. |
Reimplemented in InterfaceKinetics, and BulkKinetics.
Definition at line 692 of file Kinetics.h.
|
inlinevirtual |
Set/modify derivative settings.
settings | AnyMap containing settings determining derivative evaluation. |
Reimplemented in InterfaceKinetics, and BulkKinetics.
Definition at line 703 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for forward rate constants with respect to temperature at constant pressure, molar concentration and mole fractions.
[out] | dkfwd | Output vector of derivatives. Length: nReactions(). |
Reimplemented in BulkKinetics.
Definition at line 714 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for forward rate constants with respect to pressure at constant temperature, molar concentration and mole fractions.
[out] | dkfwd | Output vector of derivatives. Length: nReactions(). |
Reimplemented in BulkKinetics.
Definition at line 725 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for forward rate constants with respect to molar concentration at constant temperature, pressure and mole fractions.
[out] | dkfwd | Output vector of derivatives. Length: nReactions(). |
Reimplemented in BulkKinetics.
Definition at line 739 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions.
[out] | drop | Output vector of derivatives. Length: nReactions(). |
Reimplemented in BulkKinetics.
Definition at line 750 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions.
[out] | drop | Output vector of derivatives. Length: nReactions(). |
Reimplemented in BulkKinetics.
Definition at line 761 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions.
[out] | drop | Output vector of derivatives. Length: nReactions(). |
Reimplemented in BulkKinetics.
Definition at line 775 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration.
The method returns a matrix with nReactions() rows and nTotalSpecies() columns. For a derivative with respect to \( X_i \), all other \( X_j \) are held constant, rather than enforcing \( \sum X_j = 1 \).
Reimplemented in BulkKinetics.
Definition at line 792 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant temperature, pressure and remaining species concentrations.
The method returns a matrix with nReactions() rows and nTotalSpecies() columns. For a derivative with respect to \( c_i \), all other \( c_j \) are held constant.
Reimplemented in InterfaceKinetics, and BulkKinetics.
Definition at line 812 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions.
[out] | drop | Output vector of derivatives. Length: nReactions(). |
Reimplemented in BulkKinetics.
Definition at line 823 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions.
[out] | drop | Output vector of derivatives. Length: nReactions(). |
Reimplemented in BulkKinetics.
Definition at line 834 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions.
[out] | drop | Output vector of derivatives. Length: nReactions(). |
Reimplemented in BulkKinetics.
Definition at line 848 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration.
The method returns a matrix with nReactions() rows and nTotalSpecies() columns. For a derivative with respect to \( X_i \), all other \( X_j \) are held constant, rather than enforcing \( \sum X_j = 1 \).
Reimplemented in BulkKinetics.
Definition at line 865 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant temperature, pressure and remaining species concentrations.
The method returns a matrix with nReactions() rows and nTotalSpecies() columns. For a derivative with respect to \( c_i \), all other \( c_j \) are held constant.
Reimplemented in InterfaceKinetics, and BulkKinetics.
Definition at line 885 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions.
[out] | drop | Output vector of derivatives. Length: nReactions(). |
Reimplemented in BulkKinetics.
Definition at line 896 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions.
[out] | drop | Output vector of derivatives. Length: nReactions(). |
Reimplemented in BulkKinetics.
Definition at line 907 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for net rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions.
[out] | drop | Output vector of derivatives. Length: nReactions(). |
Reimplemented in BulkKinetics.
Definition at line 921 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for net rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration.
The method returns a matrix with nReactions() rows and nTotalSpecies() columns. For a derivative with respect to \( X_i \), all other \( X_j \) are held constant, rather than enforcing \( \sum X_j = 1 \).
Reimplemented in BulkKinetics.
Definition at line 938 of file Kinetics.h.
|
inlinevirtual |
Calculate derivatives for net rates-of-progress with respect to species concentration at constant temperature, pressure, and remaining species concentrations.
The method returns a matrix with nReactions() rows and nTotalSpecies() columns. For a derivative with respect to \( c_i \), all other \( c_j \) are held constant.
Reimplemented in InterfaceKinetics, and BulkKinetics.
Definition at line 958 of file Kinetics.h.
void getCreationRates_ddT | ( | double * | dwdot | ) |
Calculate derivatives for species creation rates with respect to temperature at constant pressure, molar concentration and mole fractions.
[out] | dwdot | Output vector of derivatives. Length: m_kk. |
Definition at line 374 of file Kinetics.cpp.
void getCreationRates_ddP | ( | double * | dwdot | ) |
Calculate derivatives for species creation rates with respect to pressure at constant temperature, molar concentration and mole fractions.
[out] | dwdot | Output vector of derivatives. Length: m_kk. |
Definition at line 386 of file Kinetics.cpp.
void getCreationRates_ddC | ( | double * | dwdot | ) |
Calculate derivatives for species creation rates with respect to molar concentration at constant temperature, pressure and mole fractions.
[out] | dwdot | Output vector of derivatives. Length: m_kk. |
Definition at line 398 of file Kinetics.cpp.
Eigen::SparseMatrix< double > creationRates_ddX | ( | ) |
Calculate derivatives for species creation rates with respect to species mole fractions at constant temperature, pressure and molar concentration.
The method returns a square matrix with nTotalSpecies() rows and columns. For a derivative with respect to \( X_i \), all other \( X_j \) are held constant, rather than enforcing \( \sum X_j = 1 \).
Definition at line 410 of file Kinetics.cpp.
Eigen::SparseMatrix< double > creationRates_ddCi | ( | ) |
Calculate derivatives for species creation rates with respect to species concentration at constant temperature, pressure, and concentration of all other species.
The method returns a square matrix with nTotalSpecies() rows and columns. For a derivative with respect to \( c_i \), all other \( c_j \) are held constant.
Definition at line 420 of file Kinetics.cpp.
void getDestructionRates_ddT | ( | double * | dwdot | ) |
Calculate derivatives for species destruction rates with respect to temperature at constant pressure, molar concentration and mole fractions.
[out] | dwdot | Output vector of derivatives. Length: m_kk. |
Definition at line 430 of file Kinetics.cpp.
void getDestructionRates_ddP | ( | double * | dwdot | ) |
Calculate derivatives for species destruction rates with respect to pressure at constant temperature, molar concentration and mole fractions.
[out] | dwdot | Output vector of derivatives. Length: m_kk. |
Definition at line 442 of file Kinetics.cpp.
void getDestructionRates_ddC | ( | double * | dwdot | ) |
Calculate derivatives for species destruction rates with respect to molar concentration at constant temperature, pressure and mole fractions.
[out] | dwdot | Output vector of derivatives. Length: m_kk. |
Definition at line 454 of file Kinetics.cpp.
Eigen::SparseMatrix< double > destructionRates_ddX | ( | ) |
Calculate derivatives for species destruction rates with respect to species mole fractions at constant temperature, pressure and molar concentration.
The method returns a square matrix with nTotalSpecies() rows and columns. For a derivative with respect to \( X_i \), all other \( X_j \) are held constant, rather than enforcing \( \sum X_j = 1 \).
Definition at line 466 of file Kinetics.cpp.
Eigen::SparseMatrix< double > destructionRates_ddCi | ( | ) |
Calculate derivatives for species destruction rates with respect to species concentration at constant temperature, pressure, and concentration of all other species.
The method returns a square matrix with nTotalSpecies() rows and columns. For a derivative with respect to \( c_i \), all other \( c_j \) are held constant.
Definition at line 476 of file Kinetics.cpp.
void getNetProductionRates_ddT | ( | double * | dwdot | ) |
Calculate derivatives for species net production rates with respect to temperature at constant pressure, molar concentration and mole fractions.
[out] | dwdot | Output vector of derivatives. Length: m_kk. |
Definition at line 486 of file Kinetics.cpp.
void getNetProductionRates_ddP | ( | double * | dwdot | ) |
Calculate derivatives for species net production rates with respect to pressure at constant temperature, molar concentration and mole fractions.
[out] | dwdot | Output vector of derivatives. Length: m_kk. |
Definition at line 494 of file Kinetics.cpp.
void getNetProductionRates_ddC | ( | double * | dwdot | ) |
Calculate derivatives for species net production rates with respect to molar concentration at constant temperature, pressure and mole fractions.
[out] | dwdot | Output vector of derivatives. Length: m_kk. |
Definition at line 502 of file Kinetics.cpp.
Eigen::SparseMatrix< double > netProductionRates_ddX | ( | ) |
Calculate derivatives for species net production rates with respect to species mole fractions at constant temperature, pressure and molar concentration.
The method returns a square matrix with nTotalSpecies() rows and columns. For a derivative with respect to \( X_i \), all other \( X_j \) are held constant, rather than enforcing \( \sum X_j = 1 \).
Definition at line 510 of file Kinetics.cpp.
Eigen::SparseMatrix< double > netProductionRates_ddCi | ( | ) |
Calculate derivatives for species net production rates with respect to species concentration at constant temperature, pressure, and concentration of all other species.
The method returns a square matrix with nTotalSpecies() rows and columns. For a derivative with respect to \( c_i \), all other \( c_j \) are held constant.
Definition at line 515 of file Kinetics.cpp.
|
inline |
Return semi-analytical Jacobian of an AdaptivePreconditioner object.
Definition at line 53 of file AdaptivePreconditioner.h.
MultiJac & jacobian | ( | ) |
Return a reference to the Jacobian evaluator of an OneDim object.
Definition at line 87 of file OneDim.cpp.
|
inlinevirtual |
Calculate the Jacobian of a specific Reactor specialization.
Reimplemented in IdealGasMoleReactor, and IdealGasConstPressureMoleReactor.