Loading...
Searching...
No Matches
HyperReduction< SnapshotsLists > Class Template Reference

Public Types

using SnapshotsListTuple = std::tuple<std::decay_t<SnapshotsLists>...>
using FieldsTuple
template<std::size_t N>
using NthFieldListType
template<std::size_t N>
using NthFieldType = typename NthFieldListType<N>::value_type

Public Member Functions

 HyperReduction (label n_modes, label n_nodes, Eigen::VectorXi initialSeeds, word problemName, SnapshotsLists &&...snapshotsLists)
 Construct HyperReduction class, interpolation-based.
 HyperReduction (label n_modes, label n_nodes, unsigned int vectorial_dim, label n_cells, Eigen::VectorXi initialSeeds, word problemName, Eigen::VectorXd volumes=Eigen::VectorXd::Constant(1, 0.0))
 Construct HyperReduction class.
void offlineGappyDEIM (Eigen::MatrixXd &snapshotsModes, Eigen::VectorXd &normalizingWeights)
 Methods implemented: 'GappyDEIM' from "DEIM, Chaturantabut, Saifon, and Danny C. Sorensen. "Nonlinear model reduction via discrete empirical interpolation.
void offlineGappyDEIM (Eigen::MatrixXd &snapshotsModes, Eigen::VectorXd &normalizingWeights, word folderMethodName)
void offlineECP (Eigen::MatrixXd &snapshotsModes, Eigen::VectorXd &normalizingWeights)
 Methods implemented: 'ECP' from "ECP, Hernandez, Joaquin Alberto, Manuel Alejandro Caicedo, and Alex Ferrer. "Dimensional hyper-reduction of nonlinear finite element models via empirical cubature.
void offlineECP (Eigen::MatrixXd &snapshotsModes, Eigen::VectorXd &normalizingWeights, word folderMethodName)
void offlineECP (List< Eigen::MatrixXd > &listSnapshotsModes, Eigen::VectorXd &weights, word folderMethodName)
 Methods implemented: close to 'SAW-ECM' from "J.R.Bravo, J.A.Hernández, S.Ares de Parga, R.Rossi. "A subspace-adaptive weights cubature method with application to the local hyperreduction of parameterized finite element models".
void initSeeds (Eigen::VectorXd mp_not_mask, std::set< label > nodePointsSet)
 TODO.
void computeLS (Eigen::MatrixXd &J, Eigen::MatrixXd &JWhole, Eigen::VectorXd &b, Eigen::VectorXd &q, bool NNLS)
 TODO.
void getModesSVD (SnapshotsListTuple &SnapshotsListTuple, Eigen::MatrixXd &modesSVD, Eigen::VectorXd &fieldWeights, bool saveModesFlag=false)
 TODO.
void getModesSVD (SnapshotsListTuple &snapshotsListTuple, Eigen::MatrixXd &modesSVD, Eigen::VectorXd &fieldWeights, Eigen::MatrixXd &modesSVDBoundary, Eigen::VectorXd &fieldWeightsBoundary, bool saveModesFlag)
 TODO.
void updateNodes (Eigen::SparseMatrix< double > &P, label &ind_max, Eigen::VectorXd &mp_not_mask)
 TODO.
void updateNodesRemove (Eigen::SparseMatrix< double > &P, Eigen::VectorXd weights)
 Remove from the nodes the weights equal to 0 (used with NNLS).
template<typename SnapshotsList>
void stackSnapshots (SnapshotsList sList, Eigen::MatrixXd &snapshotsMatrix, Eigen::VectorXd &fieldWeights)
 TODO.
template<typename SnapshotsList>
void stackSnapshotsBoundary (SnapshotsList sList, List< Eigen::MatrixXd > &snapshotsMatrixBoundary, List< Eigen::VectorXd > &fieldWeightsBoundary)
 TODO.
template<typename SnapshotsList>
void saveModes (SnapshotsList sList, Eigen::MatrixXd &snapshotsMatrix, unsigned int &rowIndex, unsigned int &modeIndex, word folder)
 TODO.
