the main module

class pNbody.main.Nbody(p_name=None, pos=None, vel=None, mass=None, num=None, tpe=None, ftype='default', status='old', byteorder='little', pio='no', local=False, unitsfile=None, skipped_io_blocks=[], ptypes=None, verbose=0, arrays=None, **kws)

This is the reference Nbody class.

This is the constructor for the Nbody object. Optional arguments are:

p_namename of the file

in case of multiple files, files must be included in a list [“file1”,”file2”]

pos : positions (3xN array) vel : positions (3xN array) mass : positions (1x array) num : id of particles (1xN array) tpe : type of particles (1xN array)

ftype : type of input file (binary,ascii)

status‘old’open an old file

‘new’ : create a new object

byteorder : ‘little’ or ‘big’ pio : parallel io : ‘yes’ or ‘no’

local : True=local object, False=global object (paralellized) Not implemeted Yet

unitsfile : define the type of units

verbose : Verbose level (0: None, 1: general info, 2: details)

ptypes : type of particle to read

arrays : arrays to read

by default this class initialize the following variables :

self.p_name : name of the file(s) to read or write

self.pos : array of positions self.vel : array of velocities self.mass : array of masses self.num : array of id self.tpe : array of types

self.ftype : type of the file self.status : object status (‘old’ or ‘new’) self.byteorder : byter order (‘little’ or ‘big’) self.pio : parallel io (‘yes’ or ‘no’)

self.ptypes : type of particle to read

self.arrays : arrays to read

# new variables

self.nbody : local number of particles self.nbody_tot : total number of particles self.mass_tot : total mass self.npart : number of particles of each type self.npart_tot : total number of particles of each type self.spec_vars : dictionary of variables specific for the format used self.spec_vect : dictionary of vector specific for the format used

Accel(x, eps)

Return the acceleration at a given position, using the softening length eps.

x : position (array vector) eps : softening

CombiMap(*arg, **kw)

Return an image in form of a matrix (nx x ny float array). Contrary to ComputeMap, CombiMap compose different output of ComputeMap.

pos : position of particles (moment 0)

sr : dispertion in r (with respect to xp) svr : dispertion in vr

vxyr : mean velocity in the plane svxyr: dispertion in vxy

vtr : mean tangential velocity in the plane svtr : dispertion in vt

szr : ratio sigma z/sigma r

ComputeDensityAndHsml(pos=None, Hsml=None, DesNumNgb=None, MaxNumNgbDeviation=None, Tree=None)

Compute Density and Hsml (for a specific place)

ComputeHisto(bins, mode, space)

Compute and histogram

ComputeMap(*arg, **kw)

Return an image in form of a matrix (nx x ny float array)

obs : position of observer x0 : eye position xp : focal position alpha : angle of the head view : ‘xy’ ‘xz’ ‘yz’

eye : ‘right’ ‘left’ dist_eye : distance between eyes

mode : mode of map space : pos or vel

persp : ‘on’ ‘off’ clip : (near,far) size : (maxx,maxy)

cut : ‘yes’ ‘no’

frsp : factor for rsp shape : shape of the map

ComputeMeanHisto(bins, mode1, space)

Compute the mean map of an observable.

ComputeMeanMap(*arg, **kw)

Compute the mean map of an observable.

ComputeObjectMap(*arg, **kw)
      • IN DEVELOPPEMENT : allow to draw an object like a box, a grid… * * *

Return an image in form of a matrix (nx x ny float array)

obs : position of observer x0 : eye position xp : focal position alpha : angle of the head view : ‘xy’ ‘xz’ ‘yz’

eye : ‘right’ ‘left’ dist_eye : distance between eyes

mode : mode of map space : pos or vel

persp : ‘on’ ‘off’ clip : (near,far) size : (maxx,maxy)

cut : ‘yes’ ‘no’

frsp : factor for rsp shape : shape of the map

ComputeRsp(Nngb=5)

Compute a value for self.rsp base on the distance to the Nngb neighbours This function is not paralellised.

ComputeSigmaHisto(bins, mode1, mode2, space)

Compute the histogram of an observable.

ComputeSigmaMap(*arg, **kw)

Compute the sigma map of an observable.

ComputeSph(DesNumNgb=None, MaxNumNgbDeviation=None, Tree=None)

Compute self.Density and self.Hsml using sph approximation

Ekin()

Return the total kinetic energy

Epot(eps)

Return the total potential energy using the softening length eps.

eps : softening

WARNING : THIS FUNCTION DO NOT WORK IN MPI MODE

ExchangeParticles()

Exchange particles betwee procs, using peano-hilbert decomposition computed in ptree

Get_Velocities_From_AdaptativeSpherical_Grid(select=None, eps=0.1, n=1000, UseTree=True, Tree=None, phi=None, ErrTolTheta=0.5)

Computes velocities using the jeans equation in spherical coordinates. An adaptative grid is set automatically.

Get_Velocities_From_Cylindrical_Grid(select='disk', disk=('gas', 'disk'), eps=0.1, nR=32, nz=32, nt=2, Rmax=100, zmin=- 10, zmax=10, params=[None, None, None], UseTree=True, Tree=None, Phi=None, ErrTolTheta=0.5, AdaptativeSoftenning=False, g=None, gm=None, NoDispertion=False)

Computes velocities using the jeans equation in cylindrical coordinates.

Get_Velocities_From_Spherical_Grid(select=None, eps=0.1, nr=128, rmax=100.0, beta=0, UseTree=True, Tree=None, phi=None, ErrTolTheta=0.5, g=None, gm=None, NoDispertion=False, omega=None, vmax_limit=True)

Computes velocities using the jeans equation in spherical coordinates.

Get_Velocities_From_Virial_Approximation(select=None, vf=1.0, eps=0.1, UseTree=True, Tree=None, ErrTolTheta=0.5)

This routine does not work ! Do not use it, or check !

IntegrateUsingRK(tstart=0, dt=1, dt0=1e-05, epsx=1e-13, epsv=1e-13)

Integrate the equation of motion using RK78 integrator.

tstart : initial time dt : interval time dt0 : inital dt epsx : position precision epsv : velocity precision

tmin,tmax,dt,dtout,epsx,epsv,filename

L()

Return the angular momentum in x,y,z of all particles. The output is an 3xn float array.

Ltot()

Return the total angular momentum. The output is an 3x1 float array.

Map(*arg, **kw)

Return 2 final images (float and int)

Mr_Spherical(nr=25, rmin=0, rmax=50)

Return the mass inside radius r (supposing a spherical density distribution). The output is 2 n x 1 float arrays.

nr : number of bins (size of the output) rmin : minimal radius (this must be zero, instead it is wrong…) rmax : maximal radius

Pot(x, eps)

Return the potential at a given position, using the softening length eps.

x : position (array vector) eps : softening

R()

Return a 1xn float array that corresponds to the projected distance from the center of each particle.

SendAllToAll()

Send all particles to all nodes at the end of the day, all nodes have the same nbody object

SphEvaluate(val, pos=None, vel=None, hsml=None, DesNumNgb=None, MaxNumNgbDeviation=None, Tree=None)

Return an sph evaluation of the variable var

TreeAccel(pos, eps, Tree=None)

Return the acceleration at a given position, using the softening length eps and using a tree.

pos : position (array vector) eps : softening Tree: gravitational tree if already computed

WARNING : this function do not work in parallel

TreePot(pos, eps, Tree=None)

Return the potential at a given position, using the softening length eps and using a tree.

pos : position (array vector) eps : softening Tree: gravitational tree if already computed

WARNING : this function do not work in parallel

VR()

Return the norm of the radial (cylindrical) velocities of particles The output is an 1xn float array.

VT()

Return the norm of the tangiancial (cylindrical) velocities of particles The output is an 1xn float array.

Vphi()

Return the norm of the phi component of the velocities of particles The output is an 1xn float array.

Vr()

Return the norm of the radial (spherical) velocities of particles The output is an 1xn float array.

Vrxy()

Return the radial velocities (in the z=0 plane) of particles The output is an 3xn float array.

Vt()

Return the norm of the tangential velocities of particles The output is an 1xn float array.

Vtheta()

Return the norm of the theta component of the velocities of particles The output is an 1xn float array.

Vtxy()

Return the tangential velocities (in the z=0 plane) of particles The output is an 3xn float array.

Vz()

Return a 1xn float array containing z velocity

align(axis, mode='a', sgn='+', fact=None)

Rotate the object in order to align the axis ‘axis’ with the z axis.

axis : [x,y,z] mode : ‘p’ : only position

‘v’ : only velocities ‘a’ : both (default)

sgn‘+’normal rotation

‘-’ : reverse sense of rotation

fact : int : factor to increase the angle

align2(axis1=[1, 0, 0], axis2=[0, 0, 1], point=[0, 0, 0])

Rotate the object in order to align the axis ‘axis’ with the z axis.

axis1 : [x,y,z] axis2 : [x,y,z] point : [x,y,z]

align_with_main_axis(mode='a')

Rotate the object in order to align it’s major axis with the axis of its inertial tensor.

mode‘p’only position

‘v’ : only velocities ‘a’ : both (default)

append(solf, do_not_sort=False, do_init_num=True, do_not_change_num=False)

Add to the current N-body object, particles form the N-body object “new”.

solf : Nbody object

cart2sph(pos=None)

Transform carthesian coodinates x,y,z into spherical coordinates r,p,t Return a 3xn float array.

check_arrays()

check if the array contains special values like NaN or Inf

check_ftype()

check the file format

clone_empty(local=False)

create a new empty object using the same class

cm()

Return the mass center of the model. The output is an 3x1 float array.

cmcenter()

Move the N-body object in order to center the mass center at the origin.

cv()

Return the center of the velocities of the model. The output is an 3x1 float array.

cvcenter()

Center the center of velocities at the origin.

debug(msg)

print warning

dens(r=None, nb=25, rm=50)

Return the number density at radius r (supposing a spherical density distribution). If r is not specified, it is computed with nb and rm. The output is an n x 1 float array.

!!! This routine do not use masses !!!

r : radius nb : number of bins (size of the output) rm : maximal radius

display(*arg, **kw)

Display the model

dmodes(nr=32, nm=16, rm=32)

Compute the density modes of a model

nm = 16 : number of modes nr = 32 : number of radius rm = 50 : max radius

return

r : the radius used m : the modes computed m1 : the matrix of the amplitude m2 : the matrix of the phases

dump(aname=None)

Dump the array aname to a file

dumpToHDF5(aname=None)

Dump the array aname to a file

dv_mean()

Return the average relative speed between particles.

dx_mean()

Return the average distance between particles.

ekin()

Return the total specific kinetic energy

epot(eps)

Return the total specific potential energy using the softening length eps.

eps : softening

WARNING : THIS FUNCTION DO NOT WORK IN MPI MODE

expose(obs, eye=None, dist_eye=None, foc=None, space='pos', pos=None, vel=None)

Rotate and translate the object in order to be seen as if the observer was in x0, looking at a point in xp.

obs : observer matrix eye : ‘right’ or ‘left’ dist_eye : distance between eyes (separation = angle) space : pos or vel foc : focal

extend_format()

Extend format with format file (e.g. gh5) and extensions (config/extension)

find_format(default)

Test the default format and if not good, test a few. :returns: format name

find_vars()

This function return a list of variables defined in the current object

gather_mass()

Gather in a unique array all mass of all nodes.

gather_num()

Gather in a unique array all num of all nodes.

gather_pos()

Gather in a unique array all positions of all nodes.

gather_vec(vec)

Gather in a unique array all vectors vec of all nodes.

gather_vel()

Gather in a unique array all velocites of all nodes.

getAccelerationInCylindricalGrid(eps, z, Rmax, nr=32, nt=32, UseTree=False)

Compute the Acceleration in cells of a cylindrical grid

getNumberParticlesInCylindricalGrid(Rmax, nr=32, nt=32)

Compute the number of particles in cells of a cylindrical grid

getPotentialInCylindricalGrid(eps, z, Rmax, nr=32, nt=32, UseTree=False)

Compute the potential in cells of a cylindrical grid

getRadialVelocityDispersionInCylindricalGrid(Rmax, nr=32, nt=32)

Compute the radial velocity dispersion in cells of a cylindrical grid

