.. _star_analysis:

Star Particle Analysis
======================
.. sectionauthor:: Stephen Skory <sskory@physics.ucsd.edu>
.. versionadded:: 1.6

This document describes tools in yt for analyzing star particles.
The Star Formation Rate tool bins stars by time to produce star formation
statistics over several metrics.
A synthetic flux spectrum and a spectral energy density plot can be calculated
with the Spectrum tool.

.. _star_formation_rate:

Star Formation Rate
-------------------

This tool can calculate various star formation statistics binned over time.
As input it can accept either a yt ``data_source``, such as a region or 
sphere (see :ref:`available-objects`), or arrays containing the data for
the stars you wish to analyze.

This example will analyze all the stars in the volume:

.. code-block:: python

  from yt.mods import *
  from yt.extensions.StarAnalysis import *
  pf = load("data0030")
  all = pf.h.all()
  sfr = StarFormationRate(pf, data_source=all)

or just a small part of the volume:

.. code-block:: python

  from yt.mods import *
  from yt.extensions.StarAnalysis import *
  pf = load("data0030")
  sp = p.h.sphere([0.5,0.5,0.5], 0.05)
  sfr = StarFormationRate(pf, data_source=all)

If the stars to be analyzed cannot be defined by a data_source, arrays can be
passed. In this case, the units for the ``star_mass`` must be in Msun,
the ``star_creation_time`` in code units, and the volume must be specified
in mpc as a float
(but it doesn't have to be correct depending on which statistic is important).

.. code-block:: python

  from yt.mods import *
  from yt.extensions.StarAnalysis import *
  pf = load("data0030")
  re = pf.h.region([0.5,0.5,0.5], [0.4,0.5,0.6], [0.5,0.6,0.7])
  # This puts the particle data for *all* the particles in the region re
  # into the arrays sm and ct.
  sm = re["ParticleMassMsun"]
  ct = re["creation_time"]
  # First pick out only stars.
  stars = (ct > 0)
  ct = ct[stars]
  sm = sm[stars]
  # Pick out only old stars using Numpy array fancy indexing.
  # 100 is a time in code units.
  sm_old = sm[ct < 100]
  ct_old = ct[ct < 100]
  sfr = StarFormationRate(pf, star_mass=sm_old, star_creation_time=ct_old,
  volume=re.volume('mpc'))

To output the data to a text file, use the command ``.write_out``:

.. code-block:: python

  sfr.write_out(name="StarFormationRate.out")

In the file ``StarFormationRate.out``, there are seven columns of data:

  1. Time (yrs)
  2. Look-back time (yrs)
  3. Redshift
  4. Star formation rate in this bin per year (Msol/yr)
  5. Star formation rate in this bin per year per Mpc**3 (Msol/yr/Mpc**3)
  6. Stars formed in this time bin (Msol)
  7. Cumulative stars formed up to this time bin (Msol)

The output is easily plotted. This is a plot for some test data (that may or may not 
correspond to anything physical) using columns #2 and #4 for the x and y
axes, respectively:

.. image:: SFR.png
   :width: 640
   :height: 480


.. _synthetic_spectrum:

Synthetic Spectrum Generator
----------------------------

Based on code generously provided by Kentaro Nagamine <kn@physics.unlv.edu>,
this will generate a synthetic spectrum for the stars using the publicly-available
tables of Bruzual & Charlot (hereafter B&C). Please see their `2003 paper
<http://adsabs.harvard.edu/abs/2003MNRAS.344.1000B>`_ for more information
and the `main data
distribution page <http://www.cida.ve/~bruzual/bc2003>`_ for the original data.
Based on the mass, age and metallicity of each star, a cumulative spectrum is
generated and can be output in two ways, either raw, or as a spectral
energy distribution.

This analysis toolkit reads in the B&C data from HDF5 files that have been
converted from the original ASCII files (available at the link above). The
HDF5 files are one-quarter the size of the ASCII files, and greatly reduce
the time required to read the data off disk. The HDF5 files are available from
the main yt website `here <http://yt.enzotools.org/files/bc03/>`_.
Both the Salpeter and Chabrier models have been converted,
and it is simplest to download all the files to the same location.
Please read the original B&C sources for information on the differences between
the models.

In order to analyze stars, first the Bruzual & Charlot data tables need to be
read in from disk. This is accomplished by initializing ``SpectrumBuilder`` and
specifying the location of the HDF5 files with the ``bcdir`` parameter.
The models are chosen with the ``model`` parameter, which is either
"chabrier" or "salpeter".

.. code-block:: python

  from yt.mods import *
  from yt.extensions.StarAnalysis import *
  pf = load("data0030")
  spec = SpectrumBuilder(pf, bcdir="/home/username/bc/", model="chabrier")

In order to analyze a set of stars, use the ``calculate_spectrum`` command.
It accepts either a ``data_source``, or a set of arrays with the star 
information. Continuing from the above example:

.. code-block:: python

  re = pf.h.region([0.5,0.5,0.5], [0.4,0.5,0.6], [0.5,0.6,0.7])
  spec.calculate_spectrum(data_source=pf)

If a subset of stars are desired, call it like this. ``star_mass`` is in units
of Msun, ``star_creation_time`` and ``star_metallicity_fraction`` in code
units.

.. code-block:: python

  re = pf.h.region([0.5,0.5,0.5], [0.4,0.5,0.6], [0.5,0.6,0.7])
  # This puts the particle data for *all* the particles in the region re
  # into the arrays sm, ct and metal.
  sm = re["ParticleMassMsun"]
  ct = re["creation_time"]
  metal = re["metallicity_fraction"]
  # First pick out only stars.
  stars = (ct > 0)
  ct = ct[stars]
  sm = sm[stars]
  metal = metal[stars]
  # Pick out only old stars using Numpy array fancy indexing.
  # 100 is a time in code units.
  sm_old = sm[ct < 100]
  ct_old = ct[ct < 100]
  metal_old = metal[ct < 100]
  spec.calculate_spectrum(star_mass=sm_old, star_creation_time=ct_old,
  star_metallicity_fraction=metal_old)

Alternatively, when using either a ``data_source`` or individual arrays,
the option ``star_metallicity_constant`` can be specified to force all the
stars to have the same metallicity. If arrays are being used, the
``star_metallicity_fraction`` array need not be specified.

.. code-block:: python

  # Make all the stars have solar metallicity.
  spec.calculate_spectrum(data_source=pf, star_metallicity_constant=0.02)

There are two ways to write out the data once the spectrum has been calculated.
The command ``write_out`` outputs two columns of data:

  1. Wavelength (Angstrom)
  2. Flux (Luminosity per unit wavelength, L_sun Ang^-1, L_sun = 3.826 * 10^33 ergs s^-1.)

and can be called simply, specifying the output file:

.. code-block:: python

  spec.write_out(name="spec.out")

The other way is to output a spectral energy density plot. Along with the
``name`` parameter, this command can also take the ``flux_norm`` option,
which is the wavelength in Angstroms of the flux to normalize the 
distribution to. The default is 5200 Angstroms. This command outputs the data
in two columns:

  1. Wavelength (Angstrom)
  2. Relative flux normalized to the flux at *flux_norm*.

.. code-block:: python

  spec.write_out_SED(name="SED.out", flux_norm=5200)

Below is an example of an absurd SED for universe-old stars all with 
solar metallicity at a redshift of zero. Note that even in this example,
a ``pf`` is required.

.. code-block:: python

  from yt.mods import *
  from yt.extensions.StarAnalysis import *
  pf = load("data0030")
  spec = SpectrumBuilder(pf, bcdir="/home/user/bc", model="chabrier")
  sm = na.ones(100)
  ct = na.zeros(100)
  spec.calculate_spectrum(star_mass=sm, star_creation_time=ct, star_metallicity_constant=0.02)
  spec.write_out_SED('SED.out')

And the plot:

.. image:: SED.png
   :width: 640
   :height: 480

Iterate Over a Number of Haloes
^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^

In this example below, the haloes for a dataset are found, and the SED is calculated
and written out for each.

.. code-block:: python

  from yt.mods import *
  from yt.extensions.StarAnalysis import *
  pf = load("data0030")
  # Find all the haloes, and include star particles.
  haloes = HaloFinder(pf, dm_only=False)
  # Set up the spectrum builder.
  spec = SpectrumBuilder(pf, bcdir="/home/user/bc", model="salpeter")
  # Iterate over the haloes.
  for halo in haloes:
      # Get the pertinent arrays.
      ct = halo["creation_time"]
      sm = halo["ParticleMassMsun"]
      metal = halo["metallicity_fraction"]
      # Select just the stars.
      stars = (ct > 0)
      ct = ct[stars]
      sm = sm[stars]
      metal = metal[stars]
      # Calculate the spectrum.
      spec.calculate_spectrum(star_mass=sm, star_creation_time=ct,
      star_metallicity_fraction=metal)
      # Write out the SED using the default flux normalization.
      spec.write_out_SED(name="halo%05d.out" % halo.id)

