[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
ROM_error_post_sampling.cpp
1 #include "ROM_error_post_sampling.h"
2 #include "reduced_order/pod_basis_offline.h"
3 #include "parameters/all_parameters.h"
4 #include "flow_solver/flow_solver.h"
5 #include "flow_solver/flow_solver_factory.h"
6 #include "ode_solver/ode_solver_factory.h"
7 #include "reduced_order/assemble_ECSW_residual.h"
8 #include "reduced_order/assemble_ECSW_jacobian.h"
9 #include "linear_solver/NNLS_solver.h"
10 #include "linear_solver/helper_functions.h"
11 #include "reduced_order/pod_adaptive_sampling.h"
12 #include "reduced_order/hyper_reduced_adaptive_sampling.h"
13 #include "rom_import_helper_functions.h"
14 #include <eigen/Eigen/Dense>
15 #include <iostream>
16 #include <filesystem>
17 
18 namespace PHiLiP {
19 namespace Tests {
20 
21 template <int dim, int nstate>
23  const dealii::ParameterHandler &parameter_handler_input)
24  : TestsBase::TestsBase(parameters_input)
25  , parameter_handler(parameter_handler_input)
26 {}
27 
28 template <int dim, int nstate>
30  // Copy all parameters
32 
33  parameters.reduced_order_param.path_to_search = path;
34  return parameters;
35 }
36 
37 template <int dim, int nstate>
39 {
40  pcout << "Starting error analysis for ROM..." << std::endl;
41 
42  // Create POD Petrov-Galerkin ROM from Offline POD Files
43  std::unique_ptr<FlowSolver::FlowSolver<dim,nstate>> flow_solver_petrov_galerkin = FlowSolver::FlowSolverFactory<dim,nstate>::select_flow_case(all_parameters, parameter_handler);
44  std::shared_ptr<ProperOrthogonalDecomposition::OfflinePOD<dim>> pod_petrov_galerkin = std::make_shared<ProperOrthogonalDecomposition::OfflinePOD<dim>>(flow_solver_petrov_galerkin->dg);
45 
46  // Create Instance of Adaptive Sampling to calculate the error between the FOM and ROM at the points from getROMPoints
47  std::shared_ptr<AdaptiveSampling<dim,nstate>> parameter_sampling = std::make_unique<AdaptiveSampling<dim,nstate>>(all_parameters, parameter_handler);
48  parameter_sampling->current_pod->basis = pod_petrov_galerkin->basis;
49  parameter_sampling->current_pod->referenceState = pod_petrov_galerkin->referenceState;
50  parameter_sampling->current_pod->snapshotMatrix = pod_petrov_galerkin->snapshotMatrix;
52  std::string path = all_parameters->reduced_order_param.path_to_search; //Search specified directory for files containing "solutions_table"
53  bool snap_found = getSnapshotParamsFromFile(snapshot_parameters, path);
54  if (snap_found){
55  parameter_sampling->snapshot_parameters = snapshot_parameters;
56  pcout << "snapshot_parameters" << std::endl;
57  pcout << snapshot_parameters << std::endl;
58  }
59  else{
60  pcout << "File with snapshots not found in folder" << std::endl;
61  return -1;
62  }
63  getROMPoints(rom_points, all_parameters);
64  pcout << "ROM Locations" << std::endl;
65  pcout << rom_points << std::endl;
66  parameter_sampling->trueErrorROM(rom_points);
67 
68  return 0;
69 }
70 
71 #if PHILIP_DIM==1
73 #endif
74 
75 #if PHILIP_DIM!=1
77 #endif
78 } // Tests namespace
79 } // PHiLiP namespace
int run_test() const override
Evaluate and output the "true" error at ROM Points.
ROMErrorPostSampling(const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input)
Constructor.
MatrixXd snapshot_parameters
Matrix of snapshot parameters.
Files for the baseline physics.
Definition: ADTypes.hpp:10
ReducedOrderModelParam reduced_order_param
Contains parameters for the Reduced-Order model.
std::string path_to_search
Path to search for snapshots or saved POD basis.
Main parameter class that contains the various other sub-parameter classes.
static std::unique_ptr< FlowSolver< dim, nstate > > select_flow_case(const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input)
Factory to return the correct flow solver given input file.
const Parameters::AllParameters *const all_parameters
Pointer to all parameters.
Definition: tests.h:20
Parameters::AllParameters reinit_params(std::string path) const
Reinitialize parameters.
const dealii::ParameterHandler & parameter_handler
Dummy parameter handler because flowsolver requires it.
MatrixXd rom_points
Matrix of error sampling points.
dealii::ConditionalOStream pcout
ConditionalOStream.
Definition: tests.h:45
Base class of all the tests.
Definition: tests.h:17