1 #ifndef __ADAPTIVE_SAMPLING_BASE__ 2 #define __ADAPTIVE_SAMPLING_BASE__ 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" 12 using DealiiVector = dealii::LinearAlgebra::distributed::Vector<double>;
13 using Eigen::MatrixXd;
14 using Eigen::RowVectorXd;
15 using Eigen::VectorXd;
24 template <
int dim,
int nstate>
30 const dealii::ParameterHandler ¶meter_handler_input);
44 mutable std::vector<std::unique_ptr<ProperOrthogonalDecomposition::ROMTestLocation<dim,nstate>>>
rom_locations;
47 mutable std::vector<dealii::LinearAlgebra::distributed::Vector<double>>
fom_locations;
53 std::shared_ptr<ProperOrthogonalDecomposition::OnlinePOD<dim>>
current_pod;
64 dealii::ConditionalOStream
pcout;
82 dealii::LinearAlgebra::distributed::Vector<double>
solveSnapshotFOM(
const RowVectorXd& parameter)
const;
double max_error
Maximum error.
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 ¶meter) const
Solve full-order snapshot.
Files for the baseline physics.
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 ¶meter) 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.
const int mpi_rank
MPI rank.
virtual int run_sampling() const =0
Run Sampling Procedure.
AdaptiveSamplingBase(const PHiLiP::Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler ¶meter_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.