template<typename SnapshotsList>
void saveModes (SnapshotsList sList, Eigen::MatrixXd &snapshotsMatrix, Eigen::MatrixXd &snapshotsMatrixBoundary, unsigned int &rowIndex, unsigned int &rowIndexBoundary, unsigned int &modeIndex, word folder)
 TODO.
template<typename SnapshotsList>
void stackNames (SnapshotsList sList)
 TODO.
template<typename SnapshotsList>
void stackDimensions (SnapshotsList sList)
 TODO.
template<typename SnapshotsList>
void sumDimensions (double sum, SnapshotsList sList)
 TODO.
void evaluatePinv (Eigen::SparseMatrix< double > &Projector, Eigen::MatrixXd &Modes, Eigen::VectorXd &fieldWeights)
 Compute the pseudo-inverse of the matrix M restricted with the projector P.
void evaluateWPU (Eigen::SparseMatrix< double > &Projector, Eigen::MatrixXd &Modes, Eigen::VectorXd &fieldWeights, Eigen::VectorXd &quadratureWeights)
 Compute the pseudo-inverse of the matrix M restricted with the projector P.
void generateSubmesh (label layers, const fvMesh &mesh)
 Compute the submesh common to all fields in SnapshotsLists.
List< label > global2local (List< label > &points, fvMeshSubset &submesh)
 Get local indices in the submesh from indices in the global ones.
template<typename FieldType>
autoPtr< FieldType > interpolateField (const FieldType &field)
 TODO.
template<typename Field>
constexpr unsigned int get_field_dim ()
 TODO.
template<typename LastList>
constexpr unsigned int compute_vectorial_dim (LastList x)
 TODO.
template<typename List, typename... RemainingLists>
constexpr unsigned int compute_vectorial_dim (List &&head, RemainingLists &&...tail)
 TODO.
template<typename... FieldsArgs>
void eigen2fields (Eigen::VectorXd &eFields, FieldsArgs &&... oFields)
 TODO.
void initReshapeMat (Eigen::SparseMatrix< double > &reshapeMat)
 TODO.
void createMasks (bool offlineStage=true)
 TODO.
void getSnapMatrix (Eigen::MatrixXd &snapMatrix, Eigen::VectorXd &fieldWeights)
 TODO.
void getSnapMatrix (Eigen::MatrixXd &snapMatrix, Eigen::VectorXd &fieldWeights, List< Eigen::MatrixXd > &snapMatrixBoundary, List< Eigen::VectorXd > &fieldWeightsBoundary)
 TODO.
Eigen::SparseMatrix< double > maskToOtherDim (int newDim)
 Compute the sparse mask matrix P with a vectorial dimension different from the snapshots (enlarged or reduced).

Static Public Member Functions

static double s_optimality (Eigen::MatrixXd &A)
 TODO.

Public Attributes

word methodName
unsigned int vectorial_dim
unsigned int sumFieldsDim
ITHACAparameterspara
label n_boundary_patches
List< label > n_boundary_cells_list
label n_boundary_cells
label n_modes
 The maximum number of modes to be considered.
label n_nodes
 The maximum number of modes to be considered.
label n_snapshots
 The length of the snapshots lists.
Eigen::VectorXi initialSeeds
Eigen::VectorXi nodes
word problemName
 The name of the non-linear function e.g. HR_method/residual.
SnapshotsListTuple snapshotsListTuple
 The snapshots matrix containing the nonlinear function or operator.
List< word > fieldNames
 Names of the fields.
List< unsigned int > fieldDims
 Dimensions of the fields.
label n_cells
 Int Number of Cells;.
Eigen::VectorXd volumes
 Volume of each cell.
label n_cellsSubfields
 Int Number of Cells in submeshes;.
autoPtr< IOList< label > > nodePoints
 Nodes in the case of the a nonlinear function.
autoPtr< IOList< labelList > > totalNodePoints
 List of label lists of the nodes and corresponding surrounding nodes.
autoPtr< IOList< label > > uniqueNodePoints
 List of the unique indices of the nodes that define the submesh.
word folderProblem
 Folder for the HR problem.
word folderMethod
 Folder for the selected HR method.
List< label > localNodePoints
 Indices of the local node points in the subMesh.
