1 #ifndef __HYPER_REDUCED_ADAPTIVE_SAMPLING__     2 #define __HYPER_REDUCED_ADAPTIVE_SAMPLING__     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 #include "adaptive_sampling_base.h"    14 using DealiiVector = dealii::LinearAlgebra::distributed::Vector<double>;
    15 using Eigen::MatrixXd;
    16 using Eigen::RowVectorXd;
    17 using Eigen::VectorXd;
    31 template <
int dim, 
int nstate>
    37                      const dealii::ParameterHandler ¶meter_handler_input);
    46     void trueErrorROM(
const MatrixXd& rom_points, Epetra_Vector weights) 
const;
    58     std::unique_ptr<ProperOrthogonalDecomposition::ROMSolution<dim,nstate>> 
solveSnapshotROM(
const RowVectorXd& parameter, Epetra_Vector weights) 
const;
 
void solveFunctionalHROM(const RowVectorXd ¶meter, Epetra_Vector weights) const
Solve ROM and track functional. 
std::vector< double > rom_functional
Functional value predicted by the rom at each sammpling iteration at parameter location specified in ...
Files for the baseline physics. 
double solveSnapshotROMandFOM(const RowVectorXd ¶meter, Epetra_Vector weights) const
Solve FOM and ROM, return error in functional between the models. 
Main parameter class that contains the various other sub-parameter classes. 
HyperreducedAdaptiveSampling(const PHiLiP::Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler ¶meter_handler_input)
Constructor. 
bool placeROMLocations(const MatrixXd &rom_points, Epetra_Vector weights) const
Placement of ROMs. 
std::unique_ptr< ProperOrthogonalDecomposition::ROMSolution< dim, nstate > > solveSnapshotROM(const RowVectorXd ¶meter, Epetra_Vector weights) const
Solve reduced-order solution. 
int run_sampling() const override
Run test. 
Epetra_Vector allocateVectorToSingleCore(const Epetra_Vector &b) const
Copy all elements in matrix A to all cores. 
Hyperreduced adaptive sampling. 
std::shared_ptr< Epetra_Vector > ptr_weights
Ptr vector of ECSW Weights. 
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) ...
void updateNearestExistingROMs(const RowVectorXd ¶meter, Epetra_Vector weights) const
Updates nearest ROM points to snapshot if error discrepancy is above tolerance.