Elson
- cosmic.sample.cmc.elson.M_enclosed(r, gamma, rho_0)
Compute the mass enclosed in an Elson profile at radius r with slope gamma, central concentration rho_0, and assumed scale factor a = 1
In practice, this is only used to sample the positions of stars, so rho_0 is just picked to normalize the distribution (i.e. rho_0 s.t. M_enclosed(rmax) = 1)
- cosmic.sample.cmc.elson.draw_r_vr_vt(N=100000, r_max=300, gamma=4)
Draw random velocities and positions from the Elson profile.
N = number of stars r_max = maximum number of virial radii for the farthest star gamma = steepness function of the profile
Note that gamma=4 is the Plummer profile
returns (vr,vt,r) in G=M_cluster=1 units
N.B. if you’re getting a “Expect x to not have duplicates” error with a large gamma, use a smaller r_max.
- cosmic.sample.cmc.elson.find_rmax_vir(r_max, gamma)
This is a little tricky: because the virial radius of the Elson profile depends on the maximum radius, if the profile is very flat (e.g. gamma~2) then you need a large maximum radius to get a large number of virial radius in there. This basically finds r such that r / rVir(r) = r_max. In other word, how far out do we need to go to have r_max number of virial radii in the sample.
- cosmic.sample.cmc.elson.find_sigma_sqr(r, r_max_cluster, gamma)
Find the 1D velocity dispersion at a given radius r using one of the spherial Jeans equations (and assuming velocity is isotropic)
- cosmic.sample.cmc.elson.get_positions(N, r_max_cluster, gamma)
This one’s easy: the mass enclosed function is just the CDF of the mass density, so just invert that, and you’ve got positions.
- cosmic.sample.cmc.elson.get_velocities(r, r_max_cluster, gamma)
The correct way to generate velocities: generate the distribution function from rho, then use rejection sampling.
returns (vr,vt) with same length as r
- cosmic.sample.cmc.elson.get_velocities_old(r, r_max_cluster, gamma)
Uses the spherical Jeans functions to sample the velocity dispersion for the cluster at different radii, then draws a random, isotropic velocity for each star
NOTE: this gives you correct velocity dispersion and produces a cluster in virial equilibrium, but it’s not strictly speaking correct (the distribution should have some kurtosis, in addition to getting the variance of the Gaussian correct). You can see the disagreement in the tail of the distributions when sampling a Plummer sphere. This is what mcluster does…
Use the new get_velocities function, which generates a distribution function directly from rho and samples velocities from it
This is kept around because it’s super useful to sample from when the rejection sampling gets too slow at the last X samples (e.g. gamma ~ 10)
returns (vr,vt) with same length as r
- cosmic.sample.cmc.elson.get_velocities_plummer(r, r_max_cluster)
The correct way to generate velocities
this is a special case for gamma=4, which is the plummer sphere (and has an analytic distribution function)
returns (vr,vt) with same length as r
- cosmic.sample.cmc.elson.phi_r(r, gamma, rho_0)
Compute the gravitational potential of an Elson profile at radius r with slope gamma, central concentration rho_0, and assumed scale factor a = 1
Needed to compute the escape speed (for resampling stars)
- cosmic.sample.cmc.elson.rho_r(r, gamma, rho_0)
Compute the density of the Elson profile at radius r Best to use the same normalized rho_0 from M_enclosed
- cosmic.sample.cmc.elson.scale_pos_and_vel(r, vr, vt)
Scale the positions and velocities to be in N-body units If we add binaries we’ll do this again in initialcmctable.py
takes r, vr, and vt as input
returns (r,vr,vt) scaled to Henon units
- cosmic.sample.cmc.elson.virial_radius_analytic(gamma, r_max)
Virial radius is best calculated directly, since rmax may be pretty far from infinity. Directly integrate 4*pi*r*rho*m_enclosed from 0 to rMax to get the binding energy, then just divide 0.5 by that.