Eigen::MatrixXd basisMatrix
 Orthonormal basis of HR.
Eigen::VectorXd normalizingWeights
Eigen::MatrixXd renormalizedBasisMatrix
 Renormalized basis of HR.
Eigen::MatrixXd pinvPU
Eigen::MatrixXd wPU
Eigen::VectorXd eigenValueseig
Eigen::SparseMatrix< double > P
 The P matrix of the HyperReduction method. The nodes are ordered in the order of insertion during the greedy procedure. The components related to the same node follow immediately.
Eigen::SparseMatrix< double > field2submesh
Eigen::SparseMatrix< double > submesh2nodes
Eigen::VectorXi submesh2nodesMask
Eigen::VectorXd quadratureWeights
 Quadrature weights. Ordered in the same order of matrix P.
List< Eigen::VectorXd > quadratureWeightsSubspaces
autoPtr< fvMeshSubset > submesh
 Submesh of the HyperReduction method.
autoPtr< volVectorField > submesh_field
 Submeshes.
Eigen::MatrixXd MatrixOnline
 Online Matrix.

Static Public Attributes

static constexpr auto n_fields = sizeof...(SnapshotsLists)

Detailed Description

template<typename... SnapshotsLists>
class HyperReduction< SnapshotsLists >

Definition at line 48 of file hyperReduction.H.

Member Typedef Documentation

◆ FieldsTuple

template<typename... SnapshotsLists>
using HyperReduction< SnapshotsLists >::FieldsTuple
Initial value:
std::tuple<typename std::decay_t<SnapshotsLists>::value_type ...>

Definition at line 54 of file hyperReduction.H.

◆ NthFieldListType

template<typename... SnapshotsLists>
template<std::size_t N>
using HyperReduction< SnapshotsLists >::NthFieldListType
Initial value:
typename
std::tuple_element<N, SnapshotsListTuple>::type

Definition at line 63 of file hyperReduction.H.

◆ NthFieldType

template<typename... SnapshotsLists>
template<std::size_t N>
using HyperReduction< SnapshotsLists >::NthFieldType = typename NthFieldListType<N>::value_type

Definition at line 67 of file hyperReduction.H.

◆ SnapshotsListTuple

template<typename... SnapshotsLists>
using HyperReduction< SnapshotsLists >::SnapshotsListTuple = std::tuple<std::decay_t<SnapshotsLists>...>

Definition at line 53 of file hyperReduction.H.

Constructor & Destructor Documentation

◆ HyperReduction() [1/2]

template<typename... SnapshotsLists>
HyperReduction< SnapshotsLists >::HyperReduction ( label n_modes,
label n_nodes,
Eigen::VectorXi initialSeeds,
word problemName,
SnapshotsLists &&... snapshotsLists )

Construct HyperReduction class, interpolation-based.

SnapshotsLists is a variadic argument of PtrLists of fields e.g. for incompressible Navier-Stokes PtrList<volVectorField>, PtrList<volScalarField> for velocity and pressure.

Parameters
[in]hrMethodthe chosen HR method
[in]n_modesdimension of the HR basis
[in]n_nodesnumber of nodes
[in]initialSeedsinitial nodes
[in]problemNamename of the function to be hyper-reduced
[in]snapshotsListslists of fields

Definition at line 40 of file hyperReduction.templates.H.

◆ HyperReduction() [2/2]

template<typename... SnapshotsLists>
HyperReduction< SnapshotsLists >::HyperReduction ( label n_modes,
label n_nodes,
unsigned int vectorial_dim,
label n_cells,
Eigen::VectorXi initialSeeds,
word problemName,
Eigen::VectorXd volumes = Eigen::VectorXd::Constant(1,0.0) )

Construct HyperReduction class.

Parameters
[in]n_modesDimension of the HR basis
[in]n_nodesNumber of nodes
[in]vectorial_dimField dimension (1 for volScalarField 3 for volVectorField)
[in]n_cellsNumber of cells
[in]initialSeedsInitial nodes
[in]problemNameName of the function to be hyper-reduced
[in]volumesVolume of each cell (optional, used in offlineECP)

Definition at line 110 of file hyperReduction.templates.H.

Member Function Documentation

◆ compute_vectorial_dim() [1/2]

template<typename... SnapshotsLists>
template<typename LastList>
unsigned int HyperReduction< SnapshotsLists >::compute_vectorial_dim ( LastList x)
inlineconstexpr

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 449 of file hyperReduction.H.

◆ compute_vectorial_dim() [2/2]

template<typename... SnapshotsLists>
template<typename List, typename... RemainingLists>
unsigned int HyperReduction< SnapshotsLists >::compute_vectorial_dim ( List && head,
RemainingLists &&... tail )
inlineconstexpr

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 460 of file hyperReduction.H.

◆ computeLS()

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::computeLS ( Eigen::MatrixXd & J,
Eigen::MatrixXd & JWhole,
Eigen::VectorXd & b,
Eigen::VectorXd & q,
bool NNLS )

TODO.

Definition at line 1034 of file hyperReduction.templates.H.

◆ createMasks()

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::createMasks ( bool offlineStage = true)

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 1226 of file hyperReduction.templates.H.

◆ eigen2fields()

template<typename... SnapshotsLists>
template<typename... FieldsArgs>
void HyperReduction< SnapshotsLists >::eigen2fields ( Eigen::VectorXd & eFields,
FieldsArgs &&... oFields )
inline

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 486 of file hyperReduction.H.

◆ evaluatePinv()

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::evaluatePinv ( Eigen::SparseMatrix< double > & Projector,
Eigen::MatrixXd & Modes,
Eigen::VectorXd & fieldWeights )

Compute the pseudo-inverse of the matrix M restricted with the projector P.

Parameters
[in]Projectorprojector into interpolation nodes
[in]Modesmatrix to restrict and invert

Definition at line 1122 of file hyperReduction.templates.H.

◆ evaluateWPU()

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::evaluateWPU ( Eigen::SparseMatrix< double > & Projector,
Eigen::MatrixXd & Modes,
Eigen::VectorXd & fieldWeights,
Eigen::VectorXd & quadratureWeights )

Compute the pseudo-inverse of the matrix M restricted with the projector P.

Parameters
[in]Projectorprojector into interpolation nodes
[in]Modesmatrix to restrict and invert

Definition at line 1136 of file hyperReduction.templates.H.

◆ generateSubmesh()

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::generateSubmesh ( label layers,
const fvMesh & mesh )

Compute the submesh common to all fields in SnapshotsLists.

Parameters
[in]layersprojector into interpolation nodes
[in]meshmatrix to restrict and invert

Definition at line 1161 of file hyperReduction.templates.H.

◆ get_field_dim()

template<typename... SnapshotsLists>
template<typename Field>
unsigned int HyperReduction< SnapshotsLists >::get_field_dim ( )
inlineconstexpr

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 438 of file hyperReduction.H.

◆ getModesSVD() [1/2]

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::getModesSVD ( SnapshotsListTuple & SnapshotsListTuple,
Eigen::MatrixXd & modesSVD,
Eigen::VectorXd & fieldWeights,
bool saveModesFlag = false )

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 134 of file hyperReduction.templates.H.

◆ getModesSVD() [2/2]

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::getModesSVD ( SnapshotsListTuple & snapshotsListTuple,
Eigen::MatrixXd & modesSVD,
Eigen::VectorXd & fieldWeights,
Eigen::MatrixXd & modesSVDBoundary,
Eigen::VectorXd & fieldWeightsBoundary,
bool saveModesFlag )

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 224 of file hyperReduction.templates.H.

◆ getSnapMatrix() [1/2]

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::getSnapMatrix ( Eigen::MatrixXd & snapMatrix,
Eigen::VectorXd & fieldWeights )

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 344 of file hyperReduction.templates.H.

◆ getSnapMatrix() [2/2]

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::getSnapMatrix ( Eigen::MatrixXd & snapMatrix,
Eigen::VectorXd & fieldWeights,
List< Eigen::MatrixXd > & snapMatrixBoundary,
List< Eigen::VectorXd > & fieldWeightsBoundary )

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 354 of file hyperReduction.templates.H.

◆ global2local()