getRadiusInCylindricalGrid(z, Rmax, nr=32, nt=32)

Compute the radius in cells of a cylindrical grid

getSurfaceDensityInCylindricalGrid(Rmax, nr=32, nt=32)

Compute the surface density in cells of a cylindrical grid

getTree(force_computation=False, ErrTolTheta=0.8)

Return a Tree object

get_arrays_memory_footprint()

Return the memory footprint of each array

get_default_arrays_props()

return default arrays properties for the class

get_default_spec_array()

return specific array default values for the class

get_default_spec_vars()

return specific variables default values for the class

get_excluded_extension()

Return a list of file to avoid when extending the default class.

get_ext_methods()

get the list of extensions methods together with their origin

get_format_file()

return the format file

get_ftype(ftype='binary')

get the current used format

get_histocenter(rbox=50, nb=500, center=[0, 0, 0], fromBoxSize=False)

Return the position of the higher density region in x,y,z (not good) found by the function “histocenter”.

rbox : size of the box nb : number of bins in each dimension fromBoxSize : compute the histograms from the boxsize information center : center of the box

get_histocenter2(rbox=50, nb=64)

Return the position of the higher density region in x,y,z (not good) found by the function “histocenter”.

rbox : size of the box nb : number of bins in each dimension

get_histocenter3D(rbox=50, nb=256, center=[0, 0, 0], fromBoxSize=False)

Return the position of the higher density region in x,y,z (not good) found by the function “histocenter”.

rbox : size of the box nb : number of bins in each dimension fromBoxSize : compute the histograms from the boxsize information center : center of the box

get_list_of_array()

Return the list of numpy vectors of size nbody.

get_list_of_method()

Return the list of instance methods (functions).

get_list_of_vars()

Get the list of vars that are linked to the model

get_mass_tot()

Return the total mass of system.

get_mxntpe()

Return the max number of type for this format

get_nbody()

Return the local number of particles.

get_nbody_tot()

Return the total number of particles.

get_npart()

Return the local number of particles of each types, based on the variable tpe

get_npart_all(npart_tot, NTask)

From npart_tot, the total number of particles per type, return npart_per_proc, an array where each element corresponds to the value of npart of each process.

get_npart_and_npart_all(npart)

From npart (usually read for the header of a file), compute :

npart : number of particles in each type npart_tot : total number of particles in each type npart_all : npart for each process.

get_npart_tot()

Return the total number of particles of each types.

get_ns()

Return in an array the number of particles of each node.

get_ntype()

return the number of paticles types

get_num()

Compute the num variable in order to be consistent with particles types

get_potcenter(ptype=None, eps=0.1, force=False)

Return the centre defined as the position of the particle having the minimal potential

get_read_fcts()

returns the functions needed to read a snapshot file.

get_rotation_matrix_to_align_with_main_axis()

Get the rotation matrix used to rotate the object in order to align it’s main axis with the axis of its inertial tensor.

get_rsp_approximation(DesNumNgb=None, MaxNumNgbDeviation=None, Tree=None)

Return an aproximation of rsp, based on the tree.

get_write_fcts()

returns the functions needed to write a snapshot file.

getindex(num)

Return an array of index of a particle from its specific number id. The array is empty if no particle corresponds to the specific number id.

num : Id of the particle

has_array(name)

Return true if the object pNbody has an array called self.name

has_var(name)

Return true if the object pNbody has a variable called self.name

hdcenter()

Move the N-body object in order to center the higher density region found.

histocenter(rbox=50, nb=500, fromBoxSize=False)

Move the N-body object in order to center the higher density region found near the mass center. The higher density region is determined with density histograms.

rbox : box dimension, where to compute the histograms nb : number of bins for the histograms fromBoxSize : compute the histograms from the boxsize information

histocenter2(rbox=50, nb=64)

Move the N-body object in order to center the higher density region found near the mass center. The higher density region is determined whith density histograms.

rbox : box dimension, where to compute the histograms nb : number of bins for the histograms

histocenter3D(rbox=50, nb=500, center=[0, 0, 0], fromBoxSize=False)

Move the N-body object in order to center the higher density region found near the mass center. The higher density region is determined with a 3D histograms.

