hub::orbit namespace

Contents

Documentation

Classes

struct Elliptic
Derived class of Kepler orbit. Elliptical orbit.
struct Hyperbolic
Derived class of Kepler orbit. Hyperbolic orbit.
template<typename Real>
struct KeplerOrbit
struct RandomIndicator
A place holder that indicates one of the three angles in orbital parameters will be randomly generated.

Enums

enum class OrbitType { Ellipse, Parabola, Hyperbola, None }
Enum of kepler orbit type. Possible value: Ellipse, Parabola, Hyperbola, None.
enum class Hyper { in, out }
Enum type that indicates the trajectory is hyperbolically incident in or hyperbolically eject out.

Typedefs

using Kepler = KeplerOrbit<double>
Alias of OrbitArgs<double>.

Functions

template<typename Scalar>
auto E_anomaly_to_T_anomaly(Scalar E_anomaly, Scalar e) -> Scalar
Calculate the corresponding true anomaly of the eccentric anomaly.
template<typename Scalar>
auto M_anomaly_to_E_anomaly(Scalar M_anomaly, Scalar e) -> Scalar
Calculate the corresponding eccentric anomaly of the mean anomaly.
template<typename Scalar>
auto M_anomaly_to_T_anomaly(Scalar M_anomaly, Scalar e) -> Scalar
Calculate the corresponding true anomaly of the mean anomaly.
template<typename Scalar>
auto T_anomaly_to_E_anomaly(Scalar T_anomaly, Scalar e) -> Scalar
Calculate the corresponding eccentric anomaly of the true anomaly.
template<typename Scalar>
auto E_anomaly_to_M_anomaly(Scalar E_anomaly, Scalar e) -> Scalar
Calculate the corresponding mean anomaly of the eccentric anomaly.
template<typename Scalar>
auto T_anomaly_to_M_anomaly(Scalar T_anomaly, Scalar e) -> Scalar
Calculate the corresponding mean anomaly of the true anomaly.
template<typename Vector, typename Scalar>
auto coord_to_orbit(Scalar m1, Scalar m2, const Vector& dr, const Vector& dv) -> auto
Tranfer relative position and velocity between two particles to Kepler orbit parameters.
template<typename Vector>
auto orbit_to_coord(KeplerOrbit<typename Vector::value_type> const& args) -> auto
Transfer Kepler orbit parameters to relative position and velocity between two component in orbit.
template<typename Particle, typename... Args>
auto group(Particle const& ptc1, Particle const& ptc2, Args const&... ptcs) -> auto
Create a std::ranges like(Container) from individual particles.
template<typename Cluster>
auto M_tot(Cluster const& ptc) -> auto
Calculate the total mass of a cluster of particles/single particle.
template<typename Particle, typename... Args>
auto M_tot(Particle const& ptc1, Particle const& ptc2, Args const&... args) -> auto
Calculate the total mass of particles.
template<typename Cluster>
auto COM_p(Cluster const& ptc) -> auto
Calculate the centre of mass position of a particle cluster/single particle.
template<typename Particle, typename... Args>
auto COM_p(Particle const& ptc1, Particle const& ptc2, Args const&... ptcs) -> auto
Calculate the centre of mass position of particles.
template<typename Cluster>
auto COM_v(Cluster const& ptc) -> auto
Calculate the centre of mass velocity of a particle cluster/single particle.
template<typename Particle, typename... Args>
auto COM_v(Particle const& ptc1, Particle const& ptc2, Args const&... ptcs) -> auto
Calculate the centre of mass velocity of particles.
template<typename Cluster1, typename Cluster2>
auto M_rdc(Cluster1 const& m1, Cluster2 const& m2) -> auto
The reduced mass of two clusters(cluster can also be a single particle).
template<typename Vector, typename Cluster>
void move_particles_pos(Vector const& centre_mass_pos, Cluster& ptc)
Move the centre of mass position of a cluster(can be a single particle) to a specific position.
template<typename Vector, typename Particle, typename... Args>
void move_particles_pos(Vector const& centre_mass_pos, Particle& ptc1, Particle& ptc2, Args&... ptcs)
Move the centre of mass position of particles to a specific position.
template<typename Vector, typename Cluster>
void move_particles_vel(Vector const& centre_mass_vel, Cluster& ptc)
Move the centre of mass velocity of a cluster(can be a single particle) to a specific velocity.
template<typename Vector, typename Particle, typename... Args>
void move_particles_vel(Vector const& centre_mass_vel, Particle& ptc1, Particle& ptc2, Args&... ptcs)
Move the centre of mass velocity of particles to a specific velocity.
template<typename Vector, typename Cluster>
void move_particles(Vector const& centre_mass_pos, Vector const& centre_mass_vel, Cluster& ptc)
Move the centre of mass position and velocity of a cluster(can be a single particle) to a specific position and velocity.
template<typename Vector, typename Particle, typename... Args>
void move_particles(Vector const& centre_mass_pos, Vector const& centre_mass_vel, Particle& ptc1, Particle& ptc2, Args&... ptcs)
Move the centre of mass position and velocity of particles to a specific position and velocity.
template<typename Particle, typename... Args>
void move_particles(KeplerOrbit<typename Particle::Scalar> const& orbit, Particle& ptc1, Args&... ptcs)
Move the centre of mass position and velocity of particles/a cluster of particles/single particle to the corresponding position and velocity of a Kepler orbit.
template<typename... Particle>
void move_to_COM_frame(Particle&... ptc)
Move the particles/a cluster of particles/single particle to the centre of mass frame and set the centre of mass to original point.
template<typename Cluster1, typename Cluster2>
auto calc_eccentricity(Cluster1 const& p1, Cluster2 const& p2) -> auto
Calculate the eccentricity of two clusters(cluster can also be a single particle) by regarding their centre of mass as point particle.
template<typename Cluster1, typename Cluster2>
auto calc_semi_major_axis(Cluster1 const& p1, Cluster2 const& p2) -> auto
Calculate the semi-major axis of two clusters(cluster can also be a single particle) by regarding their centre of mass as point particle.
template<typename Cluster1, typename Cluster2>
auto calc_a_e(Cluster1 const& p1, Cluster2 const& p2) -> auto
Calculate the semi-major axis and eccentricity of two clusters(cluster can also be a single particle) by regarding their centre of mass as point particle.
template<typename Cluster1, typename Cluster2>
auto period(Cluster1 const& p1, Cluster2 const& p2) -> auto
Calculate the semi-major axis of two clusters(cluster can also be a single particle) by regarding their centre of mass as point particle.
template<typename Cluster1, typename Cluster2>
auto time_to_periapsis(Cluster1 const& cluster1, Cluster2 const& cluster2) -> auto
Calculate the time to the periapsis of two clusters(cluster can also be a single particle) by regarding their centre of mass as point particle.
template<typename Cluster>
auto E_k(Cluster const& ptc) -> auto
Calculate the kinetic energy of a cluster/single particle.
template<typename Particle, typename... Args>
auto E_k(Particle const& ptc1, Particle const& ptc2, Args const&... args) -> auto
Calculate the kinetic energy of particles.
template<typename Cluster>
auto E_p(Cluster const& ptc) -> auto
Calculate the potential energy of a cluster/single particle.
template<typename Particle, typename... Args>
auto E_p(Particle const& ptc1, Particle const& ptc2, Args const&... args) -> auto
Calculate the kinetic energy of particles.
template<typename Cluster>
auto E_k_COM(Cluster const& ptc) -> auto
Calculate the kinetic energy of the centre of mass of a cluster/single particle.
template<typename Particle, typename... Args>
auto E_k_COM(Particle const& ptc1, Particle const& ptc2, Args const&... args) -> auto
Calculate the kinetic energy of the centre of mass of particles.
template<typename... Particle>
auto E_tot(Particle const&... ptc) -> auto
Calculate the total energy of particles/a cluster of particles/single particle.
template<typename... Particle>
auto E_inner(Particle const&... ptc) -> auto
Calculate the inner energy of particles/a cluster of particles/single particle.
template<typename Cluster>
auto cluster_size(Cluster const& ptc) -> auto
Calculate the size of a cluster/single particle.
template<typename Particle, typename... Args>
auto cluster_size(Particle const& ptc1, Particle const& ptc2, Args const&... args) -> auto
Calculate the size of a set of particles.
template<typename Cluster1, typename Cluster2, typename Scalar>
auto E_tid(Cluster1 const& cluster1, Cluster2 const& cluster2, Scalar R1, Scalar R2) -> auto
Calculate the tidal potential between two clusters/single particles by providing the size of clusters.
template<typename Cluster1, typename Cluster2>
auto E_tid(Cluster1 const& cluster1, Cluster2 const& cluster2) -> auto
Calculate the tidal potential between two clusters/single particles.
template<typename Cluster1, typename Cluster2, typename Scalar>
auto tidal_factor(Cluster1 const& cluster1, Cluster2 const& cluster2, Scalar R1, Scalar R2) -> auto
Calculate the tidal factors between two clusters/single particles by providing the size of clusters.
template<typename Cluster1, typename Cluster2>
auto tidal_factor(Cluster1 const& cluster1, Cluster2 const& cluster2) -> auto
Calculate the tidal factors between two clusters/single particles.
template<typename Cluster1, typename Cluster2, typename Scalar>
auto tidal_radius(Scalar tidal_factor, Cluster1 const& cluster1, Cluster2 const& cluster2, Scalar R2) -> auto
Calculate the tidal radius between two clusters.
template<typename Cluster1, typename Cluster2, typename Scalar>
auto tidal_radius(Scalar tidal_factor, Cluster1 const& cluster1, Cluster2 const& cluster2) -> auto
Calculate the tidal radius between two clusters.

Enum documentation

enum class hub::orbit::OrbitType

Enum of kepler orbit type. Possible value: Ellipse, Parabola, Hyperbola, None.

enum class hub::orbit::Hyper

Enum type that indicates the trajectory is hyperbolically incident in or hyperbolically eject out.

Typedef documentation

using hub::orbit::Kepler = KeplerOrbit<double>

Alias of OrbitArgs<double>.

Function documentation

template<typename Scalar>
Scalar hub::orbit::E_anomaly_to_T_anomaly(Scalar E_anomaly, Scalar e)

Calculate the corresponding true anomaly of the eccentric anomaly.

Template parameters
Scalar Floating point like type.
Parameters
E_anomaly Eccentric anomaly.
e Eccentricity.
Returns Scalar True anomaly.

template<typename Scalar>
Scalar hub::orbit::M_anomaly_to_E_anomaly(Scalar M_anomaly, Scalar e)

Calculate the corresponding eccentric anomaly of the mean anomaly.

Template parameters
Scalar Floating point like type.
Parameters
M_anomaly Mean anomaly.
e Eccentricity.
Returns Scalar Eccentric anomaly.

template<typename Scalar>
Scalar hub::orbit::M_anomaly_to_T_anomaly(Scalar M_anomaly, Scalar e)

Calculate the corresponding true anomaly of the mean anomaly.

Template parameters
Scalar Floating point like type.
Parameters
M_anomaly Mean anomaly.
e Eccentricity.
Returns Scalar True anomaly.

template<typename Scalar>
Scalar hub::orbit::T_anomaly_to_E_anomaly(Scalar T_anomaly, Scalar e)

Calculate the corresponding eccentric anomaly of the true anomaly.

Template parameters
Scalar Floating point like type.
Parameters
T_anomaly True anomaly.
e Eccentricity.
Returns Scalar Eccentric anomaly.

template<typename Scalar>
Scalar hub::orbit::E_anomaly_to_M_anomaly(Scalar E_anomaly, Scalar e)

Calculate the corresponding mean anomaly of the eccentric anomaly.

Template parameters
Scalar Floating point like type.
Parameters
E_anomaly Eccentric anomaly.
e Eccentricity.
Returns Scalar Mean anomaly.

template<typename Scalar>
Scalar hub::orbit::T_anomaly_to_M_anomaly(Scalar T_anomaly, Scalar e)

Calculate the corresponding mean anomaly of the true anomaly.

Template parameters
Scalar Floating point like type.
Parameters
T_anomaly True anomaly.
e Eccentricity.
Returns Scalar Mean anomaly.

template<typename Vector, typename Scalar>
auto hub::orbit::coord_to_orbit(Scalar m1, Scalar m2, const Vector& dr, const Vector& dv)

Tranfer relative position and velocity between two particles to Kepler orbit parameters.

Template parameters
Vector
Scalar
Parameters
m1
m2
dr
dv
Returns auto

template<typename Vector>
auto hub::orbit::orbit_to_coord(KeplerOrbit<typename Vector::value_type> const& args)

Transfer Kepler orbit parameters to relative position and velocity between two component in orbit.

Template parameters
Vector
Parameters
args
Returns auto

template<typename Particle, typename... Args>
auto hub::orbit::group(Particle const& ptc1, Particle const& ptc2, Args const&... ptcs)

Create a std::ranges like(Container) from individual particles.

Template parameters
Particle Type of particle.
Args Types of particle, should be same as Particle.
Parameters
ptc1 in The first particle.
ptc2 in The second particle.
ptcs in The rest if exits.
Returns auto The containers contains input particles.

template<typename Cluster>
auto hub::orbit::M_tot(Cluster const& ptc)

Calculate the total mass of a cluster of particles/single particle.

Template parameters
Cluster std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
ptc in particle container/single particle.
Returns auto The total mass of the particle cluster/single particle.

template<typename Particle, typename... Args>
auto hub::orbit::M_tot(Particle const& ptc1, Particle const& ptc2, Args const&... args)

Calculate the total mass of particles.

Template parameters
Particle Type of the particle with public member mass(Scalar), pos(Vector) and vel(Vector).
Args Type of the particles, should be same as Particle.
Parameters
ptc1 in The first particle.
ptc2 in The second particle.
args in The rest particles if exits.
Returns auto The total mass of particles

template<typename Cluster>
auto hub::orbit::COM_p(Cluster const& ptc)

Calculate the centre of mass position of a particle cluster/single particle.

Template parameters
Cluster std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
ptc in particle container/single particle.
Returns auto The centre of mass position of the particle cluster/single particle.

template<typename Particle, typename... Args>
auto hub::orbit::COM_p(Particle const& ptc1, Particle const& ptc2, Args const&... ptcs)

Calculate the centre of mass position of particles.

Template parameters
Particle Type of the particle with public member mass(Scalar), pos(Vector) and vel(Vector).
Args Type of the particles, should be same as Particle.
Parameters
ptc1 in The first particle.
ptc2 in The second particle.
ptcs in The rest particles if exist.
Returns auto The centre of mass position of particles.

template<typename Cluster>
auto hub::orbit::COM_v(Cluster const& ptc)

Calculate the centre of mass velocity of a particle cluster/single particle.

Template parameters
Cluster std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
ptc in particle container/single particle.
Returns auto The centre of mass velocity of the particle cluster/single particle.

template<typename Particle, typename... Args>
auto hub::orbit::COM_v(Particle const& ptc1, Particle const& ptc2, Args const&... ptcs)

Calculate the centre of mass velocity of particles.

Template parameters
Particle Type of the particle with public member mass(Scalar), pos(Vector) and vel(Vector).
Args Type of the particles, should be same as Particle.
Parameters
ptc1 in The first particle.
ptc2 in The second particle.
ptcs in The rest particles if exist.
Returns auto The centre of mass velocity of particles.

template<typename Cluster1, typename Cluster2>
auto hub::orbit::M_rdc(Cluster1 const& m1, Cluster2 const& m2)

The reduced mass of two clusters(cluster can also be a single particle).

Template parameters
Cluster1 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Cluster2 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
m1 in The first cluster/first single particle.
m2 in The second cluster/second single particle.
Returns auto The reduced mass of the two clusters(cluster can also be a single particle).

template<typename Vector, typename Cluster>
void hub::orbit::move_particles_pos(Vector const& centre_mass_pos, Cluster& ptc)

Move the centre of mass position of a cluster(can be a single particle) to a specific position.

Template parameters
Vector 3-D Vector type.
Cluster std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
centre_mass_pos in The target centre of mass position.
ptc in/out The cluster(can be a single particle) needs to be moved.

template<typename Vector, typename Particle, typename... Args>
void hub::orbit::move_particles_pos(Vector const& centre_mass_pos, Particle& ptc1, Particle& ptc2, Args&... ptcs)

Move the centre of mass position of particles to a specific position.

Template parameters
Vector 3-D Vector type.
Particle Type of the particle with public member mass(Scalar), pos(Vector) and vel(Vector).
Args Type of the particles, should be same as Particle.
Parameters
centre_mass_pos in The target centre of mass position.
ptc1 in/out The first particle needs to be moved.
ptc2 in/out The second particle needs to be moved.
ptcs in/out The rest particles need to be moved.

