libROM  v1.0
Data-driven physical simulation library
Interpolator.h
1 
11 // Description: Computes the Interpolator algorithm on the given snapshot matrix.
12 
13 #ifndef included_Interpolator_h
14 #define included_Interpolator_h
15 
16 #include <vector>
17 #include <string>
18 #include <memory>
19 
20 namespace CAROM {
21 
22 class Matrix;
23 class Vector;
24 
33 {
34 protected:
35 
54  Interpolator(const std::vector<Vector> & parameter_points,
55  const std::vector<std::shared_ptr<Matrix>> & rotation_matrices,
56  int ref_point,
57  std::string rbf,
58  std::string interp_method,
59  double closest_rbf_val = 0.9,
60  bool compute_gradients = false);
61 
65  int d_rank;
66 
71 
77 
82  std::string d_rbf;
83 
88  std::string d_interp_method;
89 
95  double d_epsilon;
96 
100  std::vector<Vector> d_parameter_points;
101 
105  std::vector<std::shared_ptr<Matrix>> d_rotation_matrices;
106 
110  std::unique_ptr<Matrix> d_lambda_T;
111 
117 
122  std::vector<std::shared_ptr<Vector>> d_interpolation_gradient;
123 
124 private:
125 
129  Interpolator();
130 
134  Interpolator(
135  const Interpolator& other);
136 
140  Interpolator&
141  operator = (
142  const Interpolator& rhs);
143 };
144 
159 std::vector<double> obtainRBFToTrainingPoints(
160  const std::vector<Vector> & parameter_points,
161  const std::string & interp_method, const std::string & rbf, double epsilon,
162  const Vector & point);
163 
180 std::vector<double> obtainRBFGradientToTrainingPoints(
181  const std::vector<Vector> & parameter_points,
182  const std::string & interp_method, const std::string & rbf, double epsilon,
183  const Vector & point, const int index);
184 
185 
197 double obtainRBFGradient(std::string rbf, double epsilon, const Vector & point1,
198  const Vector & point2, const int index);
199 
205 double rbfWeightedSum(std::vector<double>& rbf);
206 
216 double obtainRBF(std::string rbf, double epsilon, const Vector & point1,
217  const Vector & point2);
218 
227 double convertClosestRBFToEpsilon(const std::vector<Vector> & parameter_points,
228  const std::string & rbf, double closest_rbf_val);
229 
240 std::vector<std::shared_ptr<Matrix>> obtainRotationMatrices(
241  const std::vector<Vector> & parameter_points,
242  const std::vector<std::shared_ptr<Matrix>> & bases,
243  int ref_point);
244 
245 std::unique_ptr<std::vector<Vector>>
246  scalarsToVectors(const std::vector<double> & s);
247 }
248 #endif
std::vector< std::shared_ptr< Vector > > d_interpolation_gradient
Gradient with respect to the parameters. Only exists after interpolate has been run.
Definition: Interpolator.h:122
std::vector< std::shared_ptr< Matrix > > d_rotation_matrices
The reduced bases with compatible coordinates.
Definition: Interpolator.h:105
double d_epsilon
The RBF parameter that determines the width of influence. a small epsilon: larger influential width a...
Definition: Interpolator.h:95
int d_ref_point
The index within the vector of parameter points to the reference point.
Definition: Interpolator.h:76
int d_rank
The rank of the process this object belongs to.
Definition: Interpolator.h:65
std::unique_ptr< Matrix > d_lambda_T
The RHS of the linear solve in tangential space.
Definition: Interpolator.h:110
std::vector< Vector > d_parameter_points
The sampled parameter points.
Definition: Interpolator.h:100
int d_num_procs
The number of processors being run on.
Definition: Interpolator.h:70
bool d_compute_gradients
Flag that determines if a gradient with respect to the parameters should be computed.
Definition: Interpolator.h:116
std::string d_rbf
The RBF type (gaussian, multiquadric, inverse quadratic, inverse multiquadric)
Definition: Interpolator.h:82
std::string d_interp_method
The interpolation method (linear solve, inverse distance weighting, lagrangian polynomials)
Definition: Interpolator.h:88
double rbfWeightedSum(std::vector< double > &rbf)
Compute the sum of the RBF weights.
std::vector< std::shared_ptr< Matrix > > obtainRotationMatrices(const std::vector< Vector > &parameter_points, const std::vector< std::shared_ptr< Matrix >> &bases, int ref_point)
Obtain the rotation matrices for all the parameter points using the basis of the reference point.
std::vector< double > obtainRBFGradientToTrainingPoints(const std::vector< Vector > &parameter_points, const std::string &interp_method, const std::string &rbf, double epsilon, const Vector &point, const int index)
Compute the RBF gradient from the parameter points with the unsampled parameter point.
std::vector< double > obtainRBFToTrainingPoints(const std::vector< Vector > &parameter_points, const std::string &interp_method, const std::string &rbf, double epsilon, const Vector &point)
Compute the RBF from the parameter points with the unsampled parameter point.
double convertClosestRBFToEpsilon(const std::vector< Vector > &parameter_points, const std::string &rbf, double closest_rbf_val)
Convert closest RBF value to an epsilon value.
double obtainRBF(std::string rbf, double epsilon, const Vector &point1, const Vector &point2)
Compute the RBF between two points.
double obtainRBFGradient(std::string rbf, double epsilon, const Vector &point1, const Vector &point2, const int index)
Compute the RBF gradient between two points.