template<typename... SnapshotsLists>
List< label > HyperReduction< SnapshotsLists >::global2local ( List< label > & points,
fvMeshSubset & submesh )

Get local indices in the submesh from indices in the global ones.

Parameters
pointsThe points
submeshThe submesh
Returns
The local indices

Definition at line 1292 of file hyperReduction.templates.H.

◆ initReshapeMat()

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::initReshapeMat ( Eigen::SparseMatrix< double > & reshapeMat)

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 1056 of file hyperReduction.templates.H.

◆ initSeeds()

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::initSeeds ( Eigen::VectorXd mp_not_mask,
std::set< label > nodePointsSet )

TODO.

Definition at line 925 of file hyperReduction.templates.H.

◆ interpolateField()

template<typename... SnapshotsLists>
template<typename FieldType>
autoPtr< FieldType > HyperReduction< SnapshotsLists >::interpolateField ( const FieldType & field)
inline

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 426 of file hyperReduction.H.

◆ maskToOtherDim()

template<typename... SnapshotsLists>
Eigen::SparseMatrix< double > HyperReduction< SnapshotsLists >::maskToOtherDim ( int newDim)

Compute the sparse mask matrix P with a vectorial dimension different from the snapshots (enlarged or reduced).

Parameters
[in]newDimint, the vectorial dim for the output

Definition at line 1313 of file hyperReduction.templates.H.

◆ offlineECP() [1/3]

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::offlineECP ( Eigen::MatrixXd & snapshotsModes,
Eigen::VectorXd & normalizingWeights )
inline

Methods implemented: 'ECP' from "ECP, Hernandez, Joaquin Alberto, Manuel Alejandro Caicedo, and Alex Ferrer. "Dimensional hyper-reduction of nonlinear finite element models via empirical cubature.

" Computer methods in applied mechanics and engineering 313 (2017): 687-722.".

Definition at line 231 of file hyperReduction.H.

◆ offlineECP() [2/3]

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::offlineECP ( Eigen::MatrixXd & snapshotsModes,
Eigen::VectorXd & normalizingWeights,
word folderMethodName )
inline

Definition at line 239 of file hyperReduction.H.

◆ offlineECP() [3/3]

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::offlineECP ( List< Eigen::MatrixXd > & listSnapshotsModes,
Eigen::VectorXd & weights,
word folderMethodName )

Methods implemented: close to 'SAW-ECM' from "J.R.Bravo, J.A.Hernández, S.Ares de Parga, R.Rossi. "A subspace-adaptive weights cubature method with application to the local hyperreduction of parameterized finite element models".

(2024). and also close to the above method when vectorial_dim > 1 Snapshots are grouped into subspaces. HR nodes are shared by all subspaces. HR weights can differ.

Definition at line 676 of file hyperReduction.templates.H.

◆ offlineGappyDEIM() [1/2]

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::offlineGappyDEIM ( Eigen::MatrixXd & snapshotsModes,
Eigen::VectorXd & normalizingWeights )
inline

Methods implemented: 'GappyDEIM' from "DEIM, Chaturantabut, Saifon, and Danny C. Sorensen. "Nonlinear model reduction via discrete empirical interpolation.

" SIAM Journal on Scientific Computing 32.5 (2010): 2737-2764" and "GNAT, Carlberg, Kevin, et al. "The GNAT method for nonlinear model reduction: effective implementation and application to computational fluid dynamics and turbulent flows." Journal of Computational Physics 242 (2013): 623-647.". Quasi-optimal point selection algorithm 'SOPT' and optimal point selection 'SOPTE' from "Lauzon, Jessica T., Siu Wun Cheung, Yeonjong Shin, Youngsoo Choi, Dylan Matthew Copeland, and Kevin Huynh. "S-OPT: A Points Selection Algorithm for Hyper-Reduction in Reduced Order Models." arXiv preprint arXiv:2203.16494 (2022).""

Definition at line 217 of file hyperReduction.H.

◆ offlineGappyDEIM() [2/2]

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::offlineGappyDEIM ( Eigen::MatrixXd & snapshotsModes,
Eigen::VectorXd & normalizingWeights,
word folderMethodName )

