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

Complete LU factorization of the cell-local Jacobian for GMRES preconditioning. More...

#include <EulerJacobian.hpp>

Inheritance diagram for DNDS::Euler::JacobianLocalLU< nVarsFixed >:
[legend]
Collaboration diagram for DNDS::Euler::JacobianLocalLU< nVarsFixed >:
[legend]

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::LocalLUBase< JacobianLocalLU< nVarsFixed >, Eigen::Matrix< real, nVarsFixed, nVarsFixed >, ArrayDOFV< nVarsFixed > >
 CRTP base providing decompose/solve.
 

Public Member Functions

 JacobianLocalLU (const ssp< Direct::SerialSymLUStructure > &nMesh, int nVarsC)
 Construct and allocate the LU storage for the given adjacency structure.
 
void setZero ()
 Zero all D, L, U blocks and mark the factorization as not yet decomposed.
 
auto GetDiag (index i)
 Return the diagonal block for row i (compliant with LocalLUBase CRTP).
 
auto GetLower (index i, int iInLow)
 Return lower-triangular block iInLow for row i (compliant with LocalLUBase CRTP).
 
auto GetUpper (index i, int iInUpp)
 Return upper-triangular block iInUpp for row i (compliant with LocalLUBase 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::LocalLUBase< JacobianLocalLU< nVarsFixed >, Eigen::Matrix< real, nVarsFixed, nVarsFixed >, ArrayDOFV< nVarsFixed > >
 LocalLUBase (ssp< Direct::SerialSymLUStructure > _symLU)
 
virtual ~LocalLUBase ()=default
 
void GetDiag (index i)
 
void GetLower (index i, int ij)
 
void GetUpper (index i, int ij)
 
void InvertDiag (const Eigen::Matrix< real, nVarsFixed, nVarsFixed > &v)
 
void InPlaceDecompose ()
 
void InPlaceDecomposeV1 ()
 
void InPlaceDecomposeV2 ()
 
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, L, U blocks.
 
- Public Attributes inherited from DNDS::Direct::LocalLUBase< JacobianLocalLU< nVarsFixed >, Eigen::Matrix< real, nVarsFixed, nVarsFixed >, ArrayDOFV< nVarsFixed > >
ssp< Direct::SerialSymLUStructuresymLU
 
bool isDecomposed
 

Detailed Description

template<int nVarsFixed>
struct DNDS::Euler::JacobianLocalLU< nVarsFixed >

Complete LU factorization of the cell-local Jacobian for GMRES preconditioning.

Stores diagonal (D), lower-triangular (L), and upper-triangular (U) blocks indexed by the mesh adjacency graph (SerialSymLUStructure). Inherits the decompose/solve interface from Direct::LocalLUBase via CRTP.

Row layout in LDU for each cell iCell:

  • Index 0: diagonal block D.
  • Indices [1, 1 + nLower): lower-triangular off-diagonal blocks L.
  • Indices [1 + nLower, end): upper-triangular off-diagonal blocks U.

where nLower = symLU->lowerTriStructure[iCell].size().

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

Definition at line 268 of file EulerJacobian.hpp.

Member Typedef Documentation

◆ tBase

template<int nVarsFixed>
using DNDS::Euler::JacobianLocalLU< nVarsFixed >::tBase = Direct::LocalLUBase< JacobianLocalLU<nVarsFixed>, Eigen::Matrix<real, nVarsFixed, nVarsFixed>, ArrayDOFV<nVarsFixed> >

CRTP base providing decompose/solve.

Definition at line 277 of file EulerJacobian.hpp.

◆ tComponent

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

Single block matrix type.

Definition at line 275 of file EulerJacobian.hpp.

◆ tLocalMat

template<int nVarsFixed>
using DNDS::Euler::JacobianLocalLU< nVarsFixed >::tLocalMat = ArrayEigenUniMatrixBatch<nVarsFixed, nVarsFixed>

Batch storage for nVars×nVars blocks.

Definition at line 274 of file EulerJacobian.hpp.

◆ tVec

template<int nVarsFixed>
using DNDS::Euler::JacobianLocalLU< nVarsFixed >::tVec = ArrayDOFV<nVarsFixed>

Distributed vector type.

Definition at line 276 of file EulerJacobian.hpp.

Constructor & Destructor Documentation

◆ JacobianLocalLU()

template<int nVarsFixed>
DNDS::Euler::JacobianLocalLU< nVarsFixed >::JacobianLocalLU ( const ssp< Direct::SerialSymLUStructure > &  nMesh,
int  nVarsC 
)
inline

Construct and allocate the LU storage for the given adjacency structure.

Parameters
nMeshShared pointer to the symmetric LU sparsity structure.
nVarsCRuntime number of conservative variables.

Definition at line 297 of file EulerJacobian.hpp.

Here is the call graph for this function:

Member Function Documentation

◆ GetDiag()

template<int nVarsFixed>
auto DNDS::Euler::JacobianLocalLU< nVarsFixed >::GetDiag ( index  i)
inline

Return the diagonal block for row i (compliant with LocalLUBase CRTP).

Parameters
iCell (row) index.

Definition at line 326 of file EulerJacobian.hpp.

◆ GetLower()

template<int nVarsFixed>
auto DNDS::Euler::JacobianLocalLU< nVarsFixed >::GetLower ( index  i,
int  iInLow 
)
inline

Return lower-triangular block iInLow for row i (compliant with LocalLUBase CRTP).

Parameters
iCell (row) index.
iInLowIndex within the lower-triangular neighbor list.

Definition at line 333 of file EulerJacobian.hpp.

◆ GetUpper()

template<int nVarsFixed>
auto DNDS::Euler::JacobianLocalLU< nVarsFixed >::GetUpper ( index  i,
int  iInUpp 
)
inline

Return upper-triangular block iInUpp for row i (compliant with LocalLUBase CRTP).

Parameters
iCell (row) index.
iInUppIndex within the upper-triangular neighbor list.

Definition at line 340 of file EulerJacobian.hpp.

◆ InvertDiag()

template<int nVarsFixed>
tComponent DNDS::Euler::JacobianLocalLU< nVarsFixed >::InvertDiag ( const tComponent v)
inline

Invert a diagonal block using diagonal-preconditioned full-pivot LU.

Preconditions the matrix by scaling rows with the inverse of the diagonal, then applies Eigen's fullPivLu(). This two-step approach improves numerical stability for ill-conditioned Jacobian blocks.

Parameters
vThe nVars×nVars block matrix to invert.
Returns
The inverse matrix v^{-1}.

Definition at line 369 of file EulerJacobian.hpp.

◆ PrintLog()

template<int nVarsFixed>
void DNDS::Euler::JacobianLocalLU< nVarsFixed >::PrintLog ( )
inline

Print all non-zero entries (with diagonal inverted) to the log stream.

Definition at line 346 of file EulerJacobian.hpp.

Here is the call graph for this function:

◆ setZero()

template<int nVarsFixed>
void DNDS::Euler::JacobianLocalLU< nVarsFixed >::setZero ( )
inline

Zero all D, L, U blocks and mark the factorization as not yet decomposed.

Definition at line 313 of file EulerJacobian.hpp.

Member Data Documentation

◆ LDU

template<int nVarsFixed>
tLocalMat DNDS::Euler::JacobianLocalLU< nVarsFixed >::LDU

Sparse row storage for D, L, U blocks.

Each row i stores blocks in the order: [0] = diag, [1, 1 + symLU->lowerTriStructure[i].size()) = lower-triangular, [1 + symLU->lowerTriStructure[i].size(), end) = upper-triangular.

Definition at line 290 of file EulerJacobian.hpp.


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