[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
hyper_reduced_adaptive_sampling.h
1 #ifndef __HYPER_REDUCED_ADAPTIVE_SAMPLING__
2 #define __HYPER_REDUCED_ADAPTIVE_SAMPLING__
3 
4 #include <deal.II/numerics/vector_tools.h>
5 #include "parameters/all_parameters.h"
6 #include "pod_basis_online.h"
7 #include "rom_test_location.h"
8 #include <eigen/Eigen/Dense>
9 #include "nearest_neighbors.h"
10 #include "adaptive_sampling_base.h"
11 
12 namespace PHiLiP {
13 
14 using DealiiVector = dealii::LinearAlgebra::distributed::Vector<double>;
15 using Eigen::MatrixXd;
16 using Eigen::RowVectorXd;
17 using Eigen::VectorXd;
18 
20 
21 /*
22 Based on the work in Donovan Blais' thesis:
23 Goal-Oriented Adaptive Sampling for Projection-Based Reduced-Order Models, 2022
24 
25 and the ECSW hyperreduction technique:
26 "Mesh sampling and weighting for the hyperreduction of nonlinear Petrov–Galerkin reduced-order models with local reduced-order bases"
27 Sebastian Grimberg, Charbel Farhat, Radek Tezaur, Charbel Bou-Mosleh
28 International Journal for Numerical Methods in Engineering, 2020
29 https://onlinelibrary.wiley.com/doi/10.1002/nme.6603
30 */
31 template <int dim, int nstate>
33 {
34 public:
37  const dealii::ParameterHandler &parameter_handler_input);
38 
40  int run_sampling () const override;
41 
43  bool placeROMLocations(const MatrixXd& rom_points, Epetra_Vector weights) const;
44 
46  void trueErrorROM(const MatrixXd& rom_points, Epetra_Vector weights) const;
47 
49  double solveSnapshotROMandFOM(const RowVectorXd& parameter, Epetra_Vector weights) const;
50 
52  void solveFunctionalHROM(const RowVectorXd& parameter, Epetra_Vector weights) const;
53 
55  void updateNearestExistingROMs(const RowVectorXd& parameter, Epetra_Vector weights) const;
56 
58  std::unique_ptr<ProperOrthogonalDecomposition::ROMSolution<dim,nstate>> solveSnapshotROM(const RowVectorXd& parameter, Epetra_Vector weights) const;
59 
61  Epetra_Vector allocateVectorToSingleCore(const Epetra_Vector &b) const;
62 
64  mutable std::shared_ptr<Epetra_Vector> ptr_weights;
65 
67  mutable std::vector<double> rom_functional;
68 };
69 
70 }
71 
72 
73 #endif
void solveFunctionalHROM(const RowVectorXd &parameter, Epetra_Vector weights) const
Solve ROM and track functional.
std::vector< double > rom_functional
Functional value predicted by the rom at each sammpling iteration at parameter location specified in ...
Files for the baseline physics.
Definition: ADTypes.hpp:10
double solveSnapshotROMandFOM(const RowVectorXd &parameter, Epetra_Vector weights) const
Solve FOM and ROM, return error in functional between the models.
Main parameter class that contains the various other sub-parameter classes.
HyperreducedAdaptiveSampling(const PHiLiP::Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input)
Constructor.
bool placeROMLocations(const MatrixXd &rom_points, Epetra_Vector weights) const
Placement of ROMs.
std::unique_ptr< ProperOrthogonalDecomposition::ROMSolution< dim, nstate > > solveSnapshotROM(const RowVectorXd &parameter, Epetra_Vector weights) const
Solve reduced-order solution.
Epetra_Vector allocateVectorToSingleCore(const Epetra_Vector &b) const
Copy all elements in matrix A to all cores.
std::shared_ptr< Epetra_Vector > ptr_weights
Ptr vector of ECSW Weights.
void trueErrorROM(const MatrixXd &rom_points, Epetra_Vector weights) const
Compute true/actual error at all ROM points (error in functional between FOM and ROM solution) ...
void updateNearestExistingROMs(const RowVectorXd &parameter, Epetra_Vector weights) const
Updates nearest ROM points to snapshot if error discrepancy is above tolerance.