clustertools.cluster.cluster module¶
The StarCluster class and key internal functions
- class clustertools.cluster.cluster.StarCluster(tphys=0.0, units=None, origin=None, ctype='snapshot', projected=False, **kwargs)[source]¶
Bases:
object
A class that represents a star cluster population that ooperations and functions can be performed on
- Parameters:
- tphysfloat
Time (units not necessary) associated with the population (default: 0)
- unitsstr
Units of stellar positions and velocties. Options include ‘pckms’,’pcmyr’, ‘kpckms’,’kpcgyr’,’radec’,’nbody’,and ‘galpy’. For ‘pckms’ and ‘kpckms’, stellar velocities are assumed to be km/s. (default: None)
- originstr
Origin of coordinate system within which stellar positions and velocities are defined. Options include ‘centre’, ‘cluster’, ‘galaxy’ and ‘sky’. Note that ‘centre’ corresponds to the systems centre of density or mass (as calculated by clustertools) while ‘cluster’ corresponds to the orbital position of the cluster’s initial centre of mass given ‘tphys’. In some cases these two coordinate systems may be the same. (default: None)
- ctypestr
Code type used to generate the cluster population. Current options include ‘snapshot’, ‘nbody6’,’nbody6se’,’nemo’,’gyrfalcon’,’snaptrim’, amuse’,’clustertools’,’snapauto’. The parameter informs clustertools how to load the stellar popuation and advance to the next snapshot. (default: ‘snapshot’)
- projectedbool
return projected values instead of 3D value. (default: False)
- Other Parameters:
- sfilestr
file containing single star data
- bfilestr
file contain binary star data
- ssefilestr
file containing single star evolution data
- bsefilestr
file contain binary star evolution data
- tfilestr
name of file contain tail star data
- ofilestr
orbit file
- ofilenamestr
orbit filename if ofile is not given
- orbitstr
galpy orbit instance
- ounitsstr
{‘pckms’,’pcmyr’,’kpckms’,’kpcgyr’,’radec’,’nbody’,’galpy’} units of orbital information (else assumed equal to StarCluster.units)
- rofloat
distance to the Galactic centre (Default: solar_ro)
- vofloat
circular velocity at ro (Default: solar_vo)
- zofloat
Sun’s distance above the Galactic plane (default: solar_zo)
- solarmotionfloat
array representing U,V,W of Sun (default: solar_motion)
- nsnapint
if a specific snapshot is to be read in instead of starting from zero
- nzfillint
value for zfill when reading and writing snapshots (default: 5)
- snapbasestr
string of characters in filename before nsnap (default: ‘’)
- snapendstr
string of character in filename after nsnap (default: ‘.dat’)
- snapdirstr
string name of directory of snapshots if different than wdir (Default: ‘’)
- delimiterstr
choice of delimiter when reading ascii/csv files (default: ‘,’)
- wdirstr
working directory of snapshots if not current directory
- intializebool
initialize a galpy orbit after reading in orbital information (default: False)
- advancebool
set to True if this a snapshot that has been advanced to from an initial one? (default: False)
- centre_method: str
{None,’orthographic’,’VandeVen’} method to convert to clustercentric coordinates when units are in radec (Default: None)
- givestr
set what parameters are read in from nemo/gyrfalcon (default: ‘mxv’) Currently only accepts ‘mxvpqael’ as an alternative.
- History
- ——-
- 2018 - Written - Webb (UofT)
- add_action(JR, Jphi, Jz, OR=None, Ophi=None, Oz=None, TR=None, Tphi=None, Tz=None)[source]¶
Add action values to the cluster instance
- Parameters:
- JR,Jphi,Jzfloat
orbit actions
- OR,Ophi,Ozfloat
orbit frequencies
- TR,Tphi,Tzfloat
orbit periods
- Returns:
- None
- add_actions(JR, Jphi, Jz, OR=None, Ophi=None, Oz=None, TR=None, Tphi=None, Tz=None)[source]¶
Add action values to the cluster instance
- Parameters:
- JR,Jphi,Jzfloat
orbit actions
- OR,Ophi,Ozfloat
orbit frequencies
- TR,Tphi,Tzfloat
orbit periods
- Returns:
- None
- add_binary_stars(xb1, yb1, zb1, vxb1, vyb1, vzb1, xb2, yb2, zb2, vxb2, vyb2, vzb2, mb1=None, mb2=None, id1=None, id2=None, m01=None, m02=None, npop1=None, npop2=None, set_com=True, return_com=False)[source]¶
Individually add binary stars to StarCluster.
- Parameters:
- xb1,yb1,zb1: float
stellar positions of primary. Input is assumed to be in cartesian coordinates unless self.units==’radec’ and self.origin==’sky’, then positions are assumed to be ra,dec,dist (degrees, degrees, kpc)
- vxb1,vyb1,vzb1: float
stellar velocities of primary. Input is assumed to be in cartesian coordinates unless self.units==’radec’ and self.origin==’sky’, then positions are assumed to be pmra,pmdec,vlos (mas/yr, mas/yr, km/s)
- xb2,yb2,zb2: float
stellar positions of secondary. Input is assumed to be in cartesian coordinates unless self.units==’radec’ and self.origin==’sky’, then positions are assumed to be ra,dec,dist (degrees, degrees, kpc)
- vxb2,vyb2,vzb2: float
stellar velocities of secondary. Input is assumed to be in cartesian coordinates unless self.units==’radec’ and self.origin==’sky’, then positions are assumed to be pmra,pmdec,vlos (mas/yr, mas/yr, km/s)
- mb1. mb2float
mass of primary and secondary (default: None)
- id1,id2: int
primary and secondary star ids (default: None)
- m01,m02: float
initial primary and secondary stellar masses (default: None)
- npop1,npop2int
population number of primary and secondary star, for use with multiple populations (default: None)
- set_combool
set coordinates of binary’s centre of mass in main arrays (default: False)
- return_combool
reuturn coordinates of binary’s centre of mass (default: False)
Notes
History:
2022 - Written - Webb (UofT)
- add_bse(id1, id2, kw1, kw2, kcm, ecc, pb, semi, m1, m2, logl1, logl2, logr1, logr2, ep1=None, ep2=None, ospin1=None, ospin2=None)[source]¶
Add binary star evolution information to stars
parameters are common output variables in NBODY6
values are never adjusted during unit or coordinate changes
- Parameters:
- id1/id2int
id of star1 and star2
- kw1/kw2int
stellar evolution tags (for using with NBODY6)
- kcmint
stellar evolution tag for binary star
- eccfloat
eccentricity of binary orbit
- pb: float
period of binary orbit
- semifloat
semi major axis of binary orbit
- m1/m2float
masses of binary stars
- logl1/logl2float
luminosities of binary stars
- logr1/logr2float
radii of binary stars
- ep1/ep2float
epochs of binary stars
- ospin1/ospin2float
opsin of binary stars
- Returns:
- None
- add_energies(kin, pot, etot=None)[source]¶
Add energy information for stars
total energy and Q for the cluster are also calculated
values are never adjusted during unit or coordinate changes
- Parameters:
- kinfloat
kinetic energy
- potfloat
potentail energy
- etotfloat
total energy - calculated as kin+pot if not given
- Returns:
- None
- add_nbody6(nc=0, rc=0.0, rbar=1.0, rtide=0.0, xc=0.0, yc=0.0, zc=0.0, zmbar=1.0, vbar=1.0, tbar=1.0, rscale=1.0, ns=0.0, nb=0.0, n_p=0.0)[source]¶
Add additional information to StarCluster
parameters are common output variables in NBODY6
values are never adjusted during unit or coordinate changes
- Parameters:
- ncint
number of stars in core (default:0)
- rcint
core radius (default:0)
- rbarfloat
scaling factor between NBODY units and pc (default:1.)
- rtide
rtide set in the simulation (default:0)
- xc,yc,zcfloat
position of cluster centre (default:0)
- zmbarfloat
scaling factor between NBODY units and Msun (default:1.)
- vbarfloat
scaling factor between NBODY units and km/s (default:1.)
- tbarfloat
scaling factor between NBODY units and Myr (default:1.)
- rscalefloat
the scale radius of data (default:1.)
- nsint
number of single stars (default:0)
- nbint
number of binary stars (default:0)
- npint
number of particles (default:0)
- Returns:
- None
- add_orbit(xgc, ygc, zgc, vxgc, vygc, vzgc, ounits=None, initialize=False, from_centre=False, tphys=None)[source]¶
Add orbit properties to StarCluster
- Parameters:
- xgc,ygc,zgc: float
cluster’s galactocentric position
- vxgc,vygc,vzgc: float
cluster’s galactocentric velocity
- ounits: str
units of position and velocity. Options include ‘pckms’,’pcmyr’ ‘kpckms’,’kpcgyr’,’radec’,’nbody’,and ‘galpy’. Values will be converted to match self.units
- initialize: bool
Initialize a galpy orbit for self.orbit (default: False)
- from_centrebool
genrate orbit from cluster’s exact centre instead of its assigned galactocentric coordinates (default: False)
- tphysfloat
physical time as per the orbit file
- Returns:
- None
- add_rotation(qrot)[source]¶
Add a degree of rotation to an already generated StarCluster
- Parameters:
- clusterclass
StarCluster
- qrotfloat
fraction of stars with v_phi < 0 that are switched to vphi > 0
- Returns:
- x,y,z,vx,vy,vzfloat
stellar positions and velocities (now with rotation)
- add_sse(kw, logl, logr, ep=None, ospin=None, arg=None)[source]¶
Add stellar evolution information to stars
parameters are common output variables in NBODY6
values are never adjusted during unit or coordinate changes
- Parameters:
- kwint
stellar evolution type (based on NBODY6)
- loglfloat
log of luminosity
- logrfloat
log of stellar radius
- epfloat
epoch
- ospinfloat
ospin
- argint
array address arguments if SSE parameters do not necessarily match position and velocity arrays
- Returns:
- None
- add_stars(x, y, z, vx, vy, vz, m=None, id=None, m0=None, npop=None, nb=0, sortstars=False, analyze=True)[source]¶
Add stars to StarCluster.
- Parameters:
- x,y,z: float
stellar positions. Input is assumed to be in cartesian coordinates unless self.units==’radec’ and self.origin==’sky’, then positions are assumed to be ra,dec,dist (degrees, degrees, kpc)
- vx,vy,vz: float
atellar velocities. Input is assumed to be in cartesian coordinates unless self.units==’radec’ and self.origin==’sky’, then positions are assumed to be pmra,pmdec,vlos (mas/yr, mas/yr, km/s)
- m: float int
stellar mass
- id: int
star id
- m0: float int
initial stellar mass
- npopint
population number, for use with multiple populations
- nb: int
number of binary stars, which are assumed to be the first 2*nb stars in the array (default:0)
- sortstars: bool
order stars by radius (default: False)
- analyzebool
perform analysis on cluster after stars are added
Notes
History:
2018 - Written - Webb (UofT)
- analyse(sortstars=True, projected=True)[source]¶
Call analyze with alternative spelling
- Parameters:
- sortstarsbool
sort star by radius after coordinate change (default: True)
- projectedbool
sort projected radii as well, but do not change self.projected (default: True)
- Returns:
- None
- analyze(sortstars=True, projected=True)[source]¶
Calculate properties related to mass, radius, and velocity
- Parameters:
- sortstarsbool
sort star by radius after coordinate change (default: True)
- projectedbool
sort projected radii as well, but do not change self.projected (default: True)
- Returns:
- None
- ckin(mmin=None, mmax=None, nmass=10, rmin=None, rmax=None, vmin=None, vmax=None, emin=None, emax=None, kwmin=0, kwmax=1, indx=None, projected=False)[source]¶
NAME: Find the kinematic concentration parameter ck
see Bianchini et al. 2018, MNRAS, 475, 96
- Parameters:
- clusterclass
StarCluster instance
- mmin/mmaxfloat
specific mass range
- nmass
number of mass bins used to calculate alpha
- rmin/rmax
specific radial range
- vmin/vmaxfloat
specific velocity range
- emin/emaxfloat
specific energy range
- kwmin/kwmaxint
specific stellar evolution type range
- npopint
population number
- indxbool
specific subset of stars
- projectedbool
use projected values (default: False)
- Returns:
- ckfloat
kinematic concentration
- Other Parameters:
- kwargsstr
key words for plotting
- closest_star(projected=None, argmument=False)[source]¶
Find distance to closest star for each star - uses numba
- Parameters:
- clusterclass
positions of stars within the StarCluster
- projectedbool
use projected values (default: False)
- argumentbool
return argument of closest star as well (default: False)
- Returns:
- minimum_distancefloat
distance to closest star for each star
- if argument:
- argint
argument of closest star for each star
- core_relaxation_time(coulomb=0.4, projected=None)[source]¶
Calculate the core relaxation time (Stone & Ostriker 2015) of the cluster
Stone, N.C. & Ostriker, J.P. 2015, ApJ, 806, 28
Parameters
- clusterclass
StarCluster
- coulombfloat
Coulomb parameter (default: 0.4)
- projectedbool
use projected values (default: False)
- methodstr
choose between Stone & Ostriker 2015 and other methods (in development)
Returns
trc
- energies(specific=True, ids=None, projected=None, softening=0.0, full=True, parallel=False, **kwargs)[source]¶
Calculate kinetic and potential energy of every star
- Parameters:
- clusterclass
StarCluster instance
- specificbool
find specific energies (default: True)
- ids: boolean array or integer array
if given, find the energues of a subset of stars defined either by an array of star ids, or a boolean array that can be used to slice the cluster (default: None)
- projectedbool
use projected values (default: False)
- softeningfloat
Plummer softening length in cluster.units (default: 0.0)
- fullbool
calculate distance of full array of stars at once with numbra (default: True)
- parallelbool
calculate distances in parallel (default: False)
- Returns:
- kin,potfloat
kinetic and potential energy of every star if the i_d argument is not used. If i_d argument is used, return an arrays with potential and kinetic energy in the same shape of i_d
- History
2019 - Written - Webb (UofT) 2022 - Updated with support for multiple ids or an idexing array - Gillis (UofT)
- eta_function(mmin=None, mmax=None, nmass=10, rmin=None, rmax=None, vmin=None, vmax=None, emin=None, emax=None, kwmin=0, kwmax=1, indx=None, projected=False, plot=False, **kwargs)[source]¶
NAME: Find power slope of velocity dispersion versus mass
mass bins are set up so that there are an equal number of stars in each bin
- Parameters:
- clusterclass
StarCluster instance
- mmin/mmaxfloat
specific mass range
- nmass
number of mass bins used to calculate alpha
- rmin/rmax
specific radial range
- vmin/vmaxfloat
specific velocity range
- emin/emaxfloat
specific energy range
- kwmin/kwmaxint
specific stellar evolution type range
- npopint
population number
- indxbool
specific subset of stars
- projectedbool
use projected values (default: False)
- plotbool
plot the mass function
- Returns:
- m_meanfloat
mean mass in each bin
- sigvmfloat
velocity dispersion of stars in each bin
- etafloat
power-law slope of (sigvm ~ m^eta)
- eetafloat
error in eta
- yetafloat
y-intercept of fit to log(sigvm) vs log(m)
- eetafloat
error in yeta
- Other Parameters:
- kwargsstr
key words for plotting
- find_centre(xstart=0.0, ystart=0.0, zstart=0.0, vxstart=0.0, vystart=0.0, vzstart=0.0, indx=None, nsigma=1.0, nsphere=100, density=True, rmin=0.1, rmax=None, nmax=100, method='harfst', nneighbour=6, reset_centre=False)[source]¶
Find the cluster’s centre
The default assumes the cluster’s centre is the centre of density, calculated via the find_centre_of_density function.
For density=False, the routine first works to identify a sphere of nsphere stars around the centre in which to perform a centre of mass calculation (similar to NBODY6). Stars beyond nsigma standard deviations are removed from the calculation until only nsphere stars remain. This step prevents long tidal tails from affecting the calculation
- Parameters:
- clusterclass
StarCluster
- xstart,ystart,zstartfloat
starting position for centre (default: 0,0,0)
- vxstart,vystart,vzstartfloat
starting velocity for centre (default: 0,0,0)
- indx: bool
subset of stars to perform centre of density calculation on (default: None)
- nsphereint
number of stars in centre sphere (default:100)
- rminfloat
minimum radius of sphere around which to estimate density centre (default: 0.1 cluster.units)
- rmaxfloat
maxmimum radius of sphere around which to estimate density centre (default: None cluster.units, uses maximum r)
- nmaxfloat
maximum number of iterations (default:100)
- methodstr
method for finding the centre of density (‘harfst’ (default), ‘casertano’)
- if method==’casertano’
- nneighbourint
number of neighbours for calculation local densities
- reset_centrebool
forcibly reset cluster’s centre of mass (default: False)
- Returns:
- xc,yc,zc,vxc,vyc,vzc - coordinates of centre of mass
- find_centre_of_density(xstart=0.0, ystart=0.0, zstart=0.0, vxstart=0.0, vystart=0.0, vzstart=0.0, indx=None, nsphere=100, rmin=0.1, rmax=None, nmax=100, method='harfst', nneighbour=6, reset_centre=False)[source]¶
Find cluster’s centre of density
- Find cluster’s centre of density:
if method==’harfst’ use (Harfst, S., Gualandris, A., Merritt, D., et al. 2007, NewA, 12, 357) courtesy of Yohai Meiron if method==’casertano’ use (Casertano, S., Hut, P. 1985, ApJ, 298, 80)
- Parameters:
- clusterclass
StarCluster
- xstart,ystart,zstartfloat
starting position for centre (default: 0,0,0)
- vxstart,vystart,vzstartfloat
starting velocity for centre (default: 0,0,0)
- indx: bool
subset of stars to perform centre of density calculation on (default: None)
- nsphereint
number of stars in centre sphere (default:100)
- rminfloat
minimum radius of sphere around which to estimate density centre (default: 0.1 cluster.units)
- rmaxfloat
maxmimum radius of sphere around which to estimate density centre (default: None cluster.units, uses maximum r)
- nmaxfloat
maximum number of iterations (default:100)
- methodstr
method for finding the centre of density (‘harfst’ (default), ‘casertano’)
- if method==’casertano’
- nneighbourint
number of neighbours for calculation local densities
- reset_centrebool
forcibly reset cluster’s centre of mass (default: False)
- Returns:
- xc,yc,zc,vxc,vyc,vzcfloat
coordinates of centre of mass
- find_centre_of_mass(reset_centre=False)[source]¶
Find the centre of mass of the cluster
- Parameters:
- clusterclass
StarCluster
- Returns:
- xc,yc,zc,vxc,vyc,vzcfloat
coordinates of centre of mass
- half_mass_relaxation_time(coulomb=0.4, projected=None)[source]¶
Calculate the half-mass relaxation time (Spitzer 1987) of the cluster - Spitzer, L. 1987, Dynamical evolution of globular clusters - Need to adjust amplitude for different input units
- Parameters:
- clusterclass
StarCluster
- coulombfloat
Coulomb parameter (default: 0.4)
- projectedbool
use projected values (default: False)
- Returns:
- trhfloat
half-mass relaxation time within radius rad
- initialize_orbit(from_centre=False, ro=None, vo=None, zo=None, solarmotion=None)[source]¶
Initialize a galpy orbit instance for the cluster
- Parameters:
- clusterclass
StarCluster
- from_centrebool
intialize orbits from cluster’s exact centre instead of cluster’s position in galaxy (default :False)
- rofloat
galpy distance scale (default: None)
- vofloat
galpy velocity scale (default: None)
- zofloat
Sun’s distance above the Galactic plane (default: None)
- solarmotionfloat
array representing U,V,W of Sun (default: None)
- Returns:
- orbitclass
GALPY orbit
- initialize_orbits(ro=None, vo=None, zo=None, solarmotion=None)[source]¶
Initialize a galpy orbit for every star in the cluster
- Parameters:
- clusterclass
StarCluster
- rofloat
galpy distance scale (default: None)
- vofloat
galpy velocity scale (default: None)
- zofloat
Sun’s distance above the Galactic plane (default: None)
- solarmotionfloat
array representing U,V,W of Sun (default: None)
- Returns:
- orbitclass
GALPY orbit
- interpolate_orbit(pot=None, tfinal=None, nt=1000, from_centre=False, ro=None, vo=None, zo=None, solarmotion=None)[source]¶
Interpolate past or future position of cluster and escaped stars
- Parameters:
- clusterclass
StarCluster
- cluster_potclass
Galpy potential for host cluster that orbit is to be integrated in if None, assume a Plumme Potential
- potclass
galpy Potential that orbit is to be integrate in (default: None)
- tfinalfloat
final time (in cluster.units) to integrate orbit to (default: 12 Gyr)
- ntint
number of timesteps
- from_centrebool
intialize orbits from cluster’s exact centre instead of cluster’s position in galaxy (default :False)
- ro :float
galpy distance scale (Default: None)
- vofloat
galpy velocity scale (Default: None)
- zofloat
Sun’s distance above the Galactic plane (default: None)
- solarmotionfloat
array representing U,V,W of Sun (default: None)
- Returns:
- x,y,zfloat
interpolated positions of each star
- vx,vy,vzfloat
interpolated velocities of each star
- interpolate_orbits(pot=None, tfinal=None, nt=1000, ro=None, vo=None, zo=None, solarmotion=None)[source]¶
Interpolate past or future position of stars within the cluster
- Parameters:
- clusterclass
StarCluster
- potclass
Galpy potential for host cluster that orbit is to be integrated in if None, assume a Plumme Potential
- tfinalfloat
final time (in cluster.units) to integrate orbit to (default: 12 Gyr)
- ntint
number of timesteps
- ro :float
galpy distance scale (Default: None)
- vofloat
galpy velocity scale (Default: None)
- zofloat
Sun’s distance above the Galactic plane (default: None)
- solarmotionfloat
array representing U,V,W of Sun (default: None)
- Returns:
- x,y,zfloat
interpolated positions of each star
- vx,vy,vzfloat
interpolated velocities of each star
- key_params(do_order=True, projected=True)[source]¶
Call analyze with key_params for backwards compatibility
- Parameters:
- do_orderbool
sort star by radius after coordinate change (default: True)
- projectedbool
sort projected radii as well, but do not change self.projected (default: True)
- Returns:
- None
- mass_function(mmin=None, mmax=None, nmass=10, rmin=None, rmax=None, vmin=None, vmax=None, emin=None, emax=None, kwmin=0, kwmax=1, indx=None, mcorr=None, projected=False, plot=False, **kwargs)[source]¶
Find mass function over a given mass range
mass bins are set up so that there are an equal number of stars in each bin
- Parameters:
- clusterclass
StarCluster instance
- mmin/mmaxfloat
specific mass range
- nmass
number of mass bins used to calculate alpha
- rmin/rmax
specific radial range
- vmin/vmaxfloat
specific velocity range
- emin/emaxfloat
specific energy range
- kwmin/kwmaxint
specific stellar evolution type range
- npopint
population number
- indxbool
specific subset of stars
- projectedbool
use projected values (default: False)
- mcorrfloat
completeness correction for masses
- plotbool
plot the mass function
- Returns:
- m_meanfloat
mean mass in each bin
- m_histfloat
number of stars in each bin
- dmfloat
dN/dm of each bin
- alphafloat
power-law slope of the mass function (dN/dm ~ m^alpha)
- ealphafloat
error in alpha
- yalphafloat
y-intercept of fit to log(dN/dm) vs log(m)
- eyalphafloat
error in yalpha
- Other Parameters:
- kwargsstr
key words for plotting
- meq_function(mmin=None, mmax=None, nmass=10, rmin=None, rmax=None, vmin=None, vmax=None, emin=None, emax=None, kwmin=0, kwmax=1, indx=None, projected=False, plot=False, **kwargs)[source]¶
NAME: Find meq from velocity dispersion versus mass
mass bins are set up so that there are an equal number of stars in each bin
As per Bianchini, P. et al. 2016, MNRAS, 458, 3644, velocity dispersion versus mass is fit with the following: sigma(m)= sigma e^(-1/2 m/meq) if m<= meq
= sigma0 e^(-1/2) (m/meq)^-1/2 if m > meq
- Parameters:
- clusterclass
StarCluster instance
- mmin/mmaxfloat
specific mass range
- nmass
number of mass bins used to calculate alpha
- rmin/rmax
specific radial range
- vmin/vmaxfloat
specific velocity range
- emin/emaxfloat
specific energy range
- kwmin/kwmaxint
specific stellar evolution type range
- npopint
population number
- indxbool
specific subset of stars
- projectedbool
use projected values (default: False)
- plotbool
plot the mass function
- Returns:
- m_meanfloat
mean mass in each bin
- sigvmfloat
velocity dispersion of stars in each bin
- meqfloat
Bianchini fit to sigvm vs m
- emeqfloat
error in Bianchini fit to sigvm vs m
- sigma0float
Bianchini fit to sigvm vs m
- esigma0float
error in Bianchini fit to sigvm vs m
- Other Parameters:
- kwargsstr
key words for plotting
- orbit_interpolate(pot=None, tfinal=None, nt=1000, from_centre=False, ro=None, vo=None, zo=None, solarmotion=None)[source]¶
Interpolate past or future position of cluster and escaped stars
same as interpolate_orbit, but included for legacy purposes
- Parameters:
- clusterclass
StarCluster
- cluster_potclass
Galpy potential for host cluster that orbit is to be integrated in if None, assume a Plumme Potential
- potclass
galpy Potential that orbit is to be integrate in (default: None)
- tfinalfloat
final time (in cluster.units) to integrate orbit to (default: 12 Gyr)
- ntint
number of timesteps
- from_centrebool
intialize orbits from cluster’s exact centre instead of cluster’s position in galaxy (default :False)
- ro :float
galpy distance scale (Default: None)
- vofloat
galpy velocity scale (Default: None)
- zofloat
Sun’s distance above the Galactic plane (default: None)
- solarmotionfloat
array representing U,V,W of Sun (default: None)
- Returns:
- x,y,zfloat
interpolated positions of each star
- vx,vy,vzfloat
interpolated velocities of each star
- orbital_path(tfinal=0.1, nt=1000, pot=None, from_centre=False, skypath=False, initialize=False, ro=None, vo=None, zo=None, solarmotion=None, plot=False, **kwargs)[source]¶
Calculate the cluster’s orbital path
- Parameters:
- clusterclass
StarCluster
- tfinalfloat
final time (in cluster.units) to integrate orbit to (default: 0.1 Gyr)
- ntint
number of timesteps
- potclass
galpy Potential that orbit is to be integrate in (default: None)
- from_centrebool
genrate orbit from cluster’s exact centre instead of its assigned galactocentric coordinates (default: False)
- skypathbool
return sky coordinates instead of cartesian coordinates (default: False)
- initializebool
Initialize and return Orbit (default: False)
- ro :float
galpy distance scale (Default: None)
- vofloat
galpy velocity scale (Default: None)
- zofloat
Sun’s distance above the Galactic plane (default: None)
- solarmotionfloat
array representing U,V,W of Sun (default: None)
- plotbool
plot a snapshot of the cluster in galactocentric coordinates with the orbital path (defualt: False)
- Returns:
- tfloat
times for which path is provided
- x,y,zfloat
orbit positions
- vx,vy,vzfloat
orbit velocity
- oclass
galpy orbit (if initialize==True)
- History
- 2018 - Written - Webb (UofT)
- orbital_path_match(tfinal=0.1, nt=1000, pot=None, path=None, from_centre=False, skypath=False, to_path=False, do_full=False, ro=None, vo=None, zo=None, solarmotion=None, plot=False, projected=False, **kwargs)[source]¶
Match stars to a position along the orbital path of the cluster
- Parameters:
- clusterclass
StarCluster
- tfinalfloat
final time (in cluster.units) to integrate orbit to (default: 0.1 Gyr)
- ntint
number of timesteps
- potclass
galpy Potential that orbit is to be integrate in (default: None)
- patharray
array of (t,x,y,x,vx,vy,vz) corresponding to the tail path. If none path is calculated (default: None)
- from_centrebool
genrate orbit from cluster’s exact centre instead of its assigned galactocentric coordinates (default: False)
- skypathbool
return sky coordinates instead of cartesian coordinates (default: False) if True, projected is set to True
- to_pathbool
measure distance to the path itself instead of distance to central point along the path (default: False)
- do_fullbool
calculate dpath all at once in a single numpy array (can be memory intensive) (default:False)
- ro :float
galpy distance scale (Default: None)
- vofloat
galpy velocity scale (Default: None)
- zofloat
Sun’s distance above the Galactic plane (default: None)
- solarmotionfloat
array representing U,V,W of Sun (default: None)
- plotbool
plot a snapshot of the cluster in galactocentric coordinates with the orbital path (defualt: False)
- projectedbool
match to projected orbital path, which means matching just x and y coordinates or Ra and Dec coordinates (not z, or dist) (default:False)
- Returns:
- tstarfloat
orbital time associated with star
- dprogfloat
distance along the orbit to the progenitor
- dpath
distance to centre of the orbital path bin (Default) or the orbit path (to_path = True)
- orbits_interpolate(pot=None, tfinal=None, nt=1000, ro=None, vo=None, zo=None, solarmotion=None)[source]¶
Interpolate past or future position of stars within the cluster
same as interpolate_orbits, but kept for legacy purposes
- Parameters:
- clusterclass
StarCluster
- potclass
Galpy potential for host cluster that orbit is to be integrated in if None, assume a Plumme Potential
- tfinalfloat
final time (in cluster.units) to integrate orbit to (default: 12 Gyr)
- ntint
number of timesteps
- ro :float
galpy distance scale (Default: None)
- vofloat
galpy velocity scale (Default: None)
- zofloat
Sun’s distance above the Galactic plane (default: None)
- solarmotionfloat
array representing U,V,W of Sun (default: None)
- Returns:
- x,y,zfloat
interpolated positions of each star
- vx,vy,vzfloat
interpolated velocities of each star
- rcore(method='casertano', nneighbour=6, mfrac=0.1, projected=False, plot=False, **kwargs)[source]¶
Calculate core radius of the cluster – The default method (method=’casertano’) follows Casertano, S., Hut, P. 1985, ApJ, 298, 80 to find the core – An alternative metrhod (method==’isothermal’) assumes the cluster is an isothermal sphere the core radius is where density drops to 1/3 central value — For projected core radius, the core radius is where the surface density profile drops to 1/2 the central value — Note that the inner mass fraction of stars used to calculate central density is set by mfrac (default 0.1 = 10%)
- Parameters:
- clusterclass
StarCluster instance
- methodstring
method of calculating the core radius of a star cluster (default ‘casertano’)
- if method ==’casertano’:
- nneighbourint
number of neighbours for calculation local densities
- if method==’isothermal’:
- mfracfloat
inner mass fraction to be used to establish the central density
- projectedbool
use projected values (default: False)
- plotbool
plot the density profile and mark the core radius of the cluster (default: False)
- Returns
- ——-
- rcfloat
core radius
- Other Parameters:
- None
- relaxation_time(rad=None, coulomb=0.4, projected=None, method='spitzer')[source]¶
Calculate the relaxation time (Spitzer & Hart 1971) within a given radius of the cluster
Spitzer, L. Jr, Hart, M.H. 1971, ApJ, 164, 399 (Equation 5)
Need to adjust amplitude for different input units
Parameters¶
- clusterclass
StarCluster
- radfloat
radius within which to calculate the relaxation time (defult: cluster.rm)
- coulombfloat
Coulomb parameter (default: 0.4)
- projectedbool
use projected values (default: False)
- methodstr
choose between Spitzer & Hart 1971 and other methods (in development)
- Returns:
- trelaxfloat
relaxation time within radius rad
- reset_nbody_scale(mass=True, radii=True, rvirial=True, projected=None, **kwargs)[source]¶
Assign new conversions for real mass, size, and velocity to Nbody units
kwargs are passed to the virial_radius function. See the virial_radius documenation in functions.py
- Parameters:
- clusterclass
StarCluster instance
- massbool
find new mass conversion (default: True)
- radiibool
find new radius conversion (default: True)
- rvirialbool (default: True)
use virial radius to set conversion rate for radii as opposed to the approximation that rbar=4/3 rm
- projectedbool
use projected values to calculate radius conversion (default: False)
- Returns:
- zmbarfloat
mass conversion
- rbarfloat
radius conversion
- vbarfloat
velocity conversion
- tbarfloat
time conversion
- History:
- 2018 - Written - Webb (UofT)
- return_cluster(units0=None, origin0=None, rorder0=None, rorder_origin0=None, sortstars=False)[source]¶
return cluster to a specific combination of units and origin
- Parameters:
- clusterclass
StarCluster
- units0str
units that StarCluster will be changed to
- origin0str
origin that StarCluster will be changed to
- sortstarsbool
sort star by radius after coordinate change (default: False)
- Returns:
- None
- History:
- 2018 - Written - Webb (UofT)
- rlagrange(nlagrange=10, projected=None)[source]¶
Calculate lagrange radii of the cluster by mass
- Parameters:
- clusterclass
StarCluster
- nlagrangeint
number of lagrange radii bins (default: 10)
- mfracfloat
Exact masss fraction to calculate radius. Will supersede nlagrange if not None (default : None)
- projectedbool
calculate projected lagrange radii (default: False)
- Returns:
- rnfloat
lagrange radii
- rlimiting(pot=None, rgc=None, nrad=20, projected=False, plot=False, from_centre=False, verbose=False, **kwargs)[source]¶
Calculate limiting radius of the cluster
The limiting radius is defined to be where the cluster’s density reaches the local background density of the host galaxy
for cases where the cluster’s orbital parameters are not set, it is possible to manually set rgc which is assumed to be in kpc.
- Parameters:
- clusterclass
StarCluster
- potclass
GALPY potential used to calculate actions
- rgc
Manually set galactocentric distance in kpc at which the tidal radius is to be evaluated (default: None)
- zgcfloat
For non-spherically symmetric potentials, manually set distance in kpc above disk at which the tidal radius is to be evaluated. When set, rgc becomes radius in cylindrical coordinates (default: None)
- nradint
number of radial bins used to calculate density profile (Default: 20)
- projectedbool
use projected values (default: False)
- plotbool
plot the density profile and mark the limiting radius of the cluster (default: False)
- from_centrebool
calculate tidal radius based on location of cluster’s exact centre instead of its assigned galactocentric coordinates (default: False)
- verbosebool
- Returns
- ——-
- rlfloat
limiting radius
- Other Parameters:
- kwargsstr
key words for plotting
- rtidal(pot=None, rtiterate=0, rtconverge=0.9, indx=None, rgc=None, from_centre=False, plot=False, verbose=False, **kwargs)[source]¶
Calculate tidal radius of the cluster - The calculation uses Galpy (Bovy 2015_, which takes the formalism of Bertin & Varri 2008 to calculate the tidal radius – Bertin, G. & Varri, A.L. 2008, ApJ, 689, 1005 – Bovy J., 2015, ApJS, 216, 29 - riterate = 0 corresponds to a single calculation of the tidal radius based on the cluster’s mass (np.sum(cluster.m)) – Additional iterations take the mass within the previous iteration’s calculation of the tidal radius and calculates the tidal
radius again using the new mass until the change is less than 90%
for cases where the cluster’s orbital parameters are not set, it is possible to manually set rgc which is assumed to be in kpc.
- Parameters:
- clusterclass
StarCluster instance
- potclass
GALPY potential used to calculate tidal radius (default: None)
- rtiterateint
how many times to iterate on the calculation of r_t (default: 0)
- rtconvergefloat
criteria for tidal radius convergence within iterations (default 0.9)
- indxbool
subset of stars to use when calculate the tidal radius (default: None)
- rgcfloat
Manually set galactocentric distance in kpc at which the tidal radius is to be evaluated (default: None)
- zgcfloat
For non-spherically symmetric potentials, manually set distance in kpc above disk at which the tidal radius is to be evaluated. When set, rgc becomes radius in cylindrical coordinates (default: None)
- from_centrebool
calculate tidal radius based on location of cluster’s exact centre instead of its assigned galactocentric coordinates (default: False)
- plotbool
plot the x and y coordinates of stars and mark the tidal radius of the cluster (default: False)
- verbosebool
Print information about iterative calculation of rt
- Returns:
- rtfloat
tidal radius
- Other Parameters:
- kwargsstr
key words for plotting
- save_cluster()[source]¶
Save cluster’s units and origin
- Parameters:
- clusterclass
StarCluster
- Returns:
- cluster.units, cluster.originstr
units and origin of StarCluster
- History:
- 2018 - Written - Webb (UofT)
- sortstars(projected=True)[source]¶
sort stars based on radius
- Parameters:
- projectedbool
sort projected radii as well, but do not change self.projected (default: True)
- Returns:
- None
- subset(**kwargs)[source]¶
Generate a boolean array that corresponds to subset of star cluster members that meet a certain criteria
- Parameters:
- rmin/rmaxfloat
minimum and maximum stellar radii
- mmin/mmaxfloat
minimum and maximum stellar mass
- vmin/vmaxfloat
minimum and maximum stellar velocity
- emin/emaxfloat
minimum and maximum stellar energy
- kwmin/kwmaxint
minimum and maximum stellar type (kw)
- npopint
population number
- indxbool
user defined boolean array from which to extract the subset
- projectedbool
use projected values and constraints (default:False)
- Returns:
- indxbool
boolean array that is True for stars that meet the criteria
- tail_path(dt=0.1, no=1000, nt=100, ntail=100, pot=None, dmax=None, bintype='fix', from_centre=False, skypath=False, to_path=False, do_full=False, ro=None, vo=None, zo=None, solarmotion=None, plot=False, **kwargs)[source]¶
Calculate tail path +/- dt Gyr around the cluster
- Parameters:
- clusterclass
StarCluster
- dtfloat
timestep that StarCluster is to be moved to
- noint
number of timesteps for orbit integration (default:1000)
- ntint
number of points along the tail to set the tail spacing (default: 100)
- ntailint
number of points along the tail with roaming average (default: 1000)
- potclass
galpy Potential that orbit is to be integrate in (default: None)
- dmaxfloat
maximum distance (assumed to be same units as cluster) from orbital path to be included in generating tail path (default: None)
- bintypestr
type of binning for tail stars (default : ‘fix’)
- from_centrebool
genrate orbit from cluster’s exact centre instead of its assigned galactocentric coordinates (default: False)
- skypathbool
return sky coordinates instead of cartesian coordinates (default: False)
- to_pathbool
measure distance to the path itself instead of distance to central point along the path (default: False)
- do_fullbool
calculate dpath all at once in a single numpy array (can be memory intensive) (default:False)
- ro :float
galpy distance scale (Default: None)
- vofloat
galpy velocity scale (Default: None)
- zofloat
Sun’s distance above the Galactic plane (default: None)
- solarmotionfloat
array representing U,V,W of Sun (default: None)
- plotbool
plot a snapshot of the cluster in galactocentric coordinates with the orbital path (defualt: False)
- projectedbool
match to projected orbital path, which means matching just x and y coordinates or Ra and Dec coordinates (not z, or dist) (default:False)
- Returns:
- tfloat
times for which path is provided
- x,y,zfloat
tail path positions
- vx,vy,vzfloat
tail path velocities
- History
- 2018 - Written - Webb (UofT)
- 2019 - Implemented numpy array preallocation to minimize runtime - Nathaniel Starkman (UofT)
- tail_path_match(dt=0.1, no=1000, nt=100, ntail=100, pot=None, dmax=None, from_centre=False, to_path=False, do_full=False, ro=None, vo=None, zo=None, solarmotion=None, plot=False, **kwargs)[source]¶
Match stars to a position along the tail path of the cluster
- Parameters:
- clusterclass
StarCluster
- dtfloat
timestep that StarCluster is to be moved to
- noint
number of timesteps for orbit integration (default:1000)
- ntint
number of points along the tail to set the tail spacing (default: 100)
- ntailint
number of points along the tail with roaming average (default: 1000)
- potclass
galpy Potential that orbit is to be integrate in (default: None)
- patharray
array of (t,x,y,x,vx,vy,vz) corresponding to the tail path. If none path is calculated (default: None)
- from_centrebool
genrate orbit from cluster’s exact centre instead of its assigned galactocentric coordinates (default: False)
- skypathbool
return sky coordinates instead of cartesian coordinates (default: False) if True, projected is set to True
- to_pathbool
measure distance to the path itself instead of distance to central point along the path (default: False)
- do_fullbool
calculate dpath all at once in a single numpy array (can be memory intensive) (default:False)
- ro :float
galpy distance scale (Default: None)
- vofloat
galpy velocity scale (Default: None)
- zofloat
Sun’s distance above the Galactic plane (default: None)
- solarmotionfloat
array representing U,V,W of Sun (default: None)
- plotbool
plot a snapshot of the cluster in galactocentric coordinates with the orbital path (defualt: False)
- projectedbool
match to projected orbital path, which means matching just x and y coordinates or Ra and Dec coordinates (not z, or dist) (default:False)
- Returns:
- tstarfloat
orbital time associated with star
- dprogfloat
distance along the path to the progenitor
- dpath
distance to centre of the tail path bin (default) or the tail path (to_path = True)
- tapered_mass_function(mmin=None, mmax=None, nmass=10, rmin=None, rmax=None, vmin=None, vmax=None, emin=None, emax=None, kwmin=0, kwmax=1, npop=None, indx=None, projected=False, mcorr=None, plot=False, **kwargs)[source]¶
Find a tapered mass function over a given mass range
mass bins are set up so that there are an equal number of stars in each bin
functional form of the tapered mass function is taken from De Marchi, Paresce & Portegies Zwart 2010
Parameters¶
- clusterclass
StarCluster instance
- mmin/mmaxfloat
specific mass range
- nmass
number of mass bins used to calculate alpha
- rmin/rmax
specific radial range
- vmin/vmaxfloat
specific velocity range
- emin/emaxfloat
specific energy range
- kwmin/kwmaxint
specific stellar evolution type range
- npopint
population number
- indxbool
specific subset of stars
- projectedbool
use projected values (default: False)
- mcorrfloat
completeness correction for masses
- plotbool
plot the mass function
- Returns:
- m_meanfloat
mean mass in each bin
- m_histfloat
number of stars in each bin
- dmfloat
dN/dm of each bin
- alphafloat
power-law slope of the mass function (dN/dm ~ m^alpha)
- ealphafloat
error in alpha
- yalphafloat
y-intercept of fit to log(dN/dm) vs log(m)
- eyalphafloat
error in yalpha
- Other Parameters:
- kwargsstr
key words for plotting
- to_WDunits(analyze=True)[source]¶
Convert stellar positions/velocities, centre of mass, and orbital position and velocity to Walter Dehnen Units
- Parameters:
- clusterclass
StarCluster
- analyzebool
run analysis function (default: True)
- Returns:
- None
- History:
- 2022 - Written - Webb (UofT)
- to_amuse(analyze=True)[source]¶
Convert stellar positions/velocities, centre of mass, and orbital position and velocity to AMUSE Vector and Scalar Quantities
- Parameters:
- clusterclass
StarCluster
- analyzebool
run analysis function (default: True)
- Returns:
- None
- History:
- 2022 - Written - Webb (UofT)
- to_audays()[source]¶
Convert binary star semi-major axis and periods to solar radii and days
- Parameters:
- clusterclass
StarCluster
- Returns:
- None
- to_center(sortstars=True, centre_method=None)[source]¶
Shift coordinates such that origin is the centre of mass
this is the same to to_centre, but allows for alternate spelling
- Parameters:
- clusterclass
StarCluster
- sortstarsbool
sort star by radius after coordinate change (default: False)
- centre_methodstr
method for shifting coordinates to clustercentric coordinates (see to_cluster). (default: None)
- Returns:
- None
- History:
- 2018 - Written - Webb (UofT)
- to_centre(sortstars=True, centre_method=None)[source]¶
Shift coordinates such that origin is the centre of mass
- Parameters:
- clusterclass
StarCluster
- sortstarsbool
sort star by radius after coordinate change (default: False)
- centre_methodstr
method for shifting coordinates to clustercentric coordinates (see to_cluster). (default: None)
- Returns:
- None
- History:
- 2018 - Written - Webb (UofT)
- to_cluster(sortstars=True, centre_method=None)[source]¶
Shift coordinates to clustercentric reference frame
When units=’radec’ and origin=’sky’, the cluster will be shifted to clustercentric coordinates using either:
–centre_method=None: angular distance between each star’s RA/DEC and the RA/DEC of the cluster’s centre with proper motions directly subtracted off –centre_method=’orthographic’ , positions and velocities changed to orthnormal coordinates (Helmi et al. 2018) – centre_method=’VandeVen’ , positions and velocities changed to clustercentric coordinates using method outlined by Van de Ven et al. 2005
Note the the line of site positions and velocities will just have the galactocentric coordinates of the cluster
subtracted off. Be sure to set projected=True when making any calculations to use only x and y coordinates
- Parameters:
- clusterclass
StarCluster
- sortstarsbool
sort star by radius after coordinate change (default: False)
- centre_methodstr
method for shifting coordinates to clustercentric coordinates. (default: None)
- Returns:
- None
- History:
- 2018 - Written - Webb (UofT)
- to_galaxy(sortstars=True)[source]¶
Shift coordinates to galactocentric reference frame
- Parameters:
- clusterclass
StarCluster
- sortstarsbool
sort star by radius after coordinate change (default: False)
- Returns:
- None
- History:
- 2018 - Written - Webb (UofT)
- to_galpy(analyze=True)[source]¶
Convert stellar positions/velocities, centre of mass, and orbital position and velocity to galpy units
- Parameters:
- clusterclass
StarCluster
- analyzebool
run analysis function (default: True)
- Returns:
- None
- History:
- 2018 - Written - Webb (UofT)
- to_kpcgyr(analyze=True)[source]¶
Convert stellar positions/velocities, centre of mass, and orbital position and velocity to kpc and kpc/Gyr
- Parameters:
- clusterclass
StarCluster
- analyzebool
run analysis function (default: True)
- Returns:
- None
- History:
- 2022 - Written - Webb (UofT)
- to_kpckms(analyze=True)[source]¶
Convert stellar positions/velocities, centre of mass, and orbital position and velocity to kpc and km/s
- Parameters:
- clusterclass
StarCluster
- analyzebool
run analysis function (default: True)
- Returns:
- None
- History:
- 2018 - Written - Webb (UofT)
- to_nbody(analyze=True)[source]¶
Convert stellar positions/velocities, centre of mass, and orbital position and velocity to Nbody units
requires that cluster.zmbar, cluster.rbar, cluster.vbar are set (defaults are 1)
- Parameters:
- clusterclass
StarCluster
- analyzebool
run analysis function (default: True)
- Returns:
- None
- History:
- 2018 - Written - Webb (UofT)
- to_origin(origin, sortstars=True, centre_method=None)[source]¶
Shift cluster to origin as defined by keyword
- Parameters:
- clusterclass
StarCluster
- originstr
origin of coordinate system to be shifted to (centre, cluster, galaxy, sky)
- sortstarsbool
sort star by radius after coordinate change (default: False)
- centre_methodstr
method for shifting coordinates to clustercentric coordinates (see to_cluster). (default: None)
- Returns:
- None
- History:
- 2018 - Written - Webb (UofT)
- to_pckms(analyze=True)[source]¶
Convert stellar positions/velocities, centre of mass, and orbital position and velocity to pc and km/s
- Parameters:
- clusterclass
StarCluster
- analyzebool
run analysis function (default: True)
- Returns:
- None
- to_pcmyr(analyze=True)[source]¶
Convert stellar positions/velocities, centre of mass, and orbital position and velocity to pc and pc/Myr
- Parameters:
- clusterclass
StarCluster
- analyzebool
run analysis function (default: True)
- Returns:
- None
- History:
- 2022 - Written - Webb (UofT)
- to_radec(sortstars=True, centre_method=None, analyze=True)[source]¶
Convert to on-sky position, proper motion, and radial velocity of cluster
- Parameters:
- clusterclass
StarCluster
- sortstarsbool
sort star by radius after coordinate change (default: False)
- centre_methodstr
method for shifting coordinates to clustercentric coordinates (see to_cluster). (default: None)
- analyzebool
run analysis function (default: True)
- Returns:
- None
- History:
- 2018 - Written - Webb (UofT)
- to_sky(sortstars=True, centre_method=None)[source]¶
Calculate on-sky position, proper motion, and radial velocity of cluster
Also changes units to radec
- Parameters:
- clusterclass
StarCluster
- sortstarsbool
sort star by radius after coordinate change (default: False)
- centre_methodstr
method for shifting coordinates to clustercentric coordinates. (default: None)
- Returns:
- None
- History:
- 2018 - Written - Webb (UofT)
- to_sudays()[source]¶
- Convert binary star semi-major axis and periods to solar radii and days
Note: Masses will be converted to solar masses
- Parameters:
- clusterclass
StarCluster
- Returns:
- None
- to_tail()[source]¶
- Calculate positions and velocities of stars when rotated such that clusters velocity vector
points along x-axis
no change to coordinates in StarCluster
- Parameters:
- clusterclass
StarCluster
- Returns:
- x_tail,y_tail,z_tail,vx_tail,vy_tail,vz_tailfloat
rotated coordinates with cluster’s velocity vector point along x-axis
- History:
- 2018 - Written - Webb (UofT)
- to_units(units)[source]¶
Convert stellar positions/velocities, centre of mass, and orbital position and velocity to user defined units
- Parameters:
- clusterclass
StarCluster
- unitsstr
units to be converted to (nbody,galpy,pckms,kpckms,radec)
- Returns:
- None
- History:
- 2018 - Written - Webb (UofT)
- virial_radius(method='inverse_distance', full=True, H=70.0, Om=0.3, overdens=200.0, nrad=20, projected=None, plot=False, **kwargs)[source]¶
Calculate virial radius of the cluster - Virial radius is calculated using either: – the average inverse distance between particles, weighted by their masses (default) – the radius at which the density is equal to the critical density of the Universe at the redshift of the system, multiplied by an overdensity constant Parameters ———- cluster : class
StarCluster
- methodstr
method for calculating virial radius (default: ‘inverse_distance’, alternative: ‘critical_density’)
- Returns:
- rvfloat
virial radius
- Other Parameters:
- fullbool
Use Numba to calculate average inverse distance between stars (default:True)
- Hfloat
Hubble constant
- Omfloat
density of matter
- overdensfloat
overdensity constant
- projectedbool
calculate projected virial radius (default: False)
- plotbool
plot cluster density profile and illustrate virial radius calculation
- kwargsstr
key word arguments for plotting function
- virialize(qvir=0.5, specific=True, full=True, projected=None, softening=0.0)[source]¶
Adjust stellar velocities so cluster is in virial equilibrium
- Parameters:
- clusterclass
StarCluster
- qvirfloat
value you wish to virial parameter to be (default: 0.5)
- specificbool
find specific energies (default: True)
- full: bool
do full array of stars at once with numba (default: full_default)
- projectedbool
use projected values when calculating energies (default: False)
- softeningfloat
Plummer softening length in cluster.units (default: 0.0)
- Returns
- ——-
- qvfloat
scaling factor used to adjust velocities
- clustertools.cluster.cluster.overlap_cluster(cluster1, cluster2, tol=0.1, projected=False, return_cluster=True)[source]¶
Extract a sub population of stars from cluster1 that spatially overlaps with cluster2
- Parameters:
- cluster1StarCluster
cluster from which stars are to be extracted
- cluster2StarCluster
cluster from which overlap region is defined
- tol: float
tolerance parameter for how much regions need to overlap (default: 0.1)
- projectedbool
use projected values and constraints (default:False)
- return_cluster: bool
returns a sub_cluster if True, otherwise return the boolean array (default:True)
- Returns:
- StarCluster
- clustertools.cluster.cluster.sub_cluster(cluster, rmin=None, rmax=None, mmin=None, mmax=None, vmin=None, vmax=None, emin=None, emax=None, kwmin=0, kwmax=15, npop=None, indx=[None], projected=False, sortstars=True, reset_centre=False, reset_nbody=False, reset_nbody_mass=False, reset_nbody_radii=False, reset_rvirial=False, reset_projected=False, **kwargs)[source]¶
Extract a sub population of stars from a StarCluster
automatically moves cluster to centre of mass, so all constraints are in clustercentric coordinates and current StarCluster.units
- Parameters:
- rmin/rmaxfloat
minimum and maximum stellar radii
- mmin/mmaxfloat
minimum and maximum stellar mass
- vmin/vmaxfloat
minimum and maximum stellar velocity
- emin/emaxfloat
minimum and maximum stellar energy
- kwmin/kwmaxint
minimum and maximum stellar type (kw)
- npopint
population number
- indxbool
user defined boolean array from which to extract the subset
- projectedbool
use projected values and constraints (default:False)
- sortstars: bool
order stars by radius (default: True)
- reset_centrebool
re-calculate cluster centre after extraction (default:False)
- reset_nbodybool
reset nbody scaling factors (default:False)
- reset_nbody_massbool
find new nbody scaling mass (default:False)
- reset_nbody_radiibool
find new nbody scaling radius (default:False)
- reset_rvirialbool
use virial radius to find nbody scaling radius (default: True)
- reset_projectedbool
use projected radii to find nbody scaling radius (default: False)
- Returns:
- StarCluster