BubbleProfiler
0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
|
Class implementing the system of ODEs obeyed by profile perturbations. More...
#include <perturbations_ode_system.hpp>
Public Member Functions | |
Perturbations_ODE_system (Potential &potential_, Field_profiles &profiles_, int n_spacetime_dimensions_) | |
Constructs the system of ODEs for the given potential. More... | |
void | operator() (const Eigen::VectorXd &eps, Eigen::VectorXd &depsdr, double r) const |
Calculate the derivatives of the perturbations. More... | |
Private Member Functions | |
bool | is_finite (const Eigen::VectorXd &) const |
Eigen::VectorXd | calculate_inhomogeneities (double) const |
Eigen::MatrixXd | calculate_mass_matrix (double) const |
Private Attributes | |
int | n_fields {1} |
int | n_spacetime_dimensions {3} |
Potential * | potential {nullptr} |
Field_profiles * | profiles {nullptr} |
Class implementing the system of ODEs obeyed by profile perturbations.
This class provides for the calculation of the first-order system of ODEs obeyed by the perturbations to the field profiles, and their first derivatives, that arise when using the perturbative method.
Definition at line 36 of file perturbations_ode_system.hpp.
BubbleProfiler::Perturbations_ODE_system::Perturbations_ODE_system | ( | Potential & | potential_, |
Field_profiles & | profiles_, | ||
int | n_spacetime_dimensions_ | ||
) |
Constructs the system of ODEs for the given potential.
[in] | potential_ | the potential to calculate profile perturbations for |
[in] | profiles_ | the current estimate for the bubble profile for which perturbations are to be computed |
[in] | n_spacetime_dimensions_ | the number of spacetime dimensions |
Definition at line 28 of file perturbations_ode_system.cpp.
|
private |
Definition at line 85 of file perturbations_ode_system.cpp.
References BubbleProfiler::Field_profiles::derivative_at(), BubbleProfiler::Field_profiles::evaluate_at(), n_fields, n_spacetime_dimensions, BubbleProfiler::Potential::partial(), potential, and profiles.
Referenced by operator()().
|
private |
Definition at line 103 of file perturbations_ode_system.cpp.
References BubbleProfiler::Field_profiles::evaluate_at(), n_fields, BubbleProfiler::Potential::partial(), potential, and profiles.
Referenced by operator()().
|
private |
Definition at line 70 of file perturbations_ode_system.cpp.
Referenced by operator()().
void BubbleProfiler::Perturbations_ODE_system::operator() | ( | const Eigen::VectorXd & | eps, |
Eigen::VectorXd & | depsdr, | ||
double | r | ||
) | const |
Calculate the derivatives of the perturbations.
The state vector eps
is assumed to contain the values of the perturbations in the first n_fields
entries, and the values of their first derivatives in the remaining n_fields
entries.
[in] | eps | the values of the perturbations and their derivatives |
[out] | depsdr | the values of the derivatives of the perturbations |
[in] | r | the value of the radial coordinate |
Definition at line 38 of file perturbations_ode_system.cpp.
References calculate_inhomogeneities(), calculate_mass_matrix(), is_finite(), n_fields, and n_spacetime_dimensions.
|
private |
Number of fields in the potential
Definition at line 64 of file perturbations_ode_system.hpp.
Referenced by calculate_inhomogeneities(), calculate_mass_matrix(), and operator()().
|
private |
Number of spacetime dimensions
Definition at line 66 of file perturbations_ode_system.hpp.
Referenced by calculate_inhomogeneities(), and operator()().
|
private |
Potential for which profile is calculated
Definition at line 68 of file perturbations_ode_system.hpp.
Referenced by calculate_inhomogeneities(), and calculate_mass_matrix().
|
private |
Current estimate for bubble profile to be perturbed
Definition at line 70 of file perturbations_ode_system.hpp.
Referenced by calculate_inhomogeneities(), and calculate_mass_matrix().