template<typename Vector, typename Cluster>
void hub::orbit::move_particles_vel(Vector const& centre_mass_vel, Cluster& ptc)

Move the centre of mass velocity of a cluster(can be a single particle) to a specific velocity.

Template parameters
Vector 3-D Vector type.
Cluster std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
centre_mass_vel in The target centre of mass position.
ptc in/out The cluster(can be a single particle) needs to be moved.

template<typename Vector, typename Particle, typename... Args>
void hub::orbit::move_particles_vel(Vector const& centre_mass_vel, Particle& ptc1, Particle& ptc2, Args&... ptcs)

Move the centre of mass velocity of particles to a specific velocity.

Template parameters
Vector 3-D Vector type.
Particle Type of the particle with public member mass(Scalar), pos(Vector) and vel(Vector).
Args Type of the particles, should be same as Particle.
Parameters
centre_mass_vel in The target centre of mass position.
ptc1 in/out The first particle needs to be moved.
ptc2 in/out The first particle needs to be moved.
ptcs in/out The rest particles need to be moved.

template<typename Vector, typename Cluster>
void hub::orbit::move_particles(Vector const& centre_mass_pos, Vector const& centre_mass_vel, Cluster& ptc)

Move the centre of mass position and velocity of a cluster(can be a single particle) to a specific position and velocity.

Template parameters
Vector 3-D Vector type.
Cluster std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
centre_mass_pos in The target centre of mass position.
centre_mass_vel in The target centre of mass velocity.
ptc in/out The cluster(can be a single particle) needs to be moved.

template<typename Vector, typename Particle, typename... Args>
void hub::orbit::move_particles(Vector const& centre_mass_pos, Vector const& centre_mass_vel, Particle& ptc1, Particle& ptc2, Args&... ptcs)

Move the centre of mass position and velocity of particles to a specific position and velocity.

Template parameters
Vector Vector 3-D Vector type.
Particle Type of the particle with public member mass(Scalar), pos(Vector) and vel(Vector).
Args Type of the particles, should be same as Particle.
Parameters
centre_mass_pos in The target centre of mass position.
centre_mass_vel in The target centre of mass velocity.
ptc1 in/out The first particle needs to be moved.
ptc2 in/out The second particle needs to be moved.
ptcs in/out The rest particles need to be moved.

template<typename Particle, typename... Args>
void hub::orbit::move_particles(KeplerOrbit<typename Particle::Scalar> const& orbit, Particle& ptc1, Args&... ptcs)

Move the centre of mass position and velocity of particles/a cluster of particles/single particle to the corresponding position and velocity of a Kepler orbit.

Template parameters
Particle Type of the first particle/std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)..
Args Type of the particles if exits, should be same as Particle.
Parameters
orbit in The Kepler orbit.
ptc1 in/out The first particle/The cluster/single particle needs to be moved.
ptcs in/out The rest particles need to be moved.

template<typename... Particle>
void hub::orbit::move_to_COM_frame(Particle&... ptc)

Move the particles/a cluster of particles/single particle to the centre of mass frame and set the centre of mass to original point.

Template parameters
Particle Type of the first particle/std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector).
Parameters
ptc in/out The particles/The cluster/single particle needs to be moved.

template<typename Cluster1, typename Cluster2>
auto hub::orbit::calc_eccentricity(Cluster1 const& p1, Cluster2 const& p2)

Calculate the eccentricity of two clusters(cluster can also be a single particle) by regarding their centre of mass as point particle.

Template parameters
Cluster1 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Cluster2 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
p1 in The first cluster/first single particle.
p2 in The second cluster/first single particle.
Returns auto The eccentricity.

template<typename Cluster1, typename Cluster2>
auto hub::orbit::calc_semi_major_axis(Cluster1 const& p1, Cluster2 const& p2)

Calculate the semi-major axis of two clusters(cluster can also be a single particle) by regarding their centre of mass as point particle.

Template parameters
Cluster1 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Cluster2 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
p1 in The first cluster/first single particle.
p2 in The second cluster/first single particle.
Returns auto The semi-major axis.

