[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
hrom_test_location.h
1 #ifndef __HROM_TEST_LOCATION__
2 #define __HROM_TEST_LOCATION__
3 
4 #include "parameters/all_parameters.h"
5 #include "pod_basis_base.h"
6 #include "reduced_order_solution.h"
7 #include <eigen/Eigen/Dense>
8 #include "test_location_base.h"
9 
10 namespace PHiLiP {
11 namespace ProperOrthogonalDecomposition {
12 using Eigen::RowVectorXd;
13 
15 // Based very closely on the ROMTestLocation class
16 
17 /*
18 Based on the work in Donovan Blais' thesis:
19 Goal-Oriented Adaptive Sampling for Projection-Based Reduced-Order Models, 2022
20 
21 Details on the ROM points/errors can be found in sections 5 and 6
22 
23 Derivation of the new error indicator will likely be detailed in Calista Biondic's thesis
24 */
25 template <int dim, int nstate>
26 class HROMTestLocation: public TestLocationBase<dim,nstate>
27 {
28 public:
30  HROMTestLocation(const RowVectorXd& parameter, std::unique_ptr<ROMSolution < dim, nstate>> rom_solution, std::shared_ptr< DGBase<dim, double> > dg_input, Epetra_Vector weights);
31 
34 
36  std::shared_ptr<Epetra_CrsMatrix> generate_hyper_reduced_jacobian(const Epetra_CrsMatrix &system_matrix);
37 
39  std::shared_ptr<Epetra_CrsMatrix> generate_test_basis(Epetra_CrsMatrix &epetra_system_matrix, const Epetra_CrsMatrix &pod_basis);
40 
42  std::shared_ptr<DGBase<dim,double>> dg;
43 
45  Epetra_Vector ECSW_weights;
46 
47 };
48 
49 }
50 }
51 
52 
53 #endif
Base class for a ROM/HROM point, differences would be in the second DWR error indicator.
std::shared_ptr< DGBase< dim, double > > dg
Smart pointer to DGBase.
Files for the baseline physics.
Definition: ADTypes.hpp:10
std::shared_ptr< Epetra_CrsMatrix > generate_hyper_reduced_jacobian(const Epetra_CrsMatrix &system_matrix)
Generate hyper-reduced jacobian matrix.
Class to compute and store adjoint-based error estimates with hyperreduction.
Epetra_Vector ECSW_weights
ECSW hyper-reduction weights.
HROMTestLocation(const RowVectorXd &parameter, std::unique_ptr< ROMSolution< dim, nstate >> rom_solution, std::shared_ptr< DGBase< dim, double > > dg_input, Epetra_Vector weights)
Constructor.
std::unique_ptr< ROMSolution< dim, nstate > > rom_solution
ROM solution.
void compute_initial_rom_to_final_rom_error(std::shared_ptr< ProperOrthogonalDecomposition::PODBase< dim >> pod_updated) override
Compute error between initial ROM and final ROM.
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
Class to hold information about the reduced-order solution.
std::shared_ptr< Epetra_CrsMatrix > generate_test_basis(Epetra_CrsMatrix &epetra_system_matrix, const Epetra_CrsMatrix &pod_basis)
Generate test basis.