DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
DNDS::Euler::JacobianDiagBlock< nVarsFixed > Class Template Reference

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.
 

Detailed Description

template<int nVarsFixed>
class DNDS::Euler::JacobianDiagBlock< nVarsFixed >

Per-cell diagonal or block-diagonal Jacobian storage for implicit time stepping.

Operates in one of two modes selected by SetModeAndInit():

  • Scalar-diagonal (mode 0): stores an nVars-length diagonal vector per cell. Matrix-vector products are element-wise multiplications.
  • Matrix-block (mode 1): stores a full nVars×nVars matrix per cell. Matrix-vector products are dense matrix multiplies.

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.

Template Parameters
nVarsFixedCompile-time number of conservative variables (or Eigen::Dynamic).

Definition at line 42 of file EulerJacobian.hpp.

Member Typedef Documentation

◆ tComponent

template<int nVarsFixed>
using DNDS::Euler::JacobianDiagBlock< nVarsFixed >::tComponent = Eigen::Matrix<real, nVarsFixed, nVarsFixed>

Full block matrix type.

Definition at line 46 of file EulerJacobian.hpp.

◆ tComponentDiag

template<int nVarsFixed>
using DNDS::Euler::JacobianDiagBlock< nVarsFixed >::tComponentDiag = Eigen::Vector<real, nVarsFixed>

Diagonal vector type.

Definition at line 47 of file EulerJacobian.hpp.

◆ TU

template<int nVarsFixed>
using DNDS::Euler::JacobianDiagBlock< nVarsFixed >::TU = Eigen::Vector<real, nVarsFixed>

State vector type.

Definition at line 45 of file EulerJacobian.hpp.

Constructor & Destructor Documentation

◆ JacobianDiagBlock()

template<int nVarsFixed>
DNDS::Euler::JacobianDiagBlock< nVarsFixed >::JacobianDiagBlock ( )
inline

Default constructor; mode and storage are uninitialized until SetModeAndInit().

Definition at line 57 of file EulerJacobian.hpp.

Member Function Documentation

◆ clearValues()

template<int nVarsFixed>
void DNDS::Euler::JacobianDiagBlock< nVarsFixed >::clearValues ( )
inline

Zero all stored values and invalidate the cached inverse.

Definition at line 224 of file EulerJacobian.hpp.

Here is the call graph for this function:

◆ getBlock()

template<int nVarsFixed>
auto DNDS::Euler::JacobianDiagBlock< nVarsFixed >::getBlock ( index  iCell)
inline

Access the full nVars×nVars block matrix for cell iCell (block mode only).

Parameters
iCellLocal cell index.
Returns
Eigen map/ref to the stored block matrix.

Definition at line 101 of file EulerJacobian.hpp.

Here is the call graph for this function:

◆ getDiag()

template<int nVarsFixed>
auto DNDS::Euler::JacobianDiagBlock< nVarsFixed >::getDiag ( index  iCell)
inline

Access the diagonal vector for cell iCell (scalar-diagonal mode only).

Parameters
iCellLocal cell index.
Returns
Eigen map/ref to the stored diagonal vector.

Definition at line 112 of file EulerJacobian.hpp.

Here is the call graph for this function:

◆ GetInvert()

template<int nVarsFixed>
void DNDS::Euler::JacobianDiagBlock< nVarsFixed >::GetInvert ( )
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.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ getValue()

template<int nVarsFixed>
tComponent DNDS::Euler::JacobianDiagBlock< nVarsFixed >::getValue ( index  iCell) const
inline

Return the Jacobian as a full matrix for cell iCell (either mode).

In scalar-diagonal mode the returned matrix is constructed as asDiagonal().

Parameters
iCellLocal cell index.
Returns
nVars×nVars matrix (by value).

Definition at line 126 of file EulerJacobian.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ isBlock()

template<int nVarsFixed>
bool DNDS::Euler::JacobianDiagBlock< nVarsFixed >::isBlock ( ) const
inline

Return true if using matrix-block mode, false for scalar-diagonal.

Definition at line 94 of file EulerJacobian.hpp.

Here is the caller graph for this function:

◆ MatVecLeft()

template<int nVarsFixed>
template<class TV >
TU DNDS::Euler::JacobianDiagBlock< nVarsFixed >::MatVecLeft ( index  iCell,
TV  v 
)
inline

Left-multiply a vector by the Jacobian block for cell iCell.

Template Parameters
TVVector type compatible with Eigen matrix-vector product.
Parameters
iCellLocal cell index.
vInput vector of length nVars.
Returns
Result of J * v.

Definition at line 196 of file EulerJacobian.hpp.

Here is the call graph for this function:

◆ MatVecLeftInvert()

template<int nVarsFixed>
template<class TV >
TU DNDS::Euler::JacobianDiagBlock< nVarsFixed >::MatVecLeftInvert ( index  iCell,
TV  v 
)
inline

Left-multiply a vector by the inverse Jacobian block for cell iCell.

GetInvert() must have been called first; otherwise behavior is undefined.

Template Parameters
TVVector type compatible with Eigen matrix-vector product.
Parameters
iCellLocal cell index.
vInput vector of length nVars.
Returns
Result of J^{-1} * v.

Definition at line 215 of file EulerJacobian.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetModeAndInit()

template<int nVarsFixed>
void DNDS::Euler::JacobianDiagBlock< nVarsFixed >::SetModeAndInit ( int  mode,
int  nVarsC,
ArrayDOFV< nVarsFixed > &  mock 
)
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.

Parameters
modeStorage mode: 0 = scalar-diagonal, nonzero = matrix-block.
nVarsCRuntime number of variables (must equal nVarsFixed when fixed).
mockA reference ArrayDOFV whose MPI layout and sizes are mimicked.

Definition at line 70 of file EulerJacobian.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Size()

template<int nVarsFixed>
index DNDS::Euler::JacobianDiagBlock< nVarsFixed >::Size ( )
inline

Return the number of cells in the local partition.

Definition at line 135 of file EulerJacobian.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

The documentation for this class was generated from the following file: