Particle Attributes

Particle attributes are defined in the amuse.datamodel.particle_attributes. These attributes can be accessed in two ways, on the particle(s) or as a function imported from the module. When accessed on the particles, the first parameter (usually called `particles`) must not be given.

To access these functions directly do:

from amuse.datamodel.particle_attributes import *
from amuse.lab import new_plummer_sphere

particles = new_plummer_sphere(100)
print kinetic_enery(particles)

To access these functions as an attribute do:

from amuse.lab import new_plummer_sphere

particles = new_plummer_sphere(100)

print particles.kinetic_energy()
amuse.datamodel.particle_attributes.LagrangianRadii(stars, unit_converter='auto', mf=[0.01, 0.02, 0.05, 0.1, 0.2, 0.5, 0.75, 0.9, 1], cm=None, number_of_neighbours=7, reuse_hop=False, hop=<amuse.datamodel.particle_attributes.HopContainer object>)

Calculate lagrangian radii. Output is radii, mass fraction

>>> import numpy
>>> from amuse.ic.plummer import new_plummer_sphere
>>> numpy.random.seed(1234)
>>> parts=new_plummer_sphere(100)
>>> lr,mf=parts.LagrangianRadii()
>>> print lr[5]
0.856966667972 length
class amuse.datamodel.particle_attributes.MassSegregationRatioResults(mass_segregation_ratio, uncertainty)
property mass_segregation_ratio

Alias for field number 0

property uncertainty

Alias for field number 1

amuse.datamodel.particle_attributes.Qparameter(parts, distfunc=None)

Calculates the minimum spanning tree Q parameter (Cartwright & Whitworth 2004) for a projection of the particle set.

Parameters

distfunc – distfunc is the distance function which can be used to select

the projection plane.

amuse.datamodel.particle_attributes.bound_subset(particles, tidal_radius=None, unit_converter=None, density_weighting_power=2, smoothing_length_squared=quantity<zero>, G=quantity<6.67428e-11 m**3 * kg**-1 * s**-2>, core=None, reuse_hop=False, hop=<amuse.datamodel.particle_attributes.HopContainer object>, gravity_code=None)

find the particles bound to the cluster. Returns a subset of bound particles.

Parameters
  • tidal_radius – particles beyond this are considered not bound

  • unit_converter – Required if the particles are in SI units

  • density_weighting_power – Particle properties are weighted by density to this power

  • smooting_length_squared – the smoothing length for gravity.

  • G – gravitational constant, need to be changed for particles in different units systems

  • core – (optional) core of the cluster

>>> from amuse.ic.plummer import new_plummer_model
>>> from amuse.units import nbody_system
>>> plum=new_plummer_model(100)
>>> print len(plum.bound_subset(G=nbody_system.G))
100
>>> plum[0].velocity*=100  
>>> plum[0].position*=100  
>>> print len(plum.bound_subset(G=nbody_system.G))
99
amuse.datamodel.particle_attributes.box_counting_dimension(particles)

Computes the box-counting dimension, a measure of the fractal dimension of a set of points. The measure is based on counting the number of boxes required to cover the set, within a regular grid of cubic, equal-size boxes, for varying box sizes.

amuse.datamodel.particle_attributes.center_of_mass(particles)

Returns the center of mass of the particles set. The center of mass is defined as the average of the positions of the particles, weighted by their masses.

>>> from amuse.datamodel import Particles
>>> particles = Particles(2)
>>> particles.x = [-1.0, 1.0] | units.m
>>> particles.y = [0.0, 0.0] | units.m
>>> particles.z = [0.0, 0.0] | units.m
>>> particles.mass = [1.0, 1.0] | units.kg
>>> particles.center_of_mass()
quantity<[0.0, 0.0, 0.0] m>
amuse.datamodel.particle_attributes.center_of_mass_velocity(particles)

Returns the center of mass velocity of the particles set. The center of mass velocity is defined as the average of the velocities of the particles, weighted by their masses.

>>> from amuse.datamodel import Particles
>>> particles = Particles(2)
>>> particles.vx = [-1.0, 1.0] | units.ms
>>> particles.vy = [0.0, 0.0] | units.ms
>>> particles.vz = [0.0, 0.0] | units.ms
>>> particles.mass = [1.0, 1.0] | units.kg
>>> particles.center_of_mass_velocity()
quantity<[0.0, 0.0, 0.0] m * s**-1>
amuse.datamodel.particle_attributes.connected_components(parts, threshold=None, distfunc=None, verbose=False)

return a list of connected component subsets of particles, connected if the distfunc is smaller than the threshold.

Parameters
  • threshold – value of the threshold. Must have consistent units with distfunc

  • distfunc – distance or weight function. Must have consistent units with threshold

amuse.datamodel.particle_attributes.correlation_dimension(particles, max_array_length=10000000)

Computes the correlation dimension, a measure of the fractal dimension of a set of points. The measure is based on counting the number of pairs with a mutual distance less than ‘eps’, for varying values of ‘eps’.

amuse.datamodel.particle_attributes.densitycentre_coreradius_coredens(particles, unit_converter=None, number_of_neighbours=7, reuse_hop=False, hop=<amuse.datamodel.particle_attributes.HopContainer object>)

calculate position of the density centre, coreradius and coredensity

>>> import numpy
>>> from amuse.ic.plummer import new_plummer_sphere
>>> numpy.random.seed(1234)
>>> particles=new_plummer_sphere(100)
>>> pos,coreradius,coredens=particles.densitycentre_coreradius_coredens()
>>> print coreradius
0.404120092331 length
amuse.datamodel.particle_attributes.distances_squared(particles, other_particles)

Returns the distance squared from each particle in this set to each of the particles in the other set.

Parameters

other_particles – the particles to which the distance squared is calculated

>>> from amuse.datamodel import Particles
>>> field_particles = Particles(2)
>>> field_particles.x = [0.0, 2.0] | units.m
>>> field_particles.y = [0.0, 0.0] | units.m
>>> field_particles.z = [0.0, 0.0] | units.m
>>> particles = Particles(3)
>>> particles.x = [1.0, 3.0, 4.0] | units.m
>>> particles.y = [0.0, 0.0, 0.0] | units.m
>>> particles.z = [0.0, 0.0, 0.0] | units.m
>>> particles.distances_squared(field_particles)
quantity<[[1.0, 1.0], [9.0, 1.0], [16.0, 4.0]] m**2>
amuse.datamodel.particle_attributes.dynamical_timescale(particles, mass_fraction=None, G=quantity<6.67428e-11 m**3 * kg**-1 * s**-2>)

Compute the dynamical (i.e. free-fall) timescale of the particles set. This is the time it would take for a pressureless homogeneous sphere of this size and average density to collapse. If ‘mass_fraction’ is supplied, only the inner particles are considered in the computation of the size of the sphere. For example, ‘mass_fraction=0.95’ ignores the positions of the outer particles comprising 5% by mass (useful for density profiles with long tails).

amuse.datamodel.particle_attributes.find_closest_particle_to(particles, x, y, z)

return closest particle to x,y,z position

>>> from amuse.datamodel import Particles
>>> particles = Particles(2)
>>> particles.x = [0.0, 2.0] | units.m
>>> particles.y = [0.0, 0.0] | units.m
>>> particles.z = [0.0, 0.0] | units.m
>>> print particles.find_closest_particle_to( -1 | units.m,0.| units.m,0.| units.m).x
0.0 m
amuse.datamodel.particle_attributes.get_binaries(particles, hardness=10, G=quantity<6.67428e-11 m**3 * kg**-1 * s**-2>)

returns the binaries in a particleset. binaries are selected according to a hardness criterion [hardness=10] This function returns the binaries as a list of i,j particles. Triple detection is not done.

>>> from amuse import datamodel
>>> m = [1,1,1] | units.MSun
>>> x = [-1,1,0] | units.AU
>>> y = [0,0,1000] | units.AU
>>> z = [0,0,0] | units.AU
>>> vx = [0,0,0] | units.kms
>>> vy = [1.,-1.,0] | units.kms
>>> vz = [0,0,0] | units.kms
>>> particles = datamodel.create_particle_set( mass=m,x=x,y=y,z=z,vx=vx,vy=vy,vz=vz )
>>> binaries = particles.get_binaries()
>>> print len(binaries)
1
amuse.datamodel.particle_attributes.kinetic_energy(particles)

Returns the total kinetic energy of the particles in the particles set.

>>> from amuse.datamodel import Particles
>>> particles = Particles(2)
>>> particles.vx = [-1.0, 1.0] | units.ms
>>> particles.vy = [0.0, 0.0] | units.ms
>>> particles.vz = [0.0, 0.0] | units.ms
>>> particles.mass = [1.0, 1.0] | units.kg
>>> particles.kinetic_energy()
quantity<1.0 m**2 * kg * s**-2>
amuse.datamodel.particle_attributes.mass_segregation_Gini_coefficient(particles, unit_converter=None, density_weighting_power=2, core=None, reuse_hop=False, hop=<amuse.datamodel.particle_attributes.HopContainer object>)

Converse & Stahler 2008 Gini coefficient for cluster.

Parameters
  • unit_converter – Required if the particles are in SI units

  • density_weighting_power – Particle properties are weighted by density to this power

  • core – (optional) core of the cluster

>>> import numpy
>>> from amuse.ic.plummer import new_plummer_model
>>> from amuse.units import nbody_system
>>> plum=new_plummer_model(100)
>>> index=plum.position.lengths_squared().argmin()
>>> plum.mass=0|nbody_system.mass
>>> plum[index].mass=1|nbody_system.mass
>>> print plum.mass_segregation_Gini_coefficient()
1.0
amuse.datamodel.particle_attributes.mass_segregation_from_nearest_neighbour(particles, number_of_particles=None, fraction_of_particles=0.01, number_of_random_sets=None, also_compute_uncertainty=False)

Calculates the mass segregation ratio based on the average inverse distance to nearest neighbours.

  1. Determine average inverse distances to the nearest neighbour for ‘number_of_random_sets’ sets of ‘number_of_particles’ random stars; mean_idnn[number_of_random_sets]

  2. Determine the average inverse distance to the nearest neighbour of the ‘number_of_particles’ most massive stars; mean_idnn_massive

  3. Determine with what statistical significance l_massive differs from l_norm: MSR = mean_idnn_massive / <mean_idnn> +/ - sigma(mean_idnn) / <mean_idnn>

Parameters
  • number_of_particles – the number of particles in each (random and massive) sample

  • number_of_random_sets – the number of randomly selected subsets for which the average inverse distances to the nearest neighbour are calculated

  • also_compute_uncertainty – if True, a namedtuple is returned with (MSR, sigma)

amuse.datamodel.particle_attributes.mass_segregation_ratio(particles, number_of_particles=20, number_of_random_sets=50, also_compute_uncertainty=False)

Calculates the mass segregation ratio (Allison et al. 2009, MNRAS 395 1449).

  1. Determine the length of the minimum spanning tree (MST) of the ‘number_of_particles’ most massive stars; l_massive

  2. Determine the average length of the MST of ‘number_of_random_sets’ sets of ‘number_of_particles’ random stars; l_norm

  3. Determine with what statistical significance l_massive differs from l_norm: MSR = (l_norm / l_massive) +/- (sigma_norm / l_massive)

Parameters
  • number_of_particles – the number of most massive stars for the MST for l_massive

  • number_of_random_sets – the number of randomly selected subsets for which the MST is calculated to determine l_norm

  • also_compute_uncertainty – if True, a namedtuple is returned with (MSR, sigma)

amuse.datamodel.particle_attributes.minimum_spanning_tree_length(particles)

Calculates the length of the minimum spanning tree (MST) of a set of particles using David Eppstein’s Python implemention of Kruskal’s algorithm.

amuse.datamodel.particle_attributes.moment_of_inertia(particles)

Returns the total moment of inertia (about the Z axis) of the particle set.

amuse.datamodel.particle_attributes.move_to_center(particles)

Shift positions and velocities of the particles such that their center of mass (velocity) is centered on the origin.

Implemented as:

particles.position -= particles.center_of_mass() particles.velocity -= particles.center_of_mass_velocity()

amuse.datamodel.particle_attributes.nearest_neighbour(particles, neighbours=None, max_array_length=10000000)

Returns the nearest neighbour of each particle in this set. If the ‘neighbours’ particle set is supplied, the search is performed on the neighbours set, for each particle in the orignal set. Otherwise the nearest neighbour in the same set is searched.

Parameters

neighbours – the particle set in which to search for the nearest neighbour (optional)

>>> from amuse.datamodel import Particles
>>> particles = Particles(3)
>>> particles.x = [1.0, 3.0, 4.0] | units.m
>>> particles.y = [0.0, 0.0, 0.0] | units.m
>>> particles.z = [0.0, 0.0, 0.0] | units.m
>>> particles.nearest_neighbour().x
quantity<[3.0, 4.0, 3.0] m>
>>> field_particles = Particles(2)
>>> field_particles.x = [0.0, 2.5] | units.m
>>> field_particles.y = [0.0, 0.0] | units.m
>>> field_particles.z = [0.0, 0.0] | units.m
>>> particles.nearest_neighbour(field_particles).x
quantity<[0.0, 2.5, 2.5] m>
amuse.datamodel.particle_attributes.new_particle_from_cluster_core(particles, unit_converter=None, density_weighting_power=2, cm=None, reuse_hop=False, hop=<amuse.datamodel.particle_attributes.HopContainer object>)

Uses Hop to find the density centre (core) of a particle distribution and stores the properties of this core on a particle: position, velocity, (core) radius and (core) density.

Particles are assigned weights that depend on the density (as determined by Hop) to a certain power. The default weighting power is 2, which is most commonly used. Set density_weighting_power to 1 in order to get the original weighting of Casertano & Hut (1985, ApJ, 298, 80).

Parameters
  • unit_converter – Required if the particles are in SI units

  • density_weighting_power – Particle properties are weighted by density to this power

amuse.datamodel.particle_attributes.particle_potential(set, particle, smoothing_length_squared=quantity<zero>, G=quantity<6.67428e-11 m**3 * kg**-1 * s**-2>)

Returns the potential at the position of the particle.

Parameters
  • smooting_length_squared – gravitational softening, added to every distance**2.

  • G – gravitational constant, need to be changed for particles in different units systems

>>> from amuse.datamodel import Particles
>>> particles = Particles(2)
>>> particles.x = [0.0, 1.0] | units.m
>>> particles.y = [0.0, 0.0] | units.m
>>> particles.z = [0.0, 0.0] | units.m
>>> particles.mass = [1.0, 1.0] | units.kg
>>> particles[1].potential()
quantity<-6.67428e-11 m**2 * s**-2>
amuse.datamodel.particle_attributes.particle_specific_kinetic_energy(set, particle)

Returns the specific kinetic energy of the particle.

>>> from amuse.datamodel import Particles
>>> particles = Particles(2)
>>> particles.vx = [0.0, 1.0] | units.ms
>>> particles.vy = [0.0, 0.0] | units.ms
>>> particles.vz = [0.0, 0.0] | units.ms
>>> particles.mass = [1.0, 1.0] | units.kg
>>> particles[1].specific_kinetic_energy()
quantity<0.5 m**2 * s**-2>
amuse.datamodel.particle_attributes.particleset_potential(particles, smoothing_length_squared=quantity<zero>, G=quantity<6.67428e-11 m**3 * kg**-1 * s**-2>, gravity_code=None, block_size=0)

Returns the potential at the position of each particle in the set.

Parameters
  • smooting_length_squared – gravitational softening, added to every distance**2.

  • G – gravitational constant, need to be changed for particles in different units systems

>>> from amuse.datamodel import Particles
>>> particles = Particles(2)
>>> particles.x = [0.0, 1.0] | units.m
>>> particles.y = [0.0, 0.0] | units.m
>>> particles.z = [0.0, 0.0] | units.m
>>> particles.mass = [1.0, 1.0] | units.kg
>>> particles.potential()
quantity<[-6.67428e-11, -6.67428e-11] m**2 * s**-2>
amuse.datamodel.particle_attributes.potential_energy(particles, smoothing_length_squared=quantity<zero>, G=quantity<6.67428e-11 m**3 * kg**-1 * s**-2>)

Returns the total potential energy of the particles in the particles set.

Parameters
  • smooting_length_squared – gravitational softening, added to every distance**2.

  • G – gravitational constant, need to be changed for particles in different units systems

>>> from amuse.datamodel import Particles
>>> particles = Particles(2)
>>> particles.x = [0.0, 1.0] | units.m
>>> particles.y = [0.0, 0.0] | units.m
>>> particles.z = [0.0, 0.0] | units.m
>>> particles.mass = [1.0, 1.0] | units.kg
>>> particles.potential_energy()
quantity<-6.67428e-11 m**2 * kg * s**-2>
amuse.datamodel.particle_attributes.potential_energy_in_field(particles, field_particles, smoothing_length_squared=quantity<zero>, G=quantity<6.67428e-11 m**3 * kg**-1 * s**-2>, just_potential=False)

Returns the total potential energy of the particles associated with an external gravitational field, which is represented by the field_particles.

Parameters
  • field_particles – the external field consists of these (i.e. potential energy is calculated relative to the field particles)

  • smooting_length_squared – gravitational softening, added to every distance**2.

  • G – gravitational constant, need to be changed for particles in different units systems

>>> from amuse.datamodel import Particles
>>> field_particles = Particles(2)
>>> field_particles.x = [0.0, 2.0] | units.m
>>> field_particles.y = [0.0, 0.0] | units.m
>>> field_particles.z = [0.0, 0.0] | units.m
>>> field_particles.mass = [1.0, 1.0] | units.kg
>>> particles = Particles(2)
>>> particles.x = [1.0, 3.0] | units.m
>>> particles.y = [0.0, 0.0] | units.m
>>> particles.z = [0.0, 0.0] | units.m
>>> particles.mass = [1.0, 1.0] | units.kg
>>> particles.potential_energy_in_field(field_particles)
quantity<-2.22476e-10 m**2 * kg * s**-2>
amuse.datamodel.particle_attributes.scale_to_standard(particles, convert_nbody=None, smoothing_length_squared=quantity<zero>, virial_ratio=0.5)

Scale the particles to a standard NBODY model with G=1, total_mass=1, and virial_radius=1 (or potential_energy=-0.5). In virial equilibrium (virial_ratio=0.5, default) the kinetic_energy=0.25 and the velocity_dispersion=1/sqrt(2).

Parameters
  • convert_nbody – the scaling is in nbody units, when the particles are in si units a convert_nbody is needed

  • smoothing_length_squared – needed for calculating the potential energy correctly.

  • virial_ratio – scale velocities to Q=K/|U|, (kinetic/potential energy); Q = virial_ratio > 0.5: supervirial, will expand Q = virial_ratio < 0.5: subvirial, will collapse

amuse.datamodel.particle_attributes.specific_kinetic_energy(particles)

Returns the specific kinetic energy of each particle in the set.

>>> from amuse.datamodel import Particles
>>> particles = Particles(2)
>>> particles.vx = [1.0, 1.0] | units.ms
>>> particles.vy = [0.0, 0.0] | units.ms
>>> particles.vz = [0.0, 0.0] | units.ms
>>> particles.mass = [1.0, 1.0] | units.kg
>>> particles.specific_kinetic_energy()
quantity<[0.5, 0.5] m**2 * s**-2>
amuse.datamodel.particle_attributes.thermal_energy(particles)