Definition at line 417 of file hyperReduction.templates.H.

◆ s_optimality()

template<typename... SnapshotsLists>
double HyperReduction< SnapshotsLists >::s_optimality ( Eigen::MatrixXd & A)
inlinestatic

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 472 of file hyperReduction.H.

◆ saveModes() [1/2]

template<typename... SnapshotsLists>
template<typename SnapshotsList>
void HyperReduction< SnapshotsLists >::saveModes ( SnapshotsList sList,
Eigen::MatrixXd & snapshotsMatrix,
Eigen::MatrixXd & snapshotsMatrixBoundary,
unsigned int & rowIndex,
unsigned int & rowIndexBoundary,
unsigned int & modeIndex,
word folder )

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 389 of file hyperReduction.templates.H.

◆ saveModes() [2/2]

template<typename... SnapshotsLists>
template<typename SnapshotsList>
void HyperReduction< SnapshotsLists >::saveModes ( SnapshotsList sList,
Eigen::MatrixXd & snapshotsMatrix,
unsigned int & rowIndex,
unsigned int & modeIndex,
word folder )

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 373 of file hyperReduction.templates.H.

◆ stackDimensions()

template<typename... SnapshotsLists>
template<typename SnapshotsList>
void HyperReduction< SnapshotsLists >::stackDimensions ( SnapshotsList sList)
inline

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 366 of file hyperReduction.H.

◆ stackNames()

template<typename... SnapshotsLists>
template<typename SnapshotsList>
void HyperReduction< SnapshotsLists >::stackNames ( SnapshotsList sList)
inline

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 354 of file hyperReduction.H.

◆ stackSnapshots()

template<typename... SnapshotsLists>
template<typename SnapshotsList>
void HyperReduction< SnapshotsLists >::stackSnapshots ( SnapshotsList sList,
Eigen::MatrixXd & snapshotsMatrix,
Eigen::VectorXd & fieldWeights )

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 1076 of file hyperReduction.templates.H.

◆ stackSnapshotsBoundary()

template<typename... SnapshotsLists>
template<typename SnapshotsList>
void HyperReduction< SnapshotsLists >::stackSnapshotsBoundary ( SnapshotsList sList,
List< Eigen::MatrixXd > & snapshotsMatrixBoundary,
List< Eigen::VectorXd > & fieldWeightsBoundary )

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 1094 of file hyperReduction.templates.H.

◆ sumDimensions()

template<typename... SnapshotsLists>
template<typename SnapshotsList>
void HyperReduction< SnapshotsLists >::sumDimensions ( double sum,
SnapshotsList sList )
inline

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 378 of file hyperReduction.H.

◆ updateNodes()

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::updateNodes ( Eigen::SparseMatrix< double > & P,
label & ind_max,
Eigen::VectorXd & mp_not_mask )

TODO.

Template Parameters
SnapshotsListPtrList of volScalarField or volVectorField
Parameters
[in]sListThe list of snapshots

Definition at line 963 of file hyperReduction.templates.H.

◆ updateNodesRemove()

template<typename... SnapshotsLists>
void HyperReduction< SnapshotsLists >::updateNodesRemove ( Eigen::SparseMatrix< double > & P,
Eigen::VectorXd weights )

Remove from the nodes the weights equal to 0 (used with NNLS).

Template Parameters
PThe mask matrix of the nodes, updated with the call
Parameters
[in]weightsA Eigen::VectorXd of quadrature weights (for comparison with 0)

Definition at line 999 of file hyperReduction.templates.H.

Member Data Documentation

◆ basisMatrix

template<typename... SnapshotsLists>
Eigen::MatrixXd HyperReduction< SnapshotsLists >::basisMatrix

Orthonormal basis of HR.

Definition at line 172 of file hyperReduction.H.

◆ eigenValueseig

template<typename... SnapshotsLists>
Eigen::VectorXd HyperReduction< SnapshotsLists >::eigenValueseig

Definition at line 187 of file hyperReduction.H.

◆ field2submesh

template<typename... SnapshotsLists>
Eigen::SparseMatrix<double> HyperReduction< SnapshotsLists >::field2submesh

Definition at line 193 of file hyperReduction.H.

