[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
HROM_error_post_sampling.h
1 #ifndef __HROM_ERROR_POST_SAMPLING_H__
2 #define __HROM_ERROR_POST_SAMPLING_H__
3 
4 #include "tests.h"
5 #include "parameters/all_parameters.h"
6 #include "dg/dg_base.hpp"
7 #include <eigen/Eigen/Dense>
8 #include <Epetra_CrsMatrix.h>
9 #include <Epetra_Map.h>
10 #include <Epetra_Vector.h>
11 
12 namespace PHiLiP {
13 namespace Tests {
14 
15 using Eigen::MatrixXd;
16 using Eigen::VectorXd;
17 
24 template <int dim, int nstate>
26 {
27 public:
29  HROMErrorPostSampling(const Parameters::AllParameters *const parameters_input,
30  const dealii::ParameterHandler &parameter_handler_input);
31 
33  Parameters::AllParameters reinit_params(std::string path) const;
34 
36  bool getWeightsFromFile(std::shared_ptr<DGBase<dim,double>> &dg) const;
37 
39  int run_test () const override;
40 
42  const dealii::ParameterHandler &parameter_handler;
43 
45  mutable MatrixXd snapshot_parameters;
46 
48  mutable MatrixXd rom_points;
49 
51  mutable std::shared_ptr<Epetra_Vector> ptr_weights;
52 };
53 } // End of Tests namespace
54 } // End of PHiLiP namespace
55 
56 #endif
HROMErrorPostSampling(const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input)
Constructor.
Files for the baseline physics.
Definition: ADTypes.hpp:10
Main parameter class that contains the various other sub-parameter classes.
const dealii::ParameterHandler & parameter_handler
Dummy parameter handler because flowsolver requires it.
MatrixXd snapshot_parameters
Matrix of snapshot parameters.
MatrixXd rom_points
Matrix of error sampling points.
bool getWeightsFromFile(std::shared_ptr< DGBase< dim, double >> &dg) const
Read ECSW weights from the text file.
int run_test() const override
Evaluate and output the "true" error at ROM Points.
Parameters::AllParameters reinit_params(std::string path) const
Reinitialize parameters.
std::shared_ptr< Epetra_Vector > ptr_weights
Ptr vector of ECSW Weights.
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
Base class of all the tests.
Definition: tests.h:17