|
DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
|
Per-cell diagonal or block-diagonal Jacobian storage for implicit time stepping. More...
#include <EulerJacobian.hpp>
Public Types | |
| using | TU = Eigen::Vector< real, nVarsFixed > |
| State vector type. | |
| using | tComponent = Eigen::Matrix< real, nVarsFixed, nVarsFixed > |
| Full block matrix type. | |
| using | tComponentDiag = Eigen::Vector< real, nVarsFixed > |
| Diagonal vector type. | |
Public Member Functions | |
| JacobianDiagBlock () | |
| Default constructor; mode and storage are uninitialized until SetModeAndInit(). | |
| void | SetModeAndInit (int mode, int nVarsC, ArrayDOFV< nVarsFixed > &mock) |
Select the storage mode and allocate arrays matching mock's layout. | |
| bool | isBlock () const |
| Return true if using matrix-block mode, false for scalar-diagonal. | |
| auto | getBlock (index iCell) |
Access the full nVars×nVars block matrix for cell iCell (block mode only). | |
| auto | getDiag (index iCell) |
Access the diagonal vector for cell iCell (scalar-diagonal mode only). | |
| tComponent | getValue (index iCell) const |
Return the Jacobian as a full matrix for cell iCell (either mode). | |
| index | Size () |
| Return the number of cells in the local partition. | |
| void | GetInvert () |
| Compute and cache the inverse of every cell's Jacobian block. | |
| template<class TV > | |
| TU | MatVecLeft (index iCell, TV v) |
Left-multiply a vector by the Jacobian block for cell iCell. | |
| template<class TV > | |
| TU | MatVecLeftInvert (index iCell, TV v) |
Left-multiply a vector by the inverse Jacobian block for cell iCell. | |
| void | clearValues () |
| Zero all stored values and invalidate the cached inverse. | |
Per-cell diagonal or block-diagonal Jacobian storage for implicit time stepping.
Operates in one of two modes selected by SetModeAndInit():
Both modes support lazy inversion via GetInvert(). The scalar mode inverts element-wise; the matrix mode uses diagonal-preconditioned full-pivot LU (InvertDiag pattern) for robustness on ill-conditioned blocks.
| nVarsFixed | Compile-time number of conservative variables (or Eigen::Dynamic). |
Definition at line 42 of file EulerJacobian.hpp.
| using DNDS::Euler::JacobianDiagBlock< nVarsFixed >::tComponent = Eigen::Matrix<real, nVarsFixed, nVarsFixed> |
Full block matrix type.
Definition at line 46 of file EulerJacobian.hpp.
| using DNDS::Euler::JacobianDiagBlock< nVarsFixed >::tComponentDiag = Eigen::Vector<real, nVarsFixed> |
Diagonal vector type.
Definition at line 47 of file EulerJacobian.hpp.
| using DNDS::Euler::JacobianDiagBlock< nVarsFixed >::TU = Eigen::Vector<real, nVarsFixed> |
State vector type.
Definition at line 45 of file EulerJacobian.hpp.
|
inline |
Default constructor; mode and storage are uninitialized until SetModeAndInit().
Definition at line 57 of file EulerJacobian.hpp.
|
inline |
Zero all stored values and invalidate the cached inverse.
Definition at line 224 of file EulerJacobian.hpp.
|
inline |
Access the full nVars×nVars block matrix for cell iCell (block mode only).
| iCell | Local cell index. |
Definition at line 101 of file EulerJacobian.hpp.
|
inline |
Access the diagonal vector for cell iCell (scalar-diagonal mode only).
| iCell | Local cell index. |
Definition at line 112 of file EulerJacobian.hpp.
|
inline |
Compute and cache the inverse of every cell's Jacobian block.
In matrix-block mode, uses diagonal-preconditioned full-pivot LU (the InvertDiag pattern) for robustness on ill-conditioned blocks. In scalar-diagonal mode, inverts element-wise.
This is a lazy operation: subsequent calls are no-ops until clearValues() invalidates the cache.
Definition at line 153 of file EulerJacobian.hpp.
|
inline |
Return the Jacobian as a full matrix for cell iCell (either mode).
In scalar-diagonal mode the returned matrix is constructed as asDiagonal().
| iCell | Local cell index. |
Definition at line 126 of file EulerJacobian.hpp.
|
inline |
Return true if using matrix-block mode, false for scalar-diagonal.
Definition at line 94 of file EulerJacobian.hpp.
|
inline |
Left-multiply a vector by the Jacobian block for cell iCell.
| TV | Vector type compatible with Eigen matrix-vector product. |
| iCell | Local cell index. |
| v | Input vector of length nVars. |
Definition at line 196 of file EulerJacobian.hpp.
|
inline |
Left-multiply a vector by the inverse Jacobian block for cell iCell.
GetInvert() must have been called first; otherwise behavior is undefined.
| TV | Vector type compatible with Eigen matrix-vector product. |
| iCell | Local cell index. |
| v | Input vector of length nVars. |
Definition at line 215 of file EulerJacobian.hpp.
|
inline |
Select the storage mode and allocate arrays matching mock's layout.
In block mode (mode != 0), allocates nVarsC×nVarsC matrices per cell. In scalar-diagonal mode (mode == 0), allocates nVarsC-length vectors. Ghost (son) arrays are allocated with zero size.
| mode | Storage mode: 0 = scalar-diagonal, nonzero = matrix-block. |
| nVarsC | Runtime number of variables (must equal nVarsFixed when fixed). |
| mock | A reference ArrayDOFV whose MPI layout and sizes are mimicked. |
Definition at line 70 of file EulerJacobian.hpp.
|
inline |
Return the number of cells in the local partition.
Definition at line 135 of file EulerJacobian.hpp.