[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
nearest_neighbors.h
1 #ifndef __NEAREST_NEIGHBORS__
2 #define __NEAREST_NEIGHBORS__
3 
4 #include <eigen/Eigen/Dense>
5 #include <deal.II/lac/la_parallel_vector.h>
6 #include "min_max_scaler.h"
7 
8 namespace PHiLiP {
9 namespace ProperOrthogonalDecomposition {
10 using Eigen::MatrixXd;
11 using Eigen::VectorXd;
12 using Eigen::RowVectorXd;
13 
16 {
17 public:
20 
22  void update_snapshots(const MatrixXd& snapshot_parameters, dealii::LinearAlgebra::distributed::Vector<double> snapshot);
23 
26 
28  MatrixXd kNearestNeighborsMidpoint(const RowVectorXd& point);
29 
30  //Given a point, return the index of nearest snapshot_parameter
31  dealii::LinearAlgebra::distributed::Vector<double> nearestNeighborMidpointSolution(const RowVectorXd& point);
32 
34  MatrixXd snapshot_params;
35 
38 
41 
43  std::vector<dealii::LinearAlgebra::distributed::Vector<double>> snapshots;
44 
45 };
46 
47 }
48 }
49 
50 #endif
void update_snapshots(const MatrixXd &snapshot_parameters, dealii::LinearAlgebra::distributed::Vector< double > snapshot)
Add snapshot.
Files for the baseline physics.
Definition: ADTypes.hpp:10
MatrixXd scaled_snapshot_params
Scaled snapshot parameters.
std::vector< dealii::LinearAlgebra::distributed::Vector< double > > snapshots
Vector containing all snapshots.
MatrixXd kPairwiseNearestNeighborsMidpoint()
Find midpoint of all snapshot locations.
MatrixXd kNearestNeighborsMidpoint(const RowVectorXd &point)
Given a point, returns midpoint between point and k nearest snapshots, where k is 1+num_parameters...