rbox : box dimension, where to compute the histograms nb : number of bins for the histograms center : center of the box fromBoxSize : compute the histograms from the boxsize information

histovel(nb=100, vmin=None, vmax=None, mode='n')

Return or plot the histrogram of the norm of velocities or of the radial velocities.

The output is a list (r,h) of 2 nx1 float arrays, where r is the radius and h the values of the histogram.

nb : number of bins (size of the output) vmax : maximum velocity vmin : minimum velocity mode : ‘n’ (norme of the velocities)

‘r’ (radial velocities)

import_check_ftype(filename)

Import check_spec_ftype from the format file.

inertial_tensor()

Return the inertial tensor.

info()

Write info

init()

Initialize normal and specific class variables

init_units()

This function is responsible for the units initialization.

It will create :

self.unitsparameters

that contains parameters like
  • the hydrogen mass fraction,

  • the metalicity ionisation flag

  • the adiabatic index

and

self.localsystem_of_units

a UnitSystem object that really defines the system of units in the Nbody object. It uses the values :

UnitLength_in_cm UnitMass_in_g UnitVelocity_in_cm_per_s

All physical values computed in pNbody should use self.localsystem_of_units to be converted in other units. self.unitsparameters is usefull if other parameters needs to be known, like the adiabatic index, etc.

l()

Return the specific angular momentum in x,y,z of all particles. The output is an 3xn float array.

load(name=None, ptypes=None, force=False)

Load array from the file. This function relays on the info stored in self._AM.arrays_props Here, we assume that npart is known and correct

name : the array name ptypes : the list of particles types to read force : for reloading the array

loadFromHDF5(name=None, ptypes=None, force=False)

Load array from the hdf5 file. This function relays on the info stored self.arrays_props Here, we assume that npart is known and correct

name : the array name ptypes : the list of particles types to read force : for reloading the array

ltot()

Return the specific total angular momentum. The output is an 3x1 float array.

make_default_vars_global()

Make specific variables global

mdens(r=None, nb=25, rm=50)

Return the density at radius r (supposing a spherical density distribution). If r is not specified, it is computed with nb and rm. The output is an n x 1 float array.

r : radius nb : number of bins (size of the output) rm : maximal radius

memory_info()

Write info on memory size of the current object (only counting arrays size)

message(msg, verbosity=1, color=None)

print message

minert()

Return the diagonal of the intertial momentum.

mr(r=None, nb=25, rm=50)

Return the mass inside radius r (supposing a spherical density distribution). If r is not specified, it is computed with nb and rm. The output is an n x 1 float array.

r : radius nb : number of bins (size of the output) rm : maximal radius

msdens(r=None, nb=25, rm=50)

Return the mass surface density at radius r. If r is not specified, it is computed with nb and rm. The output is an nx1 float array.

r : radius nb : number of bins (size of the output) rm : maximal radius

nodes_info()

Write info on nodes

object_info()

Write class(object) info

open_and_read(name, readfct)

open and read file name

name : name of the input readfct : function used to read the file

open_and_write(name, writefct)

Open and write file

name : name of the output writefct : function used to write the file

phi_xy()

Return a 1xn float array that corresponds to the azimuth in cylindrical coordinate of each particle.

phi_xyz()

Return a 1xn float array that corresponds to the azimuth in spherical coordinate of each particle.

potcenter(ptype=None, eps=0.1)

center the model according to the potential center (position of the particle having the minimal potential).

ptype : particle type on which the potential is computed eps : gravitational softening length

print_filenames()

Print files names

r(center=None)

Return a 1xn float array that corresponds to the distance from the center of each particle.

read()

Read the particle file(s)

read_arrays(ptypes=None)

Read arrays that must be loaded.

read_num(name)

Read a num file

name : name of the input

real_numngb(num)

number of particles wor wich r<h

rebox(boxsize=None, mode=None, axis='xyz')

Translate the positions of the object in order that all particles being contained in a box of size boxsize.

boxsizesize of the box

if boxsize is not defined, we try first to see if self.boxsize is defined.

modetype of reboxing

None : -> [0,boxsize] centred : -> [-boxsize/2,boxsize/2] [x,y,z] :

redistribute()

This function redistribute particles amoung all nodes in order to have a similar number of particles per nodes

reduc(n, mass=False)

Return an N-body object that contain a fraction 1/n of particles.

n : inverse of the fraction of particule to be returned

rename(p_name=None)

Rename the files

p_name : new name(s)

rotate(angle=0, axis=[1, 0, 0], point=[0, 0, 0], mode='a')

Rotate the positions and/or the velocities of the object around a specific axis defined by a vector and an point.

angle : rotation angle in radian axis : direction of the axis point : center of the rotation

mode‘p’rotate only position

‘v’ : rotate only velocities ‘a’ : rotate both (default)

rotateR(R, mode='a')

Rotate the model using the matrix R

mode‘p’only position

‘v’ : only velocities ‘a’ : both (default)

rotate_old(angle=0, mode='a', axis='x')

Rotate the positions and/or the velocities of the object around a specific axis.

angle : rotation angle in radian axis : ‘x’ : around x

‘y’ : around y ‘z’ : around z

: [x,y,z] : around this axis

mode‘p’rotate only position

‘v’ : rotate only velocities ‘a’ : rotate both (default)

rxy()

Return a 1xn float array that corresponds to the projected distance from the center of each particle.

rxyz(center=None)

Return a 1xn float array that corresponds to the distance from the center of each particle.

sdens(r=None, nb=25, rm=50)

Return the surface density at radius r. If r is not specified, it is computed with nb and rm. The output is an nx1 float array.

!!! This routine do not uses masses !!!

r : radius nb : number of bins (size of the output) rm : maximal radius

selectc(c, local=False)

Return an N-body object that contain only particles where the corresponding value in c is not zero. c is a nx1 Nbody array.

c : the condition vector local : local selection (True) or global selection (False)

selecti(i, local=False)

Return an N-body object that contain only particles having their index (not id) in i.

i : vector containing indexes local : local selection (True) or global selection (False)

selectp(lst=None, file=None, reject=False, local=False, from_num=True)

Return an N-body object that contain only particles with specific number id.

The list of id’s is given either by lst (nx1 int array) or by the name (“file”) of a file containing the list of id’s.

lst : vector list (integer)

reject : True/False : if True, reject particles in lst (default = False) local : local selection (True) or global selection (False)

frum_numif True, use self.num to select particules

if False, use arange(self.nbody)

set_filenames(p_name, pio=None)

Set the local and global names

p_name : new name(s) pio : ‘yes’ or ‘no’

set_ftype(ftype='binary')

Change the type of the file

ftype : type of the file

set_local_system_of_units(params='default', UnitLength_in_cm=None, UnitVelocity_in_cm_per_s=None, UnitMass_in_g=None, unitparameterfile=None, gadgetparameterfile=None)

Set local system of units using UnitLength_in_cm,UnitVelocity_in_cm_per_s,UnitMass_in_g

  1. if nothing is given, we use self.unitsparameters to obtain these values

2a) if UnitLength_in_cm, UnitVelocity_in_cm_per_s, UnitMass_in_g are given, we use them

2b) if UnitLength_in_cm, UnitVelocity_in_cm_per_s, UnitMass_in_g are given in a dictionary

  1. if unitparameterfile is given we read the parameters from the file (units parameter format)

  2. if gadgetparameterfile is given we read the parameters from the file (gadget param format)

set_npart(npart)

Set the local number of particles of each types. This function modifies the variable self.tpe

set_parameters(params)

Set parameters for the class

set_pio(pio)

Set parallel input/output or not io

pio : ‘yes’ or ‘no’

set_tpe(tpe)

Set all particles to the type tpe

set_unitsparameters(unitsparams)

Set units parameters for the class.

show(*arg, **kw)

Display the model this is an alias to display

sigma(r=None, nb=25.0, rm=50.0)

Return the 3 velocity dispersion (in cylindrical coordinates) and the mean azimuthal velocity curve. If r is not specified, it is computed with nb and rm.

The output is a list (r,sr,st,sz,mt) of 5 $n imes 1$ float arrays, where r is the radius, sr the radial velocity dispersion, st, the azimuthal velocity dispersion, sz, the vertical velocity dispersion and mt, the mean azimuthal velocity curve.

!!! This routine works only if particles have equal masses !!!

r : radius where to compute the values nb : number of bins (size of the output) rm : maximal radius

return : r,sr,st,sz,mt

sigma_vz(r=None, nb=25, rm=50)

Return the vertical dispertion in z at radius r. If r is not specified, it is computed with nb and rm. The output is an nx1 float array.

r : radius nb : number of bins (size of the output) rm : maximal radius

sigma_z(r=None, nb=25, rm=50)

Return the vertical dispertion in z at radius r. If r is not specified, it is computed with nb and rm. The output is an nx1 float array.

r : radius nb : number of bins (size of the output) rm : maximal radius

size()

Estimate the model size, using the inertial momentum

sort(vec=None, local=False)

sort particles according to the vector vec

vec : vector on which to sort (default=self.num)

sort_type(local=False)

Contrary to sort, this fonction sort particles respecting their type.

spec_info()

Write specific info

sph2cart(pos=None)

Transform spherical coordinates r,p,t into carthesian coodinates x,y,z Return a 3xn float array.

sphericalvel2pos()

Replace the position by velocities in spherical coordinates

spin(omega=None, L=None, j=None, E=None)

Spin the object with angular velocity “omega” (rigid rotation). Omega is a 1 x 3 array object

If L (total angular momentum) is explicitely given, compute Omega from L (1 x 3 array object).

omega : angular speed (array vector) L : desired angular momentum j : desired energy fraction in rotation E : Total energy (without rotation)

sub(n1=0, n2=None)

Return an N-body object that have particles whith indicies in the range [n1:n2].

n1 : number of the first particule n2 : number of the last particule

Note : the first particle is 0

take(vec=None, local=False)

extract particles according to the vector vec

vec : vector (default=self.num)

theta_xyz()

Return a 1xn float array that corresponds to the elevation angle in spherical coordinate of each particle.

tork(acc)

Return the total tork on the system due to the force acting on each particle (acc). The output is an 3xn float array.

acc : 3xn float array

translate(dx, mode='p')

Translate the positions or the velocities of the object.

dx : shift (array vector) mode : ‘p’ : translate positions

‘v’ : translate velocities

usual_numngb(num)

usual way to compute the number of neighbors

v_sigma()

Return the norm of the velocity dispersions.

vel2pos()

Replace the position by velocities in Cartesian coordinates

vel_cart2cyl()

Transform velocities in carthesian coordinates vx,vy,vz into cylindrical coodinates vr,vz,vz. Pos is the position of particles in cart. coord. Vel is the velocity in cart. coord. Return a 3xn float array.

vel_cyl2cart(pos=None, vel=None)

Transform velocities in cylindrical coordinates vr,vt,vz into carthesian coodinates vx,vy,vz. Pos is the position of particles in cart. coord. Vel is the velocity in cylindrical coord. Return a 3xn float array.

vn()

Return a 1xn float array that corresponds to the norm of velocities

vrxyz()

Return a 1xn float array that corresponds to the radial velocity in spherical system

vx()

Return a 1xn float array containing x velocity

vy()

Return a 1xn float array containing y velocity

vz()

Return a 1xn float array containing z velocity

warning(msg, verbosity=1)

print warning

weighted_numngb(num)

num = particle where to compute weighted_numngb see Springel 05

write(name=None)

Write the particle file(s)

write_num(name)

Write a num file

name : name of the output

x(center=None)

Return a 1xn float array containing x coordinate

x_sigma()

Return the norm of the position dispersions.

y(center=None)

Return a 1xn float array containing y coordinate

z(center=None)

Return a 1xn float array containing z coordinate

zmodes(nr=32, nm=16, rm=32)

Compute the vertical modes of a model

nm = 16 : number of modes nr = 32 : number of radius rm = 50 : max radius

return

r : the radius used m : the modes computed m1 : the matrix of the amplitude m2 : the matrix of the phases

zprof(z=None, r=2.5, dr=0.5, nb=25, zm=5.0)

Return the z-profile in a vector for a given radius

!!! This routine works only if particles have equal masses !!!

z : bins in z (optional) r : radius of the cut dr : width in r of the cut nb : number of bins (size of the output) zm : maximal height