◆ fieldDims

template<typename... SnapshotsLists>
List<unsigned int> HyperReduction< SnapshotsLists >::fieldDims

Dimensions of the fields.

Definition at line 141 of file hyperReduction.H.

◆ fieldNames

template<typename... SnapshotsLists>
List<word> HyperReduction< SnapshotsLists >::fieldNames

Names of the fields.

Definition at line 138 of file hyperReduction.H.

◆ folderMethod

template<typename... SnapshotsLists>
word HyperReduction< SnapshotsLists >::folderMethod

Folder for the selected HR method.

Definition at line 165 of file hyperReduction.H.

◆ folderProblem

template<typename... SnapshotsLists>
word HyperReduction< SnapshotsLists >::folderProblem

Folder for the HR problem.

Definition at line 162 of file hyperReduction.H.

◆ initialSeeds

template<typename... SnapshotsLists>
Eigen::VectorXi HyperReduction< SnapshotsLists >::initialSeeds

Definition at line 126 of file hyperReduction.H.

◆ localNodePoints

template<typename... SnapshotsLists>
List<label> HyperReduction< SnapshotsLists >::localNodePoints

Indices of the local node points in the subMesh.

Definition at line 168 of file hyperReduction.H.

◆ MatrixOnline

template<typename... SnapshotsLists>
Eigen::MatrixXd HyperReduction< SnapshotsLists >::MatrixOnline

Online Matrix.

Definition at line 212 of file hyperReduction.H.

◆ methodName

template<typename... SnapshotsLists>
word HyperReduction< SnapshotsLists >::methodName

Definition at line 51 of file hyperReduction.H.

◆ n_boundary_cells

template<typename... SnapshotsLists>
label HyperReduction< SnapshotsLists >::n_boundary_cells

Definition at line 114 of file hyperReduction.H.

◆ n_boundary_cells_list

template<typename... SnapshotsLists>
List<label> HyperReduction< SnapshotsLists >::n_boundary_cells_list

Definition at line 111 of file hyperReduction.H.

◆ n_boundary_patches

template<typename... SnapshotsLists>
label HyperReduction< SnapshotsLists >::n_boundary_patches

Definition at line 108 of file hyperReduction.H.

◆ n_cells

template<typename... SnapshotsLists>
label HyperReduction< SnapshotsLists >::n_cells

Int Number of Cells;.

Definition at line 144 of file hyperReduction.H.

◆ n_cellsSubfields

template<typename... SnapshotsLists>
label HyperReduction< SnapshotsLists >::n_cellsSubfields

Int Number of Cells in submeshes;.

Definition at line 150 of file hyperReduction.H.

◆ n_fields

template<typename... SnapshotsLists>
auto HyperReduction< SnapshotsLists >::n_fields = sizeof...(SnapshotsLists)
staticconstexpr

Definition at line 56 of file hyperReduction.H.

◆ n_modes

template<typename... SnapshotsLists>
label HyperReduction< SnapshotsLists >::n_modes

The maximum number of modes to be considered.

Definition at line 117 of file hyperReduction.H.

◆ n_nodes

template<typename... SnapshotsLists>
label HyperReduction< SnapshotsLists >::n_nodes

The maximum number of modes to be considered.

Definition at line 120 of file hyperReduction.H.

◆ n_snapshots

template<typename... SnapshotsLists>
label HyperReduction< SnapshotsLists >::n_snapshots

The length of the snapshots lists.

Definition at line 123 of file hyperReduction.H.

◆ nodePoints

template<typename... SnapshotsLists>
autoPtr<IOList<label > > HyperReduction< SnapshotsLists >::nodePoints

Nodes in the case of the a nonlinear function.

Definition at line 153 of file hyperReduction.H.

◆ nodes

template<typename... SnapshotsLists>
Eigen::VectorXi HyperReduction< SnapshotsLists >::nodes

Definition at line 129 of file hyperReduction.H.

◆ normalizingWeights

template<typename... SnapshotsLists>
Eigen::VectorXd HyperReduction< SnapshotsLists >::normalizingWeights

Definition at line 175 of file hyperReduction.H.

◆ P

