#include <hyperReduction.H>
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) | |
| 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 | 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) |
| 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. | |
| 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. | |
Static Public Member Functions | |
| static double | s_optimality (Eigen::MatrixXd &A) |
| TODO. | |
Public Attributes | |
| word | methodName |
| unsigned int | vectorial_dim |
| unsigned int | sumFieldsDim |
| ITHACAparameters * | para |
| 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;. | |
| 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. | |
| 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) |
Definition at line 47 of file hyperReduction.H.
| using HyperReduction< SnapshotsLists >::FieldsTuple |
Definition at line 53 of file hyperReduction.H.
| using HyperReduction< SnapshotsLists >::NthFieldListType |
Definition at line 62 of file hyperReduction.H.
| using HyperReduction< SnapshotsLists >::NthFieldType = typename NthFieldListType<N>::value_type |
Definition at line 66 of file hyperReduction.H.
| using HyperReduction< SnapshotsLists >::SnapshotsListTuple = std::tuple<std::decay_t<SnapshotsLists>...> |
Definition at line 52 of file hyperReduction.H.
| 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.
| [in] | hrMethod | the chosen HR method |
| [in] | n_modes | dimension of the HR basis |
| [in] | n_nodes | number of nodes |
| [in] | initialSeeds | initial nodes |
| [in] | problemName | name of the function to be hyper-reduced |
| [in] | snapshotsLists | lists of fields |
Definition at line 40 of file hyperReduction.templates.H.
| HyperReduction< SnapshotsLists >::HyperReduction | ( | label | n_modes, |
| label | n_nodes, | ||
| unsigned int | vectorial_dim, | ||
| label | n_cells, | ||
| Eigen::VectorXi | initialSeeds, | ||
| word | problemName ) |
Construct HyperReduction class.
| [in] | n_modes | Dimension of the HR basis |
| [in] | n_nodes | Number of nodes |
| [in] | vectorial_dim | Field dimension (1 for volScalarField 3 for volVectorField) |
| [in] | n_cells | Number of cells |
| [in] | initialSeeds | Initial nodes |
| [in] | problemName | Name of the function to be hyper-reduced |
Definition at line 111 of file hyperReduction.templates.H.
|
inlineconstexpr |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 420 of file hyperReduction.H.
|
inlineconstexpr |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 431 of file hyperReduction.H.
| void HyperReduction< SnapshotsLists >::computeLS | ( | Eigen::MatrixXd & | J, |
| Eigen::MatrixXd & | JWhole, | ||
| Eigen::VectorXd & | b, | ||
| Eigen::VectorXd & | q ) |
TODO.
Definition at line 854 of file hyperReduction.templates.H.
| void HyperReduction< SnapshotsLists >::createMasks | ( | bool | offlineStage = true | ) |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 1047 of file hyperReduction.templates.H.
|
inline |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 457 of file hyperReduction.H.
| 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.
| [in] | Projector | projector into interpolation nodes |
| [in] | Modes | matrix to restrict and invert |
Definition at line 944 of file hyperReduction.templates.H.
| 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.
| [in] | Projector | projector into interpolation nodes |
| [in] | Modes | matrix to restrict and invert |
Definition at line 958 of file hyperReduction.templates.H.
| void HyperReduction< SnapshotsLists >::generateSubmesh | ( | label | layers, |
| const fvMesh & | mesh ) |
Compute the submesh common to all fields in SnapshotsLists.
| [in] | layers | projector into interpolation nodes |
| [in] | mesh | matrix to restrict and invert |
Definition at line 983 of file hyperReduction.templates.H.
|
inlineconstexpr |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 409 of file hyperReduction.H.
| void HyperReduction< SnapshotsLists >::getModesSVD | ( | SnapshotsListTuple & | SnapshotsListTuple, |
| Eigen::MatrixXd & | modesSVD, | ||
| Eigen::VectorXd & | fieldWeights, | ||
| bool | saveModesFlag = false ) |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 133 of file hyperReduction.templates.H.
| void HyperReduction< SnapshotsLists >::getModesSVD | ( | SnapshotsListTuple & | snapshotsListTuple, |
| Eigen::MatrixXd & | modesSVD, | ||
| Eigen::VectorXd & | fieldWeights, | ||
| Eigen::MatrixXd & | modesSVDBoundary, | ||
| Eigen::VectorXd & | fieldWeightsBoundary, | ||
| bool | saveModesFlag ) |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 224 of file hyperReduction.templates.H.
| void HyperReduction< SnapshotsLists >::getSnapMatrix | ( | Eigen::MatrixXd & | snapMatrix, |
| Eigen::VectorXd & | fieldWeights ) |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 345 of file hyperReduction.templates.H.
| void HyperReduction< SnapshotsLists >::getSnapMatrix | ( | Eigen::MatrixXd & | snapMatrix, |
| Eigen::VectorXd & | fieldWeights, | ||
| List< Eigen::MatrixXd > & | snapMatrixBoundary, | ||
| List< Eigen::VectorXd > & | fieldWeightsBoundary ) |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 355 of file hyperReduction.templates.H.
| List< label > HyperReduction< SnapshotsLists >::global2local | ( | List< label > & | points, |
| fvMeshSubset & | submesh ) |
Get local indices in the submesh from indices in the global ones.
| points | The points |
| submesh | The submesh |
Definition at line 1113 of file hyperReduction.templates.H.
| void HyperReduction< SnapshotsLists >::initReshapeMat | ( | Eigen::SparseMatrix< double > & | reshapeMat | ) |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 878 of file hyperReduction.templates.H.
| void HyperReduction< SnapshotsLists >::initSeeds | ( | Eigen::VectorXd | mp_not_mask, |
| std::set< label > | nodePointsSet ) |
TODO.
Definition at line 782 of file hyperReduction.templates.H.
|
inline |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 397 of file hyperReduction.H.
|
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.". TODO implement non-negative LS.
Definition at line 224 of file hyperReduction.H.
| void HyperReduction< SnapshotsLists >::offlineECP | ( | Eigen::MatrixXd & | snapshotsModes, |
| Eigen::VectorXd & | normalizingWeights, | ||
| word | folderMethodName ) |
Definition at line 678 of file hyperReduction.templates.H.
|
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 210 of file hyperReduction.H.
| void HyperReduction< SnapshotsLists >::offlineGappyDEIM | ( | Eigen::MatrixXd & | snapshotsModes, |
| Eigen::VectorXd & | normalizingWeights, | ||
| word | folderMethodName ) |
Definition at line 419 of file hyperReduction.templates.H.
|
inlinestatic |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 443 of file hyperReduction.H.
| void HyperReduction< SnapshotsLists >::saveModes | ( | SnapshotsList | sList, |
| Eigen::MatrixXd & | snapshotsMatrix, | ||
| Eigen::MatrixXd & | snapshotsMatrixBoundary, | ||
| unsigned int & | rowIndex, | ||
| unsigned int & | rowIndexBoundary, | ||
| unsigned int & | modeIndex, | ||
| word | folder ) |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 391 of file hyperReduction.templates.H.
| void HyperReduction< SnapshotsLists >::saveModes | ( | SnapshotsList | sList, |
| Eigen::MatrixXd & | snapshotsMatrix, | ||
| unsigned int & | rowIndex, | ||
| unsigned int & | modeIndex, | ||
| word | folder ) |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 375 of file hyperReduction.templates.H.
|
inline |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 337 of file hyperReduction.H.
|
inline |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 325 of file hyperReduction.H.
| void HyperReduction< SnapshotsLists >::stackSnapshots | ( | SnapshotsList | sList, |
| Eigen::MatrixXd & | snapshotsMatrix, | ||
| Eigen::VectorXd & | fieldWeights ) |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 898 of file hyperReduction.templates.H.
| void HyperReduction< SnapshotsLists >::stackSnapshotsBoundary | ( | SnapshotsList | sList, |
| List< Eigen::MatrixXd > & | snapshotsMatrixBoundary, | ||
| List< Eigen::VectorXd > & | fieldWeightsBoundary ) |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 916 of file hyperReduction.templates.H.
|
inline |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 349 of file hyperReduction.H.
| void HyperReduction< SnapshotsLists >::updateNodes | ( | Eigen::SparseMatrix< double > & | P, |
| label & | ind_max, | ||
| Eigen::VectorXd & | mp_not_mask ) |
TODO.
| SnapshotsList | PtrList of volScalarField or volVectorField |
| [in] | sList | The list of snapshots |
Definition at line 818 of file hyperReduction.templates.H.
| Eigen::MatrixXd HyperReduction< SnapshotsLists >::basisMatrix |
Orthonormal basis of HR.
Definition at line 166 of file hyperReduction.H.
| Eigen::VectorXd HyperReduction< SnapshotsLists >::eigenValueseig |
Definition at line 181 of file hyperReduction.H.
| Eigen::SparseMatrix<double> HyperReduction< SnapshotsLists >::field2submesh |
Definition at line 187 of file hyperReduction.H.
| List<unsigned int> HyperReduction< SnapshotsLists >::fieldDims |
Dimensions of the fields.
Definition at line 138 of file hyperReduction.H.
| List<word> HyperReduction< SnapshotsLists >::fieldNames |
Names of the fields.
Definition at line 135 of file hyperReduction.H.
| word HyperReduction< SnapshotsLists >::folderMethod |
Folder for the selected HR method.
Definition at line 159 of file hyperReduction.H.
| word HyperReduction< SnapshotsLists >::folderProblem |
Folder for the HR problem.
Definition at line 156 of file hyperReduction.H.
| Eigen::VectorXi HyperReduction< SnapshotsLists >::initialSeeds |
Definition at line 123 of file hyperReduction.H.
| List<label> HyperReduction< SnapshotsLists >::localNodePoints |
Indices of the local node points in the subMesh.
Definition at line 162 of file hyperReduction.H.
| Eigen::MatrixXd HyperReduction< SnapshotsLists >::MatrixOnline |
Online Matrix.
Definition at line 205 of file hyperReduction.H.
| word HyperReduction< SnapshotsLists >::methodName |
Definition at line 50 of file hyperReduction.H.
| label HyperReduction< SnapshotsLists >::n_boundary_cells |
Definition at line 111 of file hyperReduction.H.
| List<label> HyperReduction< SnapshotsLists >::n_boundary_cells_list |
Definition at line 108 of file hyperReduction.H.
| label HyperReduction< SnapshotsLists >::n_boundary_patches |
Definition at line 105 of file hyperReduction.H.
| label HyperReduction< SnapshotsLists >::n_cells |
Int Number of Cells;.
Definition at line 141 of file hyperReduction.H.
| label HyperReduction< SnapshotsLists >::n_cellsSubfields |
Int Number of Cells in submeshes;.
Definition at line 144 of file hyperReduction.H.
|
staticconstexpr |
Definition at line 55 of file hyperReduction.H.
| label HyperReduction< SnapshotsLists >::n_modes |
The maximum number of modes to be considered.
Definition at line 114 of file hyperReduction.H.
| label HyperReduction< SnapshotsLists >::n_nodes |
The maximum number of modes to be considered.
Definition at line 117 of file hyperReduction.H.
| label HyperReduction< SnapshotsLists >::n_snapshots |
The length of the snapshots lists.
Definition at line 120 of file hyperReduction.H.
| autoPtr<IOList<label > > HyperReduction< SnapshotsLists >::nodePoints |
Nodes in the case of the a nonlinear function.
Definition at line 147 of file hyperReduction.H.
| Eigen::VectorXi HyperReduction< SnapshotsLists >::nodes |
Definition at line 126 of file hyperReduction.H.
| Eigen::VectorXd HyperReduction< SnapshotsLists >::normalizingWeights |
Definition at line 169 of file hyperReduction.H.
| 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 184 of file hyperReduction.H.
| ITHACAparameters* HyperReduction< SnapshotsLists >::para |
Definition at line 102 of file hyperReduction.H.
| Eigen::MatrixXd HyperReduction< SnapshotsLists >::pinvPU |
Definition at line 175 of file hyperReduction.H.
| word HyperReduction< SnapshotsLists >::problemName |
The name of the non-linear function e.g. HR_method/residual.
Definition at line 129 of file hyperReduction.H.
| Eigen::VectorXd HyperReduction< SnapshotsLists >::quadratureWeights |
Quadrature weights. Ordered in the same order of matrix P.
Definition at line 196 of file hyperReduction.H.
| Eigen::MatrixXd HyperReduction< SnapshotsLists >::renormalizedBasisMatrix |
Renormalized basis of HR.
Definition at line 172 of file hyperReduction.H.
| SnapshotsListTuple HyperReduction< SnapshotsLists >::snapshotsListTuple |
The snapshots matrix containing the nonlinear function or operator.
Definition at line 132 of file hyperReduction.H.
| autoPtr<fvMeshSubset> HyperReduction< SnapshotsLists >::submesh |
Submesh of the HyperReduction method.
Definition at line 199 of file hyperReduction.H.
| Eigen::SparseMatrix<double> HyperReduction< SnapshotsLists >::submesh2nodes |
Definition at line 190 of file hyperReduction.H.
| Eigen::VectorXi HyperReduction< SnapshotsLists >::submesh2nodesMask |
Definition at line 193 of file hyperReduction.H.
| autoPtr<volVectorField> HyperReduction< SnapshotsLists >::submesh_field |
Submeshes.
Definition at line 202 of file hyperReduction.H.
| unsigned int HyperReduction< SnapshotsLists >::sumFieldsDim |
Definition at line 59 of file hyperReduction.H.
| autoPtr<IOList<labelList > > HyperReduction< SnapshotsLists >::totalNodePoints |
List of label lists of the nodes and corresponding surrounding nodes.
Definition at line 150 of file hyperReduction.H.
| autoPtr<IOList<label > > HyperReduction< SnapshotsLists >::uniqueNodePoints |
List of the unique indices of the nodes that define the submesh.
Definition at line 153 of file hyperReduction.H.
| unsigned int HyperReduction< SnapshotsLists >::vectorial_dim |
Definition at line 58 of file hyperReduction.H.
| Eigen::MatrixXd HyperReduction< SnapshotsLists >::wPU |
Definition at line 178 of file hyperReduction.H.
1.13.2