hyperion.model.Model

class hyperion.model.Model(name=None)

Adding sources

add_source(source)
add_point_source(*args, **kwargs)
add_spherical_source(*args, **kwargs)
add_external_spherical_source(*args, **kwargs)
add_external_box_source(*args, **kwargs)
add_map_source(*args, **kwargs)
add_plane_parallel_source(*args, **kwargs)

Setting the grid

set_grid(grid)
set_cartesian_grid(x_wall, y_wall, z_wall)
set_cylindrical_polar_grid(w_wall, z_wall, …)
set_spherical_polar_grid(r_wall, t_wall, p_wall)
set_octree_grid(x, y, z, dx, dy, dz, refined)
set_amr_grid(description)

Setting quantities

add_density_grid(density, dust[, …]) Add a density grid to the model

Images/SEDs

add_peeled_images([sed, image]) Define a set of (peeled) images/SEDs.
add_binned_images([sed, image]) Define a set of (binned) images/SEDs.

Configuration

set_seed(seed) Set the seed for the random number generation
set_propagation_check_frequency(frequency) Set how often to check that the photon is in the right cell
set_monochromatic(monochromatic[, …]) Set whether to do the radiation transfer at specific wavelengths.
set_minimum_temperature(temperature) Set the minimum temperature for the dust
set_minimum_specific_energy(specific_energy) Set the minimum specific energy for the dust
set_specific_energy_type(specific_energy_type) Set whether to use the specific specific energy as an initial value or an additional component at each iteration.
set_n_initial_iterations(n_iter) Set the number of initial iterations for computing the specific energy in each cell.
set_raytracing(raytracing) Set whether to use raytracing for the non-scattered flux
set_max_interactions(inter_max[, warn]) Set the maximum number of interactions a photon can have.
set_max_reabsorptions(reabs_max[, warn]) Set the maximum number of successive reabsorptions by a source that a photon can have.
set_pda(pda) Set whether to use the Partial Diffusion Approximation (PDA)
set_mrw(mrw[, gamma, inter_max, warn]) Set whether to use the Modified Random Walk (MRW) approximation
set_convergence(convergence[, percentile, …]) Set whether to check for convergence over the initial iterations
set_kill_on_absorb(kill_on_absorb) Set whether to kill absorbed photons
set_forced_first_interaction(…[, …]) Set whether to ensure that photons scatter at least once before escaping the grid.
set_output_bytes(io_bytes) Set whether to output physical quantity arrays in 32-bit or 64-bit
set_sample_sources_evenly(sample_sources_evenly) If set to ‘True’, sample evenly from all sources and apply probability weight based on relative luminosities.
set_enforce_energy_range(enforce) Set how to deal with cells that have specific energy rates that are below or above that provided in the mean opacities and emissivities.
set_copy_input(copy) Set whether to copy the input data into the output file.

Running

write([filename, compression, copy, …]) Write the model input parameters to an HDF5 file
run([filename, logfile, mpi, n_processes, …]) Run the model (should be called after write()).

Using grids from files

use_grid_from_file(filename[, path, dust]) Use a grid from disk and don’t read it in.

Re-using previous models

read(filename[, only_initial]) Read in a previous model file
use_geometry(filename) Use the grid from an existing output or input file
use_quantities(filename[, quantities, …]) Use physical quantities from an existing output file
use_sources(filename) Use sources from an existing output file
use_image_config(filename) Use image configuration from an existing output or input file
use_run_config(filename) Use runtime configuration from an existing output or input file
use_output_config(filename) Use output configuration from an existing output or input file

Methods (detail)

add_source(source)
add_point_source(*args, **kwargs)
add_spherical_source(*args, **kwargs)
add_external_spherical_source(*args, **kwargs)
add_external_box_source(*args, **kwargs)
add_map_source(*args, **kwargs)
add_plane_parallel_source(*args, **kwargs)
set_grid(grid)
set_cartesian_grid(x_wall, y_wall, z_wall)
set_cylindrical_polar_grid(w_wall, z_wall, p_wall)
set_spherical_polar_grid(r_wall, t_wall, p_wall)
set_octree_grid(x, y, z, dx, dy, dz, refined)
set_amr_grid(description)
add_density_grid(density, dust, specific_energy=None, merge_if_possible=False)

Add a density grid to the model

Parameters:
density : np.ndarray or grid quantity

The density of the dust. This can be specified either as a 3-d Numpy array for cartesian, cylindrical polar, or spherical polar grids, as a 1-d array for octree grids, or as a grid quantity object for all grid types. Grid quantity objects are obtained by taking an instance of a grid class (e.g. AMRGrid, CartesianGrid, …) and specifying the quantity as an index, e.g. amr['density'] where amr is an AMRGrid object.

dust : str or dust instance

The dust properties, specified either as a string giving the filename of the dust properties, or a as an instance of a dust class (e.g. SphericalDust, IsotropicDust, …).

specific_energy : np.ndarray or grid quantity, optional

The specific energy of the density grid. Note that in order for this to be useful, the number of initial iterations should be set to zero, otherwise these values will be overwritten after the first initial iteration.

merge_if_possible : bool

Whether to merge density arrays that have the same dust type

add_peeled_images(sed=True, image=True)

Define a set of (peeled) images/SEDs.

This returns a PeeledImageConf instance that can be used to set the image/SED parameters.

Parameters:
sed : bool, optional

Whether to compute a set of SEDs

image : bool, optional

Whether to compute a set of images

Returns:
image_conf : PeeledImageConf instance

Instance that can be used to set the image/SED parameters

Examples

>>> image = m.add_peeled_images(sed=False, image=True)
>>> image.set_wavelength_range(250, 0.01, 1000.)
>>> ...
add_binned_images(sed=True, image=True)

Define a set of (binned) images/SEDs.

This returns a BinnedImageConf instance that can be used to set the image/SED parameters.

Parameters:
sed : bool, optional

Whether to compute a set of SEDs

image : bool, optional

Whether to compute a set of images

Returns:
image_conf : BinnedImageConf instance

Instance that can be used to set the image/SED parameters

Examples

>>> image = m.add_binned_images(sed=False, image=True)
>>> image.set_wavelength_range(250, 0.01, 1000.)
>>> ...
set_seed(seed)

Set the seed for the random number generation

Parameters:
seed : int

The seed with which to initialize the random number generation. This should be negative.

set_propagation_check_frequency(frequency)

Set how often to check that the photon is in the right cell

During photon propagation, it is possible that floating point issues cause a photon to end up in the wrong cell. By default, the code will randomly double check the position and cell of a photon for every 1 in 1000 cell wall crossings, but this can be adjusted with this method. Note that values higher than around 0.001 will cause the code to slow down.

Parameters:
frequency : float

How often the photon position and cell should be double-checked (1 is always, 0 is never).

set_monochromatic(monochromatic, wavelengths=None, energy_threshold=1e-10)

Set whether to do the radiation transfer at specific wavelengths.

Parameters:
monochromatic : bool

Whether to carry out radiation transfer at specific frequencies or wavelengths

wavelengths : iterable of floats, optional

The wavelengths to compute the radiation transfer for, in microns

energy_threshold : float

In this mode, the propagation of photons is done by scattering photons and decreasing their energy by the albedo each time. Once their energy has decreased by a ratio corresponding to this parameter, the photons are terminated.

If `monochromatic` is True then `wavelengths` is required
set_minimum_temperature(temperature)

Set the minimum temperature for the dust

Parameters:
temperature : float, list, tuple, or Numpy array

Notes

This method should not be used in conjunction with set_minimum_specific_energy - only one of the two should be used.

set_minimum_specific_energy(specific_energy)

Set the minimum specific energy for the dust

Parameters:
specific_energy : float, list, tuple, or Numpy array

Notes

This method should not be used in conjunction with set_minimum_temperature - only one of the two should be used.

set_specific_energy_type(specific_energy_type)

Set whether to use the specific specific energy as an initial value or an additional component at each iteration.

This only has an effect if a specific energy was specified during add_density_grid.

Parameters:
specific_energy_type : str

Can be 'initial' (use only as initial value) or 'additional' (add at every iteration)

set_n_initial_iterations(n_iter)

Set the number of initial iterations for computing the specific energy in each cell.

Parameters:
n_iter : int

The number of initial iterations

set_raytracing(raytracing)

Set whether to use raytracing for the non-scattered flux

If enabled, only scattered photons are peeled off in the iteration following the initial iterations, and an additional final iteration is carrried out, with raytracing of the remaining flux (sources and thermal and non-thermal dust emission).

Parameters:
raytracing : bool

Whether or not to use raytracing in the final iteration

set_max_interactions(inter_max, warn=True)

Set the maximum number of interactions a photon can have.

Parameters:
inter_max : int

Maximum number of interactions for a single photon. This can be used to prevent photons from getting stuck in very optically thick regions, especially if the modified random walk is not used.

warn : bool, optional

Whether to emit a warning whenever photons are killed for exceeding the maximum number of iterations.

set_max_reabsorptions(reabs_max, warn=True)

Set the maximum number of successive reabsorptions by a source that a photon can have.

Parameters:
reabs_max : int

Maximum number of reabsorptions for a single photon.

warn : bool, optional

Whether to emit a warning whenever photons are killed for exceeding the maximum number of reabsorptions.

set_pda(pda)

Set whether to use the Partial Diffusion Approximation (PDA)

If enabled, the PDA is used to compute the specific energy in cells which have seen few or no photons by formally solving the diffusion equations, using the cells with valid specific energies as boundary conditions.

Parameters:
pda : bool

Whether or not to use the PDA

References

Min et al. 2009, Astronomy and Astrophysics, 497, 155

set_mrw(mrw, gamma=1.0, inter_max=1000, warn=True)

Set whether to use the Modified Random Walk (MRW) approximation

If enabled, the MRW speeds up the propagation of photons in very optically thick regions by locally setting up a spherical diffusion region.

Parameters:
mrw : bool

Whether or not to use the MRW

gamma : float, optional

The parameter describing the starting criterion for the MRW. The MRW is carried out if the distance to the closest cell is larger than gamma times the Rosseland mean free path.

inter_max : int, optional

Maximum number of interactions during a single random walk. This can be used to prevent photons from getting stuck in the corners of cells in very optically thick regions, where the MRW stars to become inefficient itself.

warn : bool, optional

Whether to emit a warning whenever photons are killed for exceeding the maximum number of mrw steps.

References

Min et al. 2009, Astronomy and Astrophysics, 497, 155

set_convergence(convergence, percentile=100.0, absolute=0.0, relative=0.0)

Set whether to check for convergence over the initial iterations

If enabled, the code will check whether the specific energy absorbed in each cell has converged. First, the ratio between the previous and current specific energy absorbed in each cell is computed in each cell, and the value at the specified percentile (percentile) is found. Then, convergence has been achieved if this value is less than an absolute threshold (absolute), and if it changed by less than a relative threshold ratio (relative).

Parameters:
convergence : bool

Whether or not to check for convergence.

percentile : float, optional

The percentile at which to check for convergence.

absolute : float, optional

The abolute threshold below which the percentile value of the ratio has to be for convergence.

relative : float, optional

The relative threshold below which the ratio in the percentile value has to be for convergence.

set_kill_on_absorb(kill_on_absorb)

Set whether to kill absorbed photons

Parameters:
kill_on_absorb : bool

Whether to kill absorbed photons

set_forced_first_interaction(forced_first_interaction, algorithm='wr99', baes16_xi=0.5)

Set whether to ensure that photons scatter at least once before escaping the grid.

Parameters:
forced_first_interaction : bool

Whether to force at least one scattering before escaping the grid

algorithm : ‘wr99’ or ‘baes16’

Which algorithm to use for the forced first interaction. The algorithms are described in the notes below.

Notes

The ‘wr99’ algorithm refers to that described in Wood & Reynolds, 1999, The Astrophysical Journal, 525, 799. During normal un-forced photon propagation, we sample the optical depth from a probability density function (PDF) that follows exp(-tau) from tau=0 to infinity. The Wood and Reynolds algorithm modifies the PDF to be truncated at tau_escape (the optical depth for the photon to escape the grid). This ensures that all photons interact at least once before leaving the system. This algorithm is ideal for cases where the optical depths are very small and you are interested in making images. Note that this algorithm does not apply to the temperature calculation iterations since it is not needed there.

The ‘baes16’ algorithm refers to that described in Baes et al. 2019, Astronomy and Astrophysics, 590, A55. In this algorithm, the PDF is the weighted combination of a truncated decaying exponential and a constant, which ensures that interactions will occur with a reasonable probability anywhere along the photon escape path. This is useful for cases where there are shadowed regions that otherwise would not receive many photons. The relative weight of the truncated exponential versus the constant is given by baes16_xi, which should be in the range 0 to 1.

set_output_bytes(io_bytes)

Set whether to output physical quantity arrays in 32-bit or 64-bit

Parameters:
io_bytes : int

The number of bytes for the output. This should be either 4 (for 32-bit) or 8 (for 64-bit).

set_sample_sources_evenly(sample_sources_evenly)

If set to ‘True’, sample evenly from all sources and apply probability weight based on relative luminosities. Otherwise, sample equal energy photons from sources with probability given by relative luminosities.

Parameters:
sample_evenly : bool

Whether to sample different sources evenly

set_enforce_energy_range(enforce)

Set how to deal with cells that have specific energy rates that are below or above that provided in the mean opacities and emissivities.

Parameters:
enforce : bool

Whether or not to reset specific energies that are above or below the range of values used to specify the mean opacities and emissivities to the maximum or minimum value of the range. Setting this to True modifies the energy in the simulation, but ensures that the emissivities are consistent with the energy in the cells. Setting this to False means that the total energy in the grid will be correct, but that the emissivities may be inconsistent with the energy in the cells (if an energy is out of range, the code will pick the closest available one). In both cases, warnings will be displayed to notify the user whether this is happening.

set_copy_input(copy)

Set whether to copy the input data into the output file.

Parameters:
copy : bool

Whether to copy the input data into the output file (True) or whether to link to it (False)

write(filename=None, compression=True, copy=True, absolute_paths=False, wall_dtype=<class 'float'>, physics_dtype=<class 'float'>, overwrite=True)

Write the model input parameters to an HDF5 file

Parameters:
filename : str

The name of the input file to write. If no name is specified, the filename is constructed from the model name.

compression : bool

Whether to compress the datasets inside the HDF5 file.

copy : bool

Whether to copy all external content into the input file, or whether to just link to external content.

absolute_paths : bool

If copy=False, then if absolute_paths is True, absolute filenames are used in the link, otherwise the path relative to the input file is used.

wall_dtype : type

Numerical type to use for wall positions.

physics_dtype : type

Numerical type to use for physical grids.

overwrite : bool

Whether to overwrite any pre-existing file

run(filename=None, logfile=None, mpi=False, n_processes=4, overwrite=False)

Run the model (should be called after write()).

Parameters:
filename : str, optional

The output filename for the model. If not specified, then if the input file name contains .rtin, then this is replaced with .rtout, and otherwise .rtout is appended to the input filename.

logfile : str, optional

If specified, the standard output and errors will be output to this log file

mpi : bool, optional

Whether to run the model using the parallel (MPI) version of Hyperion.

n_processes : int, optional

If mpi is set to True, this can be used to specify the number of processes to run Hyperion on.

overwrite : bool, optional

If set to True, the output file is overwritten without warning.

use_grid_from_file(filename, path='/', dust=[])

Use a grid from disk and don’t read it in.

Parameters:
filename : str

The name of the file to read from

path : str, optional

The path where the grid is located inside the file

dust : list, optional

A list containing a dust file or object for each dust index in the grid

classmethod read(filename, only_initial=True)

Read in a previous model file

This can be used to read in a previous input file, or the input in an output file (which is possible because the input to a model is stored or linked in an output file).

If you are interested in re-using the final specific energy (and final density, if present) of a previously run model, you can use:

>>> m = Model.read('previous_model.rtout', only_initial=False)
Parameters:
filename : str

The name of the file to read the input data from

only_initial : bool, optional

Whether to use only the initial quantities, or whether the final specific energy (and optionally density) from the previous model can be used. By default, only the input density (and specific energy, if present) are read in.

use_geometry(filename)

Use the grid from an existing output or input file

Parameters:
filename : str

The file to read the grid from. This can be either the input or output file from a radiation transfer run.

use_quantities(filename, quantities=['density', 'specific_energy'], use_minimum_specific_energy=True, use_dust=True, copy=True, only_initial=False)

Use physical quantities from an existing output file

Parameters:
filename : str

The file to read the quantities from. This should be the output file of a radiation transfer run.

quantities : list

Which physical quantities to read in. Can include ‘density’ and ‘specific_energy’.

copy : bool

Whether to copy the quantities into the new input file, or whether to just link to them.

use_minimum_specific_energy : bool

Whether to also use the minimum specific energy values from the file specified

use_dust : bool

Whether to also use the dust properties from the file specified

copy : bool

Whether to read in a copy of the data. If set to False, then the physical quantities will only be links to the specified HDF5 file, and can therefore not be modified.

only_initial : bool, optional

Whether to use only the initial quantities, or whether the final specific energy (and optionally density) from the previous model can be used. By default, only the input density (and specific energy, if present) are read in.

use_sources(filename)

Use sources from an existing output file

Parameters:
filename : str

The file to read the sources from. This should be the input or output file of a radiation transfer run.

use_image_config(filename)

Use image configuration from an existing output or input file

Parameters:
filename : str

The file to read the parameters from. This can be either the input or output file from a radiation transfer run.

use_run_config(filename)

Use runtime configuration from an existing output or input file

Parameters:
filename : str

The file to read the parameters from. This can be either the input or output file from a radiation transfer run.

use_output_config(filename)

Use output configuration from an existing output or input file

Parameters:
filename : str

The file to read the parameters from. This can be either the input or output file from a radiation transfer run.