template<typename Cluster1, typename Cluster2>
auto hub::orbit::calc_a_e(Cluster1 const& p1, Cluster2 const& p2)

Calculate the semi-major axis and eccentricity of two clusters(cluster can also be a single particle) by regarding their centre of mass as point particle.

Template parameters
Cluster1 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Cluster2 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
p1 in The first cluster/first single particle.
p2 in The second cluster/first single particle.
Returns auto A tuple of (sem-major axis, eccentricity).

template<typename Cluster1, typename Cluster2>
auto hub::orbit::period(Cluster1 const& p1, Cluster2 const& p2)

Calculate the semi-major axis of two clusters(cluster can also be a single particle) by regarding their centre of mass as point particle.

Template parameters
Cluster1 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Cluster2 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
p1 in The first cluster/first single particle.
p2 in The second cluster/first single particle.
Returns auto The period.

template<typename Cluster1, typename Cluster2>
auto hub::orbit::time_to_periapsis(Cluster1 const& cluster1, Cluster2 const& cluster2)

Calculate the time to the periapsis of two clusters(cluster can also be a single particle) by regarding their centre of mass as point particle.

Template parameters
Cluster1 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Cluster2 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
cluster1 in first cluster/first single particle.
cluster2 in The second cluster/first single particle.
Returns auto The time to the periapsis.

template<typename Cluster>
auto hub::orbit::E_k(Cluster const& ptc)

Calculate the kinetic energy of a cluster/single particle.

Template parameters
Cluster std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
ptc in particle container/single particle.
Returns auto The kinetic energy.

template<typename Particle, typename... Args>
auto hub::orbit::E_k(Particle const& ptc1, Particle const& ptc2, Args const&... args)

Calculate the kinetic energy of particles.

Template parameters
Particle Type of the particle with public member mass(Scalar), pos(Vector) and vel(Vector).
Args Type of the particles, should be same as Particle.
Parameters
ptc1 in The first particle.
ptc2 in The second particle.
args in The rest particles if exits.
Returns auto The kinetic energy of particles.

template<typename Cluster>
auto hub::orbit::E_p(Cluster const& ptc)

Calculate the potential energy of a cluster/single particle.

Template parameters
Cluster std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
ptc in particle container/single particle.
Returns auto The potential energy.

template<typename Particle, typename... Args>
auto hub::orbit::E_p(Particle const& ptc1, Particle const& ptc2, Args const&... args)

Calculate the kinetic energy of particles.

Template parameters
Particle Type of the particle with public member mass(Scalar), pos(Vector) and vel(Vector).
Args Type of the particles, should be same as Particle.
Parameters
ptc1 in The first particle.
ptc2 in The second particle.
args in The rest particles if exits.
Returns auto The potential energy of particles.

template<typename Cluster>
auto hub::orbit::E_k_COM(Cluster const& ptc)

Calculate the kinetic energy of the centre of mass of a cluster/single particle.

Template parameters
Cluster std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
ptc in particle container/single particle.
Returns auto The kinetic energy of the centre of mass.

template<typename Particle, typename... Args>
auto hub::orbit::E_k_COM(Particle const& ptc1, Particle const& ptc2, Args const&... args)

Calculate the kinetic energy of the centre of mass of particles.

Template parameters
Particle Type of the particle with public member mass(Scalar), pos(Vector) and vel(Vector).
Args Type of the particles, should be same as Particle.
Parameters
ptc1 in The first particle.
ptc2 in The second particle.
args in The rest particles if exits.
Returns auto The kinetic energy of the centre of mass of particles.

template<typename... Particle>
auto hub::orbit::E_tot(Particle const&... ptc)

Calculate the total energy of particles/a cluster of particles/single particle.

Template parameters
Particle Type of the first particle/std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector).
Parameters
ptc in The particles/The cluster/single particle needs to be evaluated.
Returns auto The total energy.

template<typename... Particle>
auto hub::orbit::E_inner(Particle const&... ptc)

Calculate the inner energy of particles/a cluster of particles/single particle.

Template parameters
Particle Type of the first particle/std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector).
Parameters
ptc in The particles/The cluster/single particle needs to be evaluated.
Returns auto The inner energy.

template<typename Cluster>
auto hub::orbit::cluster_size(Cluster const& ptc)

Calculate the size of a cluster/single particle.

