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

Abstract base class for a generic potential. More...

#include <potential.hpp>

Inheritance diagram for BubbleProfiler::Potential:
Inheritance graph
[legend]

Public Member Functions

virtual ~Potential ()=default
 
virtual Potentialclone () const =0
 Subclasses must implement a clone method. More...
 
virtual double operator() (const Eigen::VectorXd &coords) const =0
 Evaluate potential at point. More...
 
virtual double partial (const Eigen::VectorXd &coords, int i) const =0
 Partial derivative WRT coordinate i at a point. More...
 
virtual double partial (const Eigen::VectorXd &coords, int i, int j) const =0
 Partial derivative WRT coordinates i, j at a a point. More...
 
virtual std::size_t get_number_of_fields () const =0
 
virtual void translate_origin (const Eigen::VectorXd &translation)=0
 Shift the location of the origin by a specified vector. More...
 
virtual void apply_basis_change (const Eigen::MatrixXd &cob_matrix)=0
 Apply a change of basis matrix. More...
 
virtual void add_constant_term (double offset)=0
 Add a constant offset to the potential. More...
 

Detailed Description

Abstract base class for a generic potential.

Implementations must supply methods for evaluation the potential at a point in field space, and the computation of first and second partial derivatives with respect to field space coordinates.

Additionally, implementations must provide methods for translating the origin of the potential, applying a basis change, and adding a constant term.

Definition at line 36 of file potential.hpp.

Constructor & Destructor Documentation

virtual BubbleProfiler::Potential::~Potential ( )
virtualdefault

Member Function Documentation

virtual void BubbleProfiler::Potential::add_constant_term ( double  offset)
pure virtual
virtual void BubbleProfiler::Potential::apply_basis_change ( const Eigen::MatrixXd &  cob_matrix)
pure virtual
virtual Potential* BubbleProfiler::Potential::clone ( ) const
pure virtual
virtual std::size_t BubbleProfiler::Potential::get_number_of_fields ( ) const
pure virtual
virtual double BubbleProfiler::Potential::operator() ( const Eigen::VectorXd &  coords) const
pure virtual
virtual double BubbleProfiler::Potential::partial ( const Eigen::VectorXd &  coords,
int  i 
) const
pure virtual
virtual double BubbleProfiler::Potential::partial ( const Eigen::VectorXd &  coords,
int  i,
int  j 
) const
pure virtual

Partial derivative WRT coordinates i, j at a a point.

Parameters
coordsCoordinates at which to evaluate
iIndex of first coordinate to be differentiated
jIndex of second coordinate to be differentiated
Returns
Value of specified partial at point

Implemented in BubbleProfiler::Gaussian_potential, BubbleProfiler::Algebraic_potential, BubbleProfiler::Thin_wall_potential, BubbleProfiler::Generalized_fubini_potential, BubbleProfiler::Logarithmic_potential, and BubbleProfiler::Restricted_quartic_potential.

virtual void BubbleProfiler::Potential::translate_origin ( const Eigen::VectorXd &  translation)
pure virtual

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