[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_sampling_error_updated.h
1 #ifndef __HYPER_REDUCED_SAMPLING_ERROR_UPDATED__
2 #define __HYPER_REDUCED_SAMPLING_ERROR_UPDATED__
3 
4 #include <deal.II/numerics/vector_tools.h>
5 #include "parameters/all_parameters.h"
6 #include "pod_basis_online.h"
7 #include "hrom_test_location.h"
8 #include <eigen/Eigen/Dense>
9 #include "nearest_neighbors.h"
10 #include "adaptive_sampling_base.h"
11 
12 namespace PHiLiP {
13 using DealiiVector = dealii::LinearAlgebra::distributed::Vector<double>;
14 using Eigen::MatrixXd;
15 using Eigen::RowVectorXd;
16 using Eigen::VectorXd;
17 
19 // Currently seperate from the adaptive sampling base as the pointer of ROM locations has
20 // been changed to a have type HROMTestLocation to update DWR errors
21 
22 /*
23 Based on the work in Donovan Blais' thesis:
24 Goal-Oriented Adaptive Sampling for Projection-Based Reduced-Order Models, 2022
25 
26 Details on the ROM points/errors can be found in sections 5 and 6
27 
28 Derivation of the new error indicator will likely be detailed in Calista Biondic's thesis
29 */
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  RowVectorXd getMaxErrorROM() const override;
44 
46  bool placeROMLocations(const MatrixXd& rom_points, Epetra_Vector weights) const;
47 
49  void trueErrorROM(const MatrixXd& rom_points, Epetra_Vector weights) const;
50 
52  double solveSnapshotROMandFOM(const RowVectorXd& parameter, Epetra_Vector weights) const;
53 
55  void solveFunctionalHROM(const RowVectorXd& parameter, Epetra_Vector weights) const;
56 
58  void updateNearestExistingROMs(const RowVectorXd& parameter, Epetra_Vector weights) const;
59 
61  std::unique_ptr<ProperOrthogonalDecomposition::ROMSolution<dim,nstate>> solveSnapshotROM(const RowVectorXd& parameter, Epetra_Vector weights) const;
62 
64  Epetra_Vector allocateVectorToSingleCore(const Epetra_Vector &b) const;
65 
67  void outputIterationData(std::string iteration) const override;
68 
70  mutable std::vector<std::unique_ptr<ProperOrthogonalDecomposition::HROMTestLocation<dim,nstate>>> hrom_locations;
71 
73  mutable std::shared_ptr<Epetra_Vector> ptr_weights;
74 
76  mutable std::vector<double> rom_functional;
77 
78 };
79 
80 }
81 
82 
83 #endif
std::unique_ptr< ProperOrthogonalDecomposition::ROMSolution< dim, nstate > > solveSnapshotROM(const RowVectorXd &parameter, Epetra_Vector weights) const
Solve reduced-order solution.
double solveSnapshotROMandFOM(const RowVectorXd &parameter, Epetra_Vector weights) const
Solve FOM and ROM, return error in functional between the models.
HyperreducedSamplingErrorUpdated(const PHiLiP::Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input)
Constructor.
Files for the baseline physics.
Definition: ADTypes.hpp:10
std::vector< double > rom_functional
Functional value predicted by the rom at each sammpling iteration at parameter location specified in ...
Main parameter class that contains the various other sub-parameter classes.
Hyperreduced Adaptive Sampling with the updated error indicator.
Epetra_Vector allocateVectorToSingleCore(const Epetra_Vector &b) const
Copy all elements in matrix A to all cores.
bool placeROMLocations(const MatrixXd &rom_points, Epetra_Vector weights) const
Placement of ROMs.
std::vector< std::unique_ptr< ProperOrthogonalDecomposition::HROMTestLocation< dim, nstate > > > hrom_locations
Vector of parameter-HROMTestLocation pairs.
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) ...
RowVectorXd getMaxErrorROM() const override
Compute RBF and find max error.
void updateNearestExistingROMs(const RowVectorXd &parameter, Epetra_Vector weights) const
Updates nearest ROM points to snapshot if error discrepancy is above tolerance.
void outputIterationData(std::string iteration) const override
Output for each iteration.
void solveFunctionalHROM(const RowVectorXd &parameter, Epetra_Vector weights) const
Solve HROM and track functional.
std::shared_ptr< Epetra_Vector > ptr_weights
Ptr vector of ECSW Weights.