Template parameters
Cluster std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
ptc in particle container/single particle.
Returns auto The size of the cluster.

The size of the cluster are calculated by detecting the farthest distance between two particles in this cluster.

template<typename Particle, typename... Args>
auto hub::orbit::cluster_size(Particle const& ptc1, Particle const& ptc2, Args const&... args)

Calculate the size of a set of particles.

Template parameters
Particle Type of the particle with public member mass(Scalar), pos(Vector) and vel(Vector).
Args Type of the particles, should be same as Particle.
Parameters
ptc1 in The first particle.
ptc2 in The second particle.
args in The rest particles if exits.
Returns auto The size of the cluster.

The size of the cluster are calculated by detecting the farthest distance between two particles in this cluster.

template<typename Cluster1, typename Cluster2, typename Scalar>
auto hub::orbit::E_tid(Cluster1 const& cluster1, Cluster2 const& cluster2, Scalar R1, Scalar R2)

Calculate the tidal potential between two clusters/single particles by providing the size of clusters.

Template parameters
Cluster1 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Cluster2 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Scalar Floating point like type.
Parameters
cluster1 in first cluster/first single particle.
cluster2 in The second cluster/first single particle.
R1 in The size of the first cluster.
R2 in The size of the second cluster.
Returns auto The tidal potential energy.

template<typename Cluster1, typename Cluster2>
auto hub::orbit::E_tid(Cluster1 const& cluster1, Cluster2 const& cluster2)

Calculate the tidal potential between two clusters/single particles.

Template parameters
Cluster1 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Cluster2 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
cluster1 in The first cluster/first single particle.
cluster2 in The second cluster/first single particle.
Returns auto The tidal potential energy.

The size of the clusters will be estimated by cluster_size().

template<typename Cluster1, typename Cluster2, typename Scalar>
auto hub::orbit::tidal_factor(Cluster1 const& cluster1, Cluster2 const& cluster2, Scalar R1, Scalar R2)

Calculate the tidal factors between two clusters/single particles by providing the size of clusters.

Template parameters
Cluster1 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Cluster2 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Scalar Floating point like type.
Parameters
cluster1 in The first cluster/first single particle.
cluster2 in The second cluster/first single particle.
R1 in The size of the first cluster.
R2 in The size of the second cluster.
Returns auto A tuple of (factor1, factor2), where factor1 is the tidal factor of first cluster and factor2 is for the second.

template<typename Cluster1, typename Cluster2>
auto hub::orbit::tidal_factor(Cluster1 const& cluster1, Cluster2 const& cluster2)

Calculate the tidal factors between two clusters/single particles.

Template parameters
Cluster1 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Cluster2 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Parameters
cluster1 in The first cluster/first single particle.
cluster2 in The second cluster/first single particle.
Returns auto A tuple of (factor1, factor2), where factor1 is the tidal factor of first cluster and factor2 is for the second.

The size of the clusters will be estimated by cluster_size().

template<typename Cluster1, typename Cluster2, typename Scalar>
auto hub::orbit::tidal_radius(Scalar tidal_factor, Cluster1 const& cluster1, Cluster2 const& cluster2, Scalar R2)

Calculate the tidal radius between two clusters.

Template parameters
Cluster1 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Cluster2 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Scalar Floating point like type.
Parameters
tidal_factor in The tidal factor indicates where the radius is estimated. =1 means tidal disrupted.
cluster1 in The first cluster/first single particle.
cluster2 in The second cluster/first single particle.
R2 The size of the second cluster.
Returns auto The tidal radius of the second cluster tidaled by the first cluster.

template<typename Cluster1, typename Cluster2, typename Scalar>
auto hub::orbit::tidal_radius(Scalar tidal_factor, Cluster1 const& cluster1, Cluster2 const& cluster2)

Calculate the tidal radius between two clusters.

Template parameters
Cluster1 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Cluster2 std::ranges(Container) with element type has public member mass(Scalar), pos(Vector) and vel(Vector)./Type of single particle.
Scalar Floating point like type.
Parameters
tidal_factor in The tidal factor indicates where the radius is estimated. =1 means tidal disrupted.
cluster1 in The first cluster/first single particle.
cluster2 in The second cluster/first single particle.
Returns auto The tidal radius of the second cluster tidaled by the first cluster.