[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
pod_adaptive_sampling.h
1 #ifndef __POD_ADAPTIVE_SAMPLING__
2 #define __POD_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 Based on the work in Donovan Blais' thesis:
22 Goal-Oriented Adaptive Sampling for Projection-Based Reduced-Order Models, 2022
23 */
24 template <int dim, int nstate>
25 class AdaptiveSampling: public AdaptiveSamplingBase<dim,nstate>
26 {
27 public:
29  AdaptiveSampling(const PHiLiP::Parameters::AllParameters *const parameters_input,
30  const dealii::ParameterHandler &parameter_handler_input);
31 
34 
36  int run_sampling () const override;
37 
39  bool placeROMLocations(const MatrixXd& rom_points) const;
40 
42  void trueErrorROM(const MatrixXd& rom_points) const;
43 
45  double solveSnapshotROMandFOM(const RowVectorXd& parameter) const;
46 
48  void solveFunctionalROM(const RowVectorXd& parameter) const;
49 
51  void updateNearestExistingROMs(const RowVectorXd& parameter) const;
52 
54  std::unique_ptr<ProperOrthogonalDecomposition::ROMSolution<dim,nstate>> solveSnapshotROM(const RowVectorXd& parameter) const;
55 
57  mutable std::vector<double> rom_functional;
58 
59 };
60 
61 }
62 
63 
64 #endif
void trueErrorROM(const MatrixXd &rom_points) const
Compute true/actual error at all ROM points (error in functional between FOM and ROM solution) ...
void updateNearestExistingROMs(const RowVectorXd &parameter) const
Updates nearest ROM points to snapshot if error discrepancy is above tolerance.
std::unique_ptr< ProperOrthogonalDecomposition::ROMSolution< dim, nstate > > solveSnapshotROM(const RowVectorXd &parameter) const
Solve reduced-order solution.
POD adaptive sampling.
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
bool placeROMLocations(const MatrixXd &rom_points) const
Placement of ROMs.
Main parameter class that contains the various other sub-parameter classes.
AdaptiveSampling(const PHiLiP::Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input)
Constructor.
double solveSnapshotROMandFOM(const RowVectorXd &parameter) const
Solve FOM and ROM, return error in functional between the models.
int run_sampling() const override
Run Sampling Procedure.
void solveFunctionalROM(const RowVectorXd &parameter) const
Solve ROM and track functional.