|
DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
|
Symmetric LDLT factorization of the cell-local Jacobian. More...
#include <EulerJacobian.hpp>
Public Types | |
| using | tLocalMat = ArrayEigenUniMatrixBatch< nVarsFixed, nVarsFixed > |
| Batch storage for nVars×nVars blocks. | |
| using | tComponent = Eigen::Matrix< real, nVarsFixed, nVarsFixed > |
| Single block matrix type. | |
| using | tVec = ArrayDOFV< nVarsFixed > |
| Distributed vector type. | |
| using | tBase = Direct::LocalLDLTBase< JacobianLocalLDLT< nVarsFixed >, Eigen::Matrix< real, nVarsFixed, nVarsFixed >, ArrayDOFV< nVarsFixed > > |
| CRTP base providing decompose/solve. | |
Public Member Functions | |
| JacobianLocalLDLT (const ssp< Direct::SerialSymLUStructure > &nMesh, int nVarsC) | |
| Construct and allocate the LDLT storage for the given adjacency structure. | |
| void | setZero () |
| Zero all D and L blocks and mark the factorization as not yet decomposed. | |
| auto | GetDiag (index i) |
Return the diagonal block for row i (compliant with LocalLDLTBase CRTP). | |
| auto | GetLower (index i, int iInLow) |
Return lower-triangular block iInLow for row i (compliant with LocalLDLTBase CRTP). | |
| void | PrintLog () |
| Print all non-zero entries (with diagonal inverted) to the log stream. | |
| tComponent | InvertDiag (const tComponent &v) |
| Invert a diagonal block using diagonal-preconditioned full-pivot LU. | |
Public Member Functions inherited from DNDS::Direct::LocalLDLTBase< JacobianLocalLDLT< nVarsFixed >, Eigen::Matrix< real, nVarsFixed, nVarsFixed >, ArrayDOFV< nVarsFixed > > | |
| LocalLDLTBase (ssp< Direct::SerialSymLUStructure > _symLU) | |
| virtual | ~LocalLDLTBase ()=default |
| void | GetDiag (index i) |
| void | GetLower (index i, int ij) |
| void | InvertDiag (const Eigen::Matrix< real, nVarsFixed, nVarsFixed > &v) |
| void | InPlaceDecompose () |
| void | MatMul (ArrayDOFV< nVarsFixed > &x, ArrayDOFV< nVarsFixed > &result) |
| void | Solve (ArrayDOFV< nVarsFixed > &b, ArrayDOFV< nVarsFixed > &result) |
Public Attributes | |
| tLocalMat | LDU |
| Sparse row storage for D and L blocks. | |
Public Attributes inherited from DNDS::Direct::LocalLDLTBase< JacobianLocalLDLT< nVarsFixed >, Eigen::Matrix< real, nVarsFixed, nVarsFixed >, ArrayDOFV< nVarsFixed > > | |
| ssp< Direct::SerialSymLUStructure > | symLU |
| bool | isDecomposed |
Symmetric LDLT factorization of the cell-local Jacobian.
Similar to JacobianLocalLU but exploits symmetry: only the diagonal and lower-triangular blocks are stored (no separate upper-triangular storage). Inherits the decompose/solve interface from Direct::LocalLDLTBase via CRTP.
Row layout in LDU for each cell iCell:
The upper triangle is implicitly L^T.
| nVarsFixed | Compile-time number of conservative variables (or Eigen::Dynamic). |
Definition at line 409 of file EulerJacobian.hpp.
| using DNDS::Euler::JacobianLocalLDLT< nVarsFixed >::tBase = Direct::LocalLDLTBase< JacobianLocalLDLT<nVarsFixed>, Eigen::Matrix<real, nVarsFixed, nVarsFixed>, ArrayDOFV<nVarsFixed> > |
CRTP base providing decompose/solve.
Definition at line 418 of file EulerJacobian.hpp.
| using DNDS::Euler::JacobianLocalLDLT< nVarsFixed >::tComponent = Eigen::Matrix<real, nVarsFixed, nVarsFixed> |
Single block matrix type.
Definition at line 416 of file EulerJacobian.hpp.
| using DNDS::Euler::JacobianLocalLDLT< nVarsFixed >::tLocalMat = ArrayEigenUniMatrixBatch<nVarsFixed, nVarsFixed> |
Batch storage for nVars×nVars blocks.
Definition at line 415 of file EulerJacobian.hpp.
| using DNDS::Euler::JacobianLocalLDLT< nVarsFixed >::tVec = ArrayDOFV<nVarsFixed> |
Distributed vector type.
Definition at line 417 of file EulerJacobian.hpp.
|
inline |
Construct and allocate the LDLT storage for the given adjacency structure.
Unlike JacobianLocalLU, only diagonal + lower blocks are allocated (no upper-triangular storage).
| nMesh | Shared pointer to the symmetric LU sparsity structure. |
| nVarsC | Runtime number of conservative variables. |
Definition at line 441 of file EulerJacobian.hpp.
|
inline |
Return the diagonal block for row i (compliant with LocalLDLTBase CRTP).
| i | Cell (row) index. |
Definition at line 469 of file EulerJacobian.hpp.
|
inline |
Return lower-triangular block iInLow for row i (compliant with LocalLDLTBase CRTP).
| i | Cell (row) index. |
| iInLow | Index within the lower-triangular neighbor list. |
Definition at line 476 of file EulerJacobian.hpp.
|
inline |
Invert a diagonal block using diagonal-preconditioned full-pivot LU.
Same algorithm as JacobianLocalLU::InvertDiag: preconditions with the diagonal inverse, then applies Eigen's fullPivLu().
| v | The nVars×nVars block matrix to invert. |
Definition at line 504 of file EulerJacobian.hpp.
|
inline |
Print all non-zero entries (with diagonal inverted) to the log stream.
Definition at line 482 of file EulerJacobian.hpp.
|
inline |
Zero all D and L blocks and mark the factorization as not yet decomposed.
Definition at line 456 of file EulerJacobian.hpp.
| tLocalMat DNDS::Euler::JacobianLocalLDLT< nVarsFixed >::LDU |
Sparse row storage for D and L blocks.
Each row i stores blocks in the order: [0] = diag, [1, 1 + symLU->lowerTriStructure[i].size()) = lower-triangular.
Definition at line 430 of file EulerJacobian.hpp.