BubbleProfiler
0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
|
Discretized set of field profiles. More...
#include <field_profiles.hpp>
Public Member Functions | |
Field_profiles ()=default | |
Field_profiles (const Eigen::VectorXd &gridpoints_, const Eigen::MatrixXd &profiles_, double interpolation_points_fraction_=1.0) | |
Create field profiles at given coordinates from data matrix. More... | |
Field_profiles (int num_fields_, double domain_start_, double domain_end_, int num_steps_, double interpolation_points_fraction_=1.0) | |
Create empty field profiles with specified parameters. More... | |
Field_profiles (const Eigen::MatrixXd &profiles_, double domain_start_, double domain_end_, double interpolation_points_fraction_=1.0) | |
Create field profiles from data matrix. More... | |
~Field_profiles ()=default | |
Field_profiles (const Field_profiles &)=default | |
Field_profiles & | operator= (const Field_profiles &)=default |
Field_profiles (Field_profiles &&)=default | |
Field_profiles & | operator= (Field_profiles &&)=default |
void | set_number_of_dimensions (int d) |
void | set_interpolation_points_fraction (double f) |
void | set_field_profile (int i, const Eigen::VectorXd &p) |
int | get_number_of_fields () const |
int | get_number_of_dimensions () const |
int | get_number_of_grid_points () const |
double | get_domain_start () const |
double | get_domain_end () const |
const Eigen::VectorXd & | get_spatial_grid () const |
Get a vector of the grid point coordinates. More... | |
const Eigen::MatrixXd & | get_field_profiles () const |
Get the field profile data in matrix form. More... | |
Eigen::VectorXd | evaluate_at (double rho) const |
Get all field values at a given radius. More... | |
double | evaluate_at (int field, double rho) const |
Get a specific field value at a given radius. More... | |
double | derivative_at (int field, int order, double rho) const |
Take the radial derivative of a field at a given radius. More... | |
Private Member Functions | |
void | initialize_interpolation_grid_values () |
void | initialize_interpolation_field_values (int) |
void | initialize_interpolation_field_values () |
void | build_spline_for_field (int) |
void | build_splines () |
Private Attributes | |
Eigen::VectorXd | gridpoints {} |
Eigen::MatrixXd | profiles {} |
double | n_spatial_dims {3} |
bool | must_update_coords {true} |
double | interpolation_points_fraction {1.0} |
Eigen::VectorXd | interpolation_grid_values {} |
Eigen::MatrixXd | interpolation_field_values {} |
std::vector< GSL_interpolator > | splines {} |
Discretized set of field profiles.
This class is a container for discretized field profiles, representing a solution for or correction to the bubble profile. Methods are provided to interrogate the dimensionality, size and discretization parameters of the profile set. It is also possible to evaluate a point in field space, either directly or via interpolation for off-grid points. Finally, it is possible to add one field profile to another. This is used to implement iterative corrections when calculating a profile solution.
It should be noted that the profile data is assumed to be spherically symmetric.
Definition at line 49 of file field_profiles.hpp.
|
default |
BubbleProfiler::Field_profiles::Field_profiles | ( | const Eigen::VectorXd & | gridpoints_, |
const Eigen::MatrixXd & | profiles_, | ||
double | interpolation_points_fraction_ = 1.0 |
||
) |
Create field profiles at given coordinates from data matrix.
gridpoints_ | Vector containing radial coordinate values |
profiles_ | Matrix containing field data |
interpolation_points_fraction_ | fraction of grid points to use for interpolation |
Definition at line 30 of file field_profiles.cpp.
References build_splines().
BubbleProfiler::Field_profiles::Field_profiles | ( | int | num_fields_, |
double | domain_start_, | ||
double | domain_end_, | ||
int | num_steps_, | ||
double | interpolation_points_fraction_ = 1.0 |
||
) |
Create empty field profiles with specified parameters.
num_fields_ | number of scalar fields |
domain_start_ | domain start (minimum radius) |
domain_end_ | domain end (maximum radius) |
num_steps_ | number of grid points per dimension |
interpolation_points_fraction_ | fraction of grid points to use for interpolation |
Definition at line 46 of file field_profiles.cpp.
BubbleProfiler::Field_profiles::Field_profiles | ( | const Eigen::MatrixXd & | profiles_, |
double | domain_start_, | ||
double | domain_end_, | ||
double | interpolation_points_fraction_ = 1.0 |
||
) |
Create field profiles from data matrix.
profiles_ | matrix containing field data |
domain_start_ | domain start (minimum radius) |
domain_end_ | domain end (maximum radius) |
interpolation_points_fraction_ | fraction of grid points to use for interpolation |
Definition at line 56 of file field_profiles.cpp.
|
default |
|
default |
|
default |
|
private |
Definition at line 143 of file field_profiles.cpp.
References build_splines(), initialize_interpolation_field_values(), interpolation_field_values, interpolation_grid_values, must_update_coords, and splines.
Referenced by set_field_profile().
|
private |
Definition at line 157 of file field_profiles.cpp.
References initialize_interpolation_field_values(), initialize_interpolation_grid_values(), interpolation_field_values, interpolation_grid_values, must_update_coords, gaussian_alpha_tests::n_fields, profiles, and splines.
Referenced by build_spline_for_field(), and Field_profiles().
double BubbleProfiler::Field_profiles::derivative_at | ( | int | field, |
int | order, | ||
double | rho | ||
) | const |
Take the radial derivative of a field at a given radius.
field | target field index |
order | order of the derivative (1 or 2) |
rho | target radius |
Definition at line 197 of file field_profiles.cpp.
References profiles, and splines.
Referenced by BubbleProfiler::calculate_full_action(), BubbleProfiler::Perturbations_ODE_system::calculate_inhomogeneities(), BubbleProfiler::calculate_kinetic_action(), BubbleProfiler::Perturbative_profiler< Integration_policy >::calculate_perturbation(), get_field_profiles(), and BubbleProfiler::Shooting_profile_guesser::get_large_distance_solution().
Eigen::VectorXd BubbleProfiler::Field_profiles::evaluate_at | ( | double | rho | ) | const |
Get all field values at a given radius.
rho | target radius |
Definition at line 186 of file field_profiles.cpp.
References gaussian_alpha_tests::n_fields, profiles, and splines.
Referenced by BubbleProfiler::Shooting_profile_guesser::calculate_field_profiles(), BubbleProfiler::calculate_full_action(), BubbleProfiler::Perturbations_ODE_system::calculate_inhomogeneities(), BubbleProfiler::Perturbations_ODE_system::calculate_mass_matrix(), BubbleProfiler::Perturbative_profiler< Integration_policy >::calculate_perturbation(), BubbleProfiler::calculate_potential_action(), get_field_profiles(), and BubbleProfiler::Relative_convergence_tester::is_converged().
double BubbleProfiler::Field_profiles::evaluate_at | ( | int | field, |
double | rho | ||
) | const |
Get a specific field value at a given radius.
field | target field index |
rho | target radius |
Definition at line 174 of file field_profiles.cpp.
|
inline |
Definition at line 104 of file field_profiles.hpp.
References gridpoints.
Referenced by BubbleProfiler::Perturbative_profiler< Integration_policy >::calculate_bubble_profile(), BubbleProfiler::Shooting_profile_guesser::calculate_field_profiles(), BubbleProfiler::calculate_full_action(), BubbleProfiler::calculate_kinetic_action(), BubbleProfiler::calculate_potential_action(), and BubbleProfiler::Shooting_profile_guesser::get_large_distance_solution().
|
inline |
Definition at line 103 of file field_profiles.hpp.
References gridpoints.
Referenced by BubbleProfiler::Perturbative_profiler< Integration_policy >::calculate_bubble_profile(), BubbleProfiler::Shooting_profile_guesser::calculate_field_profiles(), BubbleProfiler::calculate_full_action(), BubbleProfiler::calculate_kinetic_action(), BubbleProfiler::calculate_potential_action(), and BubbleProfiler::Relative_convergence_tester::is_converged().
|
inline |
Get the field profile data in matrix form.
Definition at line 117 of file field_profiles.hpp.
References derivative_at(), evaluate_at(), and profiles.
Referenced by BubbleProfiler::Perturbative_profiler< Integration_policy >::calculate_bubble_profile(), BubbleProfiler::solve_perturbative(), BubbleProfiler::Perturbative_profiler< Integration_policy >::update_field_profiles(), BubbleProfiler::Thin_wall_observer::write_data(), BubbleProfiler::Logarithmic_observer::write_data(), BubbleProfiler::Generalized_fubini_observer::write_data(), BubbleProfiler::Profile_observer::write_field_profiles_to_file(), BubbleProfiler::write_profile(), and BubbleProfiler::write_profiles().
|
inline |
Definition at line 100 of file field_profiles.hpp.
References n_spatial_dims.
Referenced by BubbleProfiler::calculate_full_action(), BubbleProfiler::calculate_kinetic_action(), and BubbleProfiler::calculate_potential_action().
|
inline |
Definition at line 99 of file field_profiles.hpp.
References profiles.
Referenced by BubbleProfiler::Perturbative_profiler< Integration_policy >::calculate_bubble_profile(), BubbleProfiler::calculate_full_action(), BubbleProfiler::calculate_kinetic_action(), and BubbleProfiler::Perturbative_profiler< Integration_policy >::calculate_perturbation().
|
inline |
Definition at line 102 of file field_profiles.hpp.
References gridpoints.
Referenced by BubbleProfiler::Shooting_profile_guesser::calculate_field_profiles(), BubbleProfiler::solve_perturbative(), and BubbleProfiler::write_profile().
|
inline |
Get a vector of the grid point coordinates.
Definition at line 111 of file field_profiles.hpp.
References gridpoints.
Referenced by BubbleProfiler::Perturbative_profiler< Integration_policy >::calculate_bubble_profile(), BubbleProfiler::Shooting_profile_guesser::calculate_field_profiles(), main(), BubbleProfiler::Generalized_fubini_observer::operator()(), BubbleProfiler::Logarithmic_observer::operator()(), BubbleProfiler::solve_perturbative(), BubbleProfiler::Perturbative_profiler< Integration_policy >::update_field_profiles(), BubbleProfiler::Thin_wall_observer::write_data(), BubbleProfiler::Generalized_fubini_observer::write_data(), BubbleProfiler::Logarithmic_observer::write_data(), BubbleProfiler::Profile_observer::write_field_profiles_to_file(), BubbleProfiler::write_profile(), and BubbleProfiler::write_profiles().
|
private |
Definition at line 111 of file field_profiles.cpp.
References gridpoints, interpolation_field_values, interpolation_grid_values, and profiles.
|
private |
Definition at line 125 of file field_profiles.cpp.
References gridpoints, interpolation_field_values, interpolation_grid_values, and profiles.
Referenced by build_spline_for_field(), and build_splines().
|
private |
Definition at line 87 of file field_profiles.cpp.
References gridpoints, interpolation_grid_values, and interpolation_points_fraction.
Referenced by build_splines().
|
default |
|
default |
void BubbleProfiler::Field_profiles::set_field_profile | ( | int | i, |
const Eigen::VectorXd & | p | ||
) |
Definition at line 76 of file field_profiles.cpp.
References build_spline_for_field(), and profiles.
Referenced by set_number_of_dimensions().
void BubbleProfiler::Field_profiles::set_interpolation_points_fraction | ( | double | f | ) |
Definition at line 65 of file field_profiles.cpp.
References interpolation_points_fraction, and must_update_coords.
Referenced by BubbleProfiler::Shooting_profile_guesser::calculate_field_profiles(), and set_number_of_dimensions().
|
inline |
Definition at line 95 of file field_profiles.hpp.
References n_spatial_dims, set_field_profile(), and set_interpolation_points_fraction().
Referenced by BubbleProfiler::Shooting_profile_guesser::calculate_field_profiles(), BubbleProfiler::Kink_profile_guesser::calculate_field_profiles(), BubbleProfiler::Instream_profile_guesser::get_profile_guess(), and BubbleProfiler::Perturbative_profiler< Integration_policy >::update_field_profiles().
|
private |
Definition at line 144 of file field_profiles.hpp.
Referenced by get_domain_end(), get_domain_start(), get_number_of_grid_points(), get_spatial_grid(), initialize_interpolation_field_values(), and initialize_interpolation_grid_values().
|
private |
Definition at line 150 of file field_profiles.hpp.
Referenced by build_spline_for_field(), build_splines(), and initialize_interpolation_field_values().
|
private |
Definition at line 149 of file field_profiles.hpp.
Referenced by build_spline_for_field(), build_splines(), initialize_interpolation_field_values(), and initialize_interpolation_grid_values().
|
private |
Definition at line 148 of file field_profiles.hpp.
Referenced by initialize_interpolation_grid_values(), and set_interpolation_points_fraction().
|
private |
Definition at line 147 of file field_profiles.hpp.
Referenced by build_spline_for_field(), build_splines(), and set_interpolation_points_fraction().
|
private |
Definition at line 146 of file field_profiles.hpp.
Referenced by get_number_of_dimensions(), and set_number_of_dimensions().
|
private |
Definition at line 145 of file field_profiles.hpp.
Referenced by build_splines(), derivative_at(), evaluate_at(), get_field_profiles(), get_number_of_fields(), initialize_interpolation_field_values(), and set_field_profile().
|
private |
Definition at line 151 of file field_profiles.hpp.
Referenced by build_spline_for_field(), build_splines(), derivative_at(), and evaluate_at().