BubbleProfiler  0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
Public Member Functions | Private Member Functions | Private Attributes | List of all members
BubbleProfiler::Field_profiles Class Reference

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_profilesoperator= (const Field_profiles &)=default
 
 Field_profiles (Field_profiles &&)=default
 
Field_profilesoperator= (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_interpolatorsplines {}
 

Detailed Description

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.

Examples:
general_fubini.cpp, and logarithmic.cpp.

Definition at line 49 of file field_profiles.hpp.

Constructor & Destructor Documentation

BubbleProfiler::Field_profiles::Field_profiles ( )
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.

Parameters
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.

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.

Parameters
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.

BubbleProfiler::Field_profiles::~Field_profiles ( )
default
BubbleProfiler::Field_profiles::Field_profiles ( const Field_profiles )
default
BubbleProfiler::Field_profiles::Field_profiles ( Field_profiles &&  )
default

Member Function Documentation

void BubbleProfiler::Field_profiles::build_spline_for_field ( int  field)
private
void BubbleProfiler::Field_profiles::build_splines ( )
private
double BubbleProfiler::Field_profiles::derivative_at ( int  field,
int  order,
double  rho 
) const

Take the radial derivative of a field at a given radius.

Parameters
fieldtarget field index
orderorder of the derivative (1 or 2)
rhotarget radius
Returns
radial derivative of the field at target coordinate

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
double BubbleProfiler::Field_profiles::evaluate_at ( int  field,
double  rho 
) const

Get a specific field value at a given radius.

Parameters
fieldtarget field index
rhotarget radius
Returns
value of target field

Definition at line 174 of file field_profiles.cpp.

References profiles, and splines.

double BubbleProfiler::Field_profiles::get_domain_end ( ) const
inline
double BubbleProfiler::Field_profiles::get_domain_start ( ) const
inline
const Eigen::MatrixXd& BubbleProfiler::Field_profiles::get_field_profiles ( ) const
inline
int BubbleProfiler::Field_profiles::get_number_of_dimensions ( ) const
inline
int BubbleProfiler::Field_profiles::get_number_of_fields ( ) const
inline
int BubbleProfiler::Field_profiles::get_number_of_grid_points ( ) const
inline
const Eigen::VectorXd& BubbleProfiler::Field_profiles::get_spatial_grid ( ) const
inline
void BubbleProfiler::Field_profiles::initialize_interpolation_field_values ( int  field)
private
void BubbleProfiler::Field_profiles::initialize_interpolation_field_values ( )
private
void BubbleProfiler::Field_profiles::initialize_interpolation_grid_values ( )
private

Definition at line 87 of file field_profiles.cpp.

References gridpoints, interpolation_grid_values, and interpolation_points_fraction.

Referenced by build_splines().

Field_profiles& BubbleProfiler::Field_profiles::operator= ( const Field_profiles )
default
Field_profiles& BubbleProfiler::Field_profiles::operator= ( Field_profiles &&  )
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)
void BubbleProfiler::Field_profiles::set_number_of_dimensions ( int  d)
inline

Member Data Documentation

Eigen::VectorXd BubbleProfiler::Field_profiles::gridpoints {}
private
Eigen::MatrixXd BubbleProfiler::Field_profiles::interpolation_field_values {}
private
Eigen::VectorXd BubbleProfiler::Field_profiles::interpolation_grid_values {}
private
double BubbleProfiler::Field_profiles::interpolation_points_fraction {1.0}
private
bool BubbleProfiler::Field_profiles::must_update_coords {true}
private
double BubbleProfiler::Field_profiles::n_spatial_dims {3}
private

Definition at line 146 of file field_profiles.hpp.

Referenced by get_number_of_dimensions(), and set_number_of_dimensions().

Eigen::MatrixXd BubbleProfiler::Field_profiles::profiles {}
private
std::vector<GSL_interpolator> BubbleProfiler::Field_profiles::splines {}
private

The documentation for this class was generated from the following files: