[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
adaptive_sampling_base.h
1 #ifndef __ADAPTIVE_SAMPLING_BASE__
2 #define __ADAPTIVE_SAMPLING_BASE__
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 
11 namespace PHiLiP {
12 using DealiiVector = dealii::LinearAlgebra::distributed::Vector<double>;
13 using Eigen::MatrixXd;
14 using Eigen::RowVectorXd;
15 using Eigen::VectorXd;
16 
19 
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>
26 {
27 public:
29  AdaptiveSamplingBase(const PHiLiP::Parameters::AllParameters *const parameters_input,
30  const dealii::ParameterHandler &parameter_handler_input);
31 
33  virtual ~AdaptiveSamplingBase() = default;
34 
36 
38  const dealii::ParameterHandler &parameter_handler;
39 
41  mutable MatrixXd snapshot_parameters;
42 
44  mutable std::vector<std::unique_ptr<ProperOrthogonalDecomposition::ROMTestLocation<dim,nstate>>> rom_locations;
45 
47  mutable std::vector<dealii::LinearAlgebra::distributed::Vector<double>> fom_locations;
48 
50  mutable double max_error;
51 
53  std::shared_ptr<ProperOrthogonalDecomposition::OnlinePOD<dim>> current_pod;
54 
56  std::shared_ptr<ProperOrthogonalDecomposition::NearestNeighbors> nearest_neighbors;
57 
58  const MPI_Comm mpi_communicator;
59  const int mpi_rank;
60 
62 
64  dealii::ConditionalOStream pcout;
65 
67  virtual void outputIterationData(std::string iteration) const;
68 
70  RowVectorXd readROMFunctionalPoint() const;
71 
73  virtual int run_sampling () const = 0;
74 
76  void placeInitialSnapshots() const;
77 
79  virtual RowVectorXd getMaxErrorROM() const;
80 
82  dealii::LinearAlgebra::distributed::Vector<double> solveSnapshotFOM(const RowVectorXd& parameter) const;
83 
85  Parameters::AllParameters reinit_params(const RowVectorXd& parameter) const;
86 
88  void configureInitialParameterSpace() const;
89 
90 };
91 
92 }
93 
94 
95 #endif
const dealii::ParameterHandler & parameter_handler
Parameter handler for storing the .prm file being ran.
dealii::ConditionalOStream pcout
ConditionalOStream.
RowVectorXd readROMFunctionalPoint() const
Find point to solve for functional from param file.
std::shared_ptr< ProperOrthogonalDecomposition::NearestNeighbors > nearest_neighbors
Nearest neighbors of snapshots.
dealii::LinearAlgebra::distributed::Vector< double > solveSnapshotFOM(const RowVectorXd &parameter) const
Solve full-order snapshot.
Files for the baseline physics.
Definition: ADTypes.hpp:10
void configureInitialParameterSpace() const
Set up parameter space depending on test case.
Main parameter class that contains the various other sub-parameter classes.
Parameters::AllParameters reinit_params(const RowVectorXd &parameter) const
Reinitialize parameters.
std::shared_ptr< ProperOrthogonalDecomposition::OnlinePOD< dim > > current_pod
Most up to date POD basis.
const MPI_Comm mpi_communicator
MPI communicator.
virtual void outputIterationData(std::string iteration) const
Output for each iteration.
void placeInitialSnapshots() const
Placement of initial snapshots.
virtual int run_sampling() const =0
Run Sampling Procedure.
AdaptiveSamplingBase(const PHiLiP::Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input)
Default constructor that will set the constants.
std::vector< dealii::LinearAlgebra::distributed::Vector< double > > fom_locations
Vector of parameter-ROMTestLocation pairs.
MatrixXd snapshot_parameters
Matrix of snapshot parameters.
const Parameters::AllParameters *const all_parameters
Pointer to all parameters.
virtual ~AdaptiveSamplingBase()=default
Virtual destructor.
std::vector< std::unique_ptr< ProperOrthogonalDecomposition::ROMTestLocation< dim, nstate > > > rom_locations
Vector of parameter-ROMTestLocation pairs.
virtual RowVectorXd getMaxErrorROM() const
Compute RBF and find max error.