Returns the total internal energy of the (gas) particles in the particles set.

>>> from amuse.datamodel import Particles
>>> particles = Particles(2)
>>> particles.u = [0.5, 0.5] | units.ms**2
>>> particles.mass = [1.0, 1.0] | units.kg
>>> particles.thermal_energy()
quantity<1.0 m**2 * kg * s**-2>
amuse.datamodel.particle_attributes.total_angular_momentum(particles)

Returns the total angular momentum of the particles set.

>>> from amuse.datamodel import Particles
>>> particles = Particles(2)
>>> particles.x = [-1.0, 1.0] | units.m
>>> particles.y = [0.0, 0.0] | units.m
>>> particles.z = [0.0, 0.0] | units.m
>>> particles.vx = [0.0, 0.0] | units.ms
>>> particles.vy = [-1.0, 1.0] | units.ms
>>> particles.vz = [0.0, 0.0] | units.ms
>>> particles.mass = [1.0, .5] | units.kg
>>> particles.total_angular_momentum()
quantity<[0.0, 0.0, 1.5] m**2 * kg * s**-1>
amuse.datamodel.particle_attributes.total_mass(particles)

Returns the total mass of the particles set.

>>> from amuse.datamodel import Particles
>>> particles = Particles(3)
>>> particles.mass = [1.0, 2.0, 3.0] | units.kg
>>> particles.total_mass()
quantity<6.0 kg>
amuse.datamodel.particle_attributes.total_momentum(particles)

Returns the total momentum of the particles set.

>>> from amuse.datamodel import Particles
>>> particles = Particles(2)
>>> particles.vx = [-1.0, 1.0] | units.ms
>>> particles.vy = [0.0, 0.0] | units.ms
>>> particles.vz = [0.0, 0.0] | units.ms
>>> particles.mass = [1.0, 1.0] | units.kg
>>> particles.total_momentum()
quantity<[0.0, 0.0, 0.0] m * kg * s**-1>
amuse.datamodel.particle_attributes.total_radius(particles)

Returns the total radius (maximum distance from center) of the particles set.

>>> from amuse.datamodel import Particles
>>> particles = Particles(3)
>>> particles.mass = [1.0, 2.0, 3.0] | units.kg
>>> particles.position = [0.0, 0.0, 0.0] | units.m
>>> particles.x = [0.0, 3.0, 6.0] | units.m
>>> particles.total_radius()
quantity<4.0 m>
amuse.datamodel.particle_attributes.velocity_diff_squared(particles, field_particles)

Returns the total potential energy of the particles in the particles set.

Parameters

field_particles – the external field consists of these (i.e. potential energy is calculated relative to the field particles)

>>> from amuse.datamodel import Particles
>>> field_particles = Particles(2)
>>> field_particles.vx = [0.0, 2.0] | units.m
>>> field_particles.vy = [0.0, 0.0] | units.m
>>> field_particles.vz = [0.0, 0.0] | units.m
>>> particles = Particles(3)
>>> particles.vx = [1.0, 3.0, 4] | units.m
>>> particles.vy = [0.0, 0.0, 0.0] | units.m
>>> particles.vz = [0.0, 0.0, 0.0] | units.m
>>> velocity_diff_squared(particles, field_particles)
quantity<[[1.0, 1.0], [9.0, 1.0], [16.0, 4.0]] m**2>
amuse.datamodel.particle_attributes.virial_radius(particles)

Returns the virial radius of the particles set. The virial radius is the inverse of the average inverse distance between particles, weighted by their masses.

>>> from amuse.datamodel import Particles
>>> particles = Particles(2)
>>> particles.x = [-1.0, 1.0] | units.m
>>> particles.y = [0.0, 0.0] | units.m
>>> particles.z = [0.0, 0.0] | units.m
>>> particles.mass = [1.0, 1.0] | units.kg
>>> particles.virial_radius()
quantity<4.0 m>