template<typename... SnapshotsLists>
Eigen::SparseMatrix<double> HyperReduction< SnapshotsLists >::P

The P matrix of the HyperReduction method. The nodes are ordered in the order of insertion during the greedy procedure. The components related to the same node follow immediately.

Definition at line 190 of file hyperReduction.H.

◆ para

template<typename... SnapshotsLists>
ITHACAparameters* HyperReduction< SnapshotsLists >::para

Definition at line 105 of file hyperReduction.H.

◆ pinvPU

template<typename... SnapshotsLists>
Eigen::MatrixXd HyperReduction< SnapshotsLists >::pinvPU

Definition at line 181 of file hyperReduction.H.

◆ problemName

template<typename... SnapshotsLists>
word HyperReduction< SnapshotsLists >::problemName

The name of the non-linear function e.g. HR_method/residual.

Definition at line 132 of file hyperReduction.H.

◆ quadratureWeights

template<typename... SnapshotsLists>
Eigen::VectorXd HyperReduction< SnapshotsLists >::quadratureWeights

Quadrature weights. Ordered in the same order of matrix P.

Definition at line 202 of file hyperReduction.H.

◆ quadratureWeightsSubspaces

template<typename... SnapshotsLists>
List<Eigen::VectorXd> HyperReduction< SnapshotsLists >::quadratureWeightsSubspaces

Definition at line 203 of file hyperReduction.H.

◆ renormalizedBasisMatrix

template<typename... SnapshotsLists>
Eigen::MatrixXd HyperReduction< SnapshotsLists >::renormalizedBasisMatrix

Renormalized basis of HR.

Definition at line 178 of file hyperReduction.H.

◆ snapshotsListTuple

template<typename... SnapshotsLists>
SnapshotsListTuple HyperReduction< SnapshotsLists >::snapshotsListTuple

The snapshots matrix containing the nonlinear function or operator.

Definition at line 135 of file hyperReduction.H.

◆ submesh

template<typename... SnapshotsLists>
autoPtr<fvMeshSubset> HyperReduction< SnapshotsLists >::submesh

Submesh of the HyperReduction method.

Definition at line 206 of file hyperReduction.H.

◆ submesh2nodes

template<typename... SnapshotsLists>
Eigen::SparseMatrix<double> HyperReduction< SnapshotsLists >::submesh2nodes

Definition at line 196 of file hyperReduction.H.

◆ submesh2nodesMask

template<typename... SnapshotsLists>
Eigen::VectorXi HyperReduction< SnapshotsLists >::submesh2nodesMask

Definition at line 199 of file hyperReduction.H.

◆ submesh_field

template<typename... SnapshotsLists>
autoPtr<volVectorField> HyperReduction< SnapshotsLists >::submesh_field

Submeshes.

Definition at line 209 of file hyperReduction.H.

◆ sumFieldsDim

template<typename... SnapshotsLists>
unsigned int HyperReduction< SnapshotsLists >::sumFieldsDim

Definition at line 60 of file hyperReduction.H.

◆ totalNodePoints

template<typename... SnapshotsLists>
autoPtr<IOList<labelList > > HyperReduction< SnapshotsLists >::totalNodePoints

List of label lists of the nodes and corresponding surrounding nodes.

Definition at line 156 of file hyperReduction.H.

◆ uniqueNodePoints

template<typename... SnapshotsLists>
autoPtr<IOList<label > > HyperReduction< SnapshotsLists >::uniqueNodePoints

List of the unique indices of the nodes that define the submesh.

Definition at line 159 of file hyperReduction.H.

◆ vectorial_dim

template<typename... SnapshotsLists>
unsigned int HyperReduction< SnapshotsLists >::vectorial_dim

Definition at line 59 of file hyperReduction.H.

◆ volumes

template<typename... SnapshotsLists>
Eigen::VectorXd HyperReduction< SnapshotsLists >::volumes

Volume of each cell.

Definition at line 147 of file hyperReduction.H.

◆ wPU

template<typename... SnapshotsLists>
Eigen::MatrixXd HyperReduction< SnapshotsLists >::wPU

Definition at line 184 of file hyperReduction.H.


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