BubbleProfiler  0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
field_profiles.hpp
Go to the documentation of this file.
1 /*
2  * This file is part of BubbleProfiler.
3  *
4  * BubbleProfiler is free software: you can redistribute it and/or modify
5  * it under the terms of the GNU General Public License as published by
6  * the Free Software Foundation, either version 3 of the License, or
7  * (at your option) any later version.
8  *
9  * BubbleProfiler is distributed in the hope that it will be useful,
10  * but WITHOUT ANY WARRANTY; without even the implied warranty of
11  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  * GNU General Public License for more details.
13  *
14  * You should have received a copy of the GNU General Public License
15  * along with BubbleProfiler. If not, see <http://www.gnu.org/licenses/>.
16  */
17 
18 #ifndef BUBBLEPROFILER_FIELD_PROFILES_HPP_INCLUDED
19 #define BUBBLEPROFILER_FIELD_PROFILES_HPP_INCLUDED
20 
26 #include "gsl_interpolator.hpp"
27 
28 #include <Eigen/Core>
29 
30 #include <vector>
31 
32 namespace BubbleProfiler {
33 
50 public:
51  Field_profiles() = default;
52 
60  Field_profiles(const Eigen::VectorXd& gridpoints_,
61  const Eigen::MatrixXd& profiles_,
62  double interpolation_points_fraction_ = 1.0);
63 
73  Field_profiles(int num_fields_, double domain_start_,
74  double domain_end_, int num_steps_,
75  double interpolation_points_fraction_ = 1.0);
76 
85  Field_profiles(const Eigen::MatrixXd& profiles_,
86  double domain_start_, double domain_end_,
87  double interpolation_points_fraction_ = 1.0);
88 
89  ~Field_profiles() = default;
90  Field_profiles(const Field_profiles&) = default;
91  Field_profiles& operator=(const Field_profiles&) = default;
92  Field_profiles(Field_profiles&&) = default;
94 
96  void set_interpolation_points_fraction(double f);
97  void set_field_profile(int i, const Eigen::VectorXd& p);
98 
99  int get_number_of_fields() const { return profiles.cols(); }
100  int get_number_of_dimensions() const { return n_spatial_dims; }
101 
102  int get_number_of_grid_points() const { return gridpoints.size(); }
103  double get_domain_start() const { return gridpoints(0); }
104  double get_domain_end() const { return gridpoints(gridpoints.size() - 1); }
105 
111  const Eigen::VectorXd& get_spatial_grid() const { return gridpoints; }
112 
117  const Eigen::MatrixXd& get_field_profiles() const { return profiles; }
118 
124  Eigen::VectorXd evaluate_at(double rho) const;
125 
132  double evaluate_at(int field, double rho) const;
133 
141  double derivative_at(int field, int order, double rho) const;
142 
143 private:
144  Eigen::VectorXd gridpoints{};
145  Eigen::MatrixXd profiles{};
146  double n_spatial_dims{3};
147  bool must_update_coords{true};
149  Eigen::VectorXd interpolation_grid_values{};
150  Eigen::MatrixXd interpolation_field_values{};
151  std::vector<GSL_interpolator> splines{};
152 
156  void build_spline_for_field(int);
157  void build_splines();
158 };
159 
160 } // namespace BubbleProfiler
161 
162 #endif
double derivative_at(int field, int order, double rho) const
Take the radial derivative of a field at a given radius.
const Eigen::VectorXd & get_spatial_grid() const
Get a vector of the grid point coordinates.
void set_interpolation_points_fraction(double f)
void set_field_profile(int i, const Eigen::VectorXd &p)
std::vector< GSL_interpolator > splines
Eigen::MatrixXd interpolation_field_values
Eigen::VectorXd evaluate_at(double rho) const
Get all field values at a given radius.
contains the definition of the GSL_interpolator class
Eigen::VectorXd interpolation_grid_values
Field_profiles & operator=(const Field_profiles &)=default
Discretized set of field profiles.
const Eigen::MatrixXd & get_field_profiles() const
Get the field profile data in matrix form.