API Documentation

class meshoid.Meshoid.Meshoid(pos, m=None, kernel_radius=None, des_ngb=None, boxsize=None, verbose=False, particle_mask=None, n_jobs=-1)

Meshoid object that stores particle data and kernel-weighting quantities, and implements various methods for integral and derivative operations.

Construct a meshoid instance.

Parameters:
  • pos (array_like, shape (n,k)) – The coordinates of the particle data

  • m (array_like, shape (n,), optional) – Masses of the particles for computing e.g. mass densities,

  • kernel_radius (array_like, shape (n,), optional) – Kernel support radii (e.g. SmoothingLength”) of the particles. Can be computed adaptively where needed, if not provided.

  • des_ngb (positive int, optional) – Number of neighbors to search for kernel density and weighting calculations (defaults to 4/20/32 for 1D/2D/3D)

  • boxsize (positive float, optional) – Side-length of box if periodic topology is to be assumed

  • verbose (boolean, optinal, default: False) – Whether to print what meshoid is doing to stdout

  • particle_mask (array_like, optional) – array-like of integer indices OR boolean mask of the particles you want to compute things for

  • n_jobs (int, optional) – number of cores to use for parallelization (default: -1, uses all cores available)

Returns:

Meshoid instance created from particle data

Return type:

Meshoid

ComputeDWeights(order=1, weighted=True)

Computes the weights required to compute least-squares gradient estimators on data colocated on the meshoid

Parameters:
  • order (int, optional) – 1 to compute weights for first derivatives, 2 for the Jacobian matrix (default: 1)

  • weighted (boolean, optional) – whether to kernel-weight the least-squares gradient solutions (default: True)

Curl(v)

Computes the curl of a vector field.

Parameters:

v (array_like) – shape (N,3) array containing a vector field colocated on the meshoid

Return type:

shape (N,3) array containing the curl of vector field v

D(f, order=1, weighted=True)

Computes the kernel-weighted least-squares gradient estimator of the function f.

Parameters:
  • f (array_like) – shape (N,…) array of (possibly vector- or tensor-valued) function values (N is the total number of particles)

  • order (int, optional) – desired order of precision, either 1 or 2

Returns:

  • (Nmask, …, dim) array of partial derivatives, evaluated at the

  • positions of the particles in the particle mask

D2(f, weighted=True)

Computes the kernel-weighted least-squares Jacobian estimator of the function f.

Parameters:

f (array_like) – shape (N,…) array of (possibly vector- or tensor-valued) function values (N is the total number of particles)

Returns:

  • (Nmask, …, N_derivs) array of partial second derivatives, evaluated at the positions of the particles in the particle mask

  • Here N_derivs is the number of unique second derivatives in the given number of dimensions (1 for 1D, 3 for 2D,)

  • 6 for 3D etc.

  • For 2D, the order is [xx,yy,xy]

  • for 3D, the order is [xx,yy,zz,xy,yz,zx]

Density()

Returns the mass density (or number density, if mass not provided) of each particle in a (N,) array

Return type:

self.density - (N,) array of particle densities

DepositToGrid(f, weights=None, size=None, center=None, res=128)

Deposits a conserved quantity (e.g. mass, momentum, energy) to a 3D grid and returns the density of that quantity on that grid

Parameters:
  • f (array_like) – Shape (N,) array of conserved quantity values colocated at the particle coordinates

  • weights (array_like, optional) – Shape (N,) array of weights for kernel-weighted interpolation

  • size (float, optional) – Side length of the grid - defaults to the self.L value: either 2 times the 90th percentile radius from the center, or the specified boxsize

  • center (array_like, optional) – Center of the grid - defaults to self.center value, either the box center if boxsize given, or average particle position

  • res (integer, optional) – Resolution of the grid - default 128

Return type:

Shape (res,res,res) grid of the density of the conserved quantity

Div(v)

Computes the divergence of a vector field.

Parameters:

v (array_like) – shape (N,3) array containing a vector field colocated on the meshoid

Return type:

shape (N,) array containing the divergence of v

Integrate(f)

Computes the volume integral of a quantity over the volume partition of the domain

Parameters:
  • f (array_like) – Shape (Nmask, …) function colocated on the meshoid

  • Returns

  • domain (integral of f over the) –

InterpToGrid(f, size=None, center=None, res: int = 128, order: int = 1, return_grid: bool = False)

Interpolates the quantity f defined on the meshoid to the cell centers of a 3D Cartesian grid

Parameters:
  • f (array_like) – Shape (N,) function colocated at the particle coordinates

  • size (float, optional) – Side length of the grid - defaults to the self.L value: either 2 times the 90th percentile radius from the center, or the specified boxsize

  • center (array_like, optional) – Center of the grid - defaults to self.center value, either the box center if boxsize given, or average particle position

  • res (integer, optional) – Resolution of the grid - default 128

  • order (int, optional) – Order of the reconstruction on the slice: 0, 1, or 2 (default: 1)

Return type:

Shape (res,res,res) grid of cell-centered interpolated values

KDE(grid, bandwidth=None)

Computes the kernel density estimate of the meshoid points on a 1D grid

Parameters:
  • grid (array_like) – 1D array of coordintes upon which to compute the KDE

  • bandwidth (float, optional) – constant bandwidth of the kernel, defined as the radius of support of the cubic spline kernel.

Return type:

array containing the density of particles defined at the grid points

KernelAverage(f)

Computes the kernel-weighted average of a function

Parameters:

f (array_like, shape (N ,...)) – Shape (N, …) function colocated on the meshoid

Return type:

Shape (N, …) array of kernel-averaged values of f

KernelVariance(f)

Computes the standard deviation of a quantity over all nearest neighbours

Parameters:

f (array_like, shape (N ,...)) – Shape (N, …) function colocated on the meshoid

Return type:

Shape (Nmask,…) array of standard deviations of f over the kernel

NearestNeighbors()

Returns the indices of the N_ngb nearest neighbors of each particle in a (N, N_ngb) array

Return type:

self.ngb - (N,ngb) array of the indices of the each particle’s Nngb nearest neighbors

NeighborDistance()

Returns the distances of the N_ngb nearest neighbors of each particle in a (N, N_ngb) array

Return type:

self.ngbdist - (N,Nngb) array of distances to nearest neighbors of each particle.

ProjectedAverage(f, size=None, center=None, res=128, smooth_fac=1.0)

Computes the average value of a quantity f along a Cartesian grid of sightlines from +/- infinity.

Parameters:
  • f – (N,) array containing the quantity you want the average of

  • size – the side length of the window of sightlines (default: None, will use the meshoid’s predefined side length’)

  • plane – the direction you want to project along, of x, y, or z (default: ‘z’)

  • center – (2,) or (3,) array containing the coordaintes of the center of the grid (default: None, will use the meshoid’s predefined center)

  • res – the resolution of the grid of sightlines (default: 128)

  • smooth_fac – smoothing lengths are increased by this factor (default: 1.)

Return type:

(res,res) array containing the averages along sightlines

Projection(f, size=None, center=None, res=128, smooth_fac=1.0)

Computes the integral of quantity f along a Cartesian grid of sightlines from +/- infinity. e.g. plugging in 3D density for f will return surface density.

Parameters:
  • f – (N,) array containing the quantity you want the projected integral of

  • size – the side length of the window of sightlines (default: None, will use the meshoid’s predefined side length’)

  • plane – the direction you want to project along, of x, y, or z (default: ‘z’)

  • center – (2,) or (3,) array containing the coordaintes of the center of the grid (default: None, will use the meshoid’s predefined center)

  • res – the resolution of the grid of sightlines (default: 128)

  • smooth_fac – smoothing lengths are increased by this factor (default: 1.)

Return type:

(res, res) array containing the projected values

Reconstruct(f: ndarray, target_points: ndarray, order: int = 1)

Gives the value of a function f colocated on the meshoid points reconstructed at an arbitrary set of points

Parameters:
  • f (array_like) – The quantity defined on the set of meshoid points that we wish to reconstruct at the target points, first dimension should be N_mask

  • target_points (array_like) – The shape (N_target,dim) array of points where you would like to reconstruct the function

  • order (int, optional) – The order of the reconstruction (default 1): 0 - nearest-neighbor value 1 - linear reconstruction from the nearest neighbor 2 - quadratic reconstruction from the nearest neighbor

Returns:

f_target – Values of f reconstructed at the target points

Return type:

ndarray

Slice(f, size=None, plane='z', center=None, res: int = 128, order: int = 1, return_grid: bool = False)

Gives the value of a function f deposited on a 2D Cartesian grid slicing through the 3D domain.

Parameters:
  • f (array_like) – Quantity to be reconstructed on the slice grid

  • size (optional) – Side length of the grid (default: None, will use the meshoid’s predefined side length’)

  • plane (optional) – The direction of the normal of the slicing plane, one of x, y, or z, OR can give an iterable containing 2 orthogonal basis vectors that specify the plane (default: ‘z’)

  • center (array_like, optional) – (2,) or (3,) array containing the coordaintes of the center of the grid (default: None, will use the meshoid’s predefined center)

  • res (int, optional) – the resolution of the grid (default: 128)

  • order (int, optional) – Order of the reconstruction on the slice: 0, 1, or 2 (default: 1)

  • return_grid (bool, optional) – Also return the grid coordinates corresponding to the slice

Returns:

slice – Shape (res,res,…) grid of values reconstructed on the slice

Return type:

array_like

SmoothingLength()

Returns the neighbor kernel radii of of each particle in a (N,) array

Return type:

(N,) array of particle smoothing lengths

SurfaceDensity(f=None, size=None, center=None, res=128, smooth_fac=1.0, conservative=False)

Computes the surface density of a quantity f defined on the meshoid on a grid of sightlines. e.g. if f is the particle masses, you will get mass surface density.

Parameters:
  • f – the quantity you want the surface density of (default: particle mass)

  • size – the side length of the window of sightlines (default: None, will use the meshoid’s predefined side length’)

  • plane – the direction you want to project along, of x, y, or z (default: ‘z’)

  • center – (2,) or (3,) array containing the coordaintes of the center of the grid (default: None, will use the meshoid’s predefined center)

  • res – the resolution of the grid of sightlines (default: 128)

  • smooth_fac – smoothing lengths are increased by this factor (default: 1.)

Return type:

(res,res) array containing the column densities integrated along sightlines

TreeUpdate()

Computes or updates the neighbor lists, smoothing lengths, and densities of particles.

Volume()

Returns the effective particle volume = (mass / density)^(1/3)

Return type:

self.vol - (N,) array of particle volumes

meshoid.kernel_density.HsmlIter(neighbor_dists, dim=3, error_norm=1e-06)

Performs the iteration to get smoothing lengths, according to Eq. 26 in Hopkins 2015 MNRAS 450.

Arguments: neighbor_dists: (N, Nngb) array of distances from particles to their Nngb nearest neighbors

Keyword arguments: dim - Dimensionality (default: 3) error_norm - Tolerance in the particle number density to stop iterating at (default: 1e-6)

meshoid.kernel_density.kernel2d(q)

Returns the normalized 2D spline kernel evaluated at radius q

meshoid.kernel_density.kernel3d(q)

Returns the normalized 3D spline kernel evaluated at radius q

meshoid.grid_deposition.GridAverage(f, x, h, center, size, res=100, box_size=-1)

Computes the number density-weighted average of a function f, integrated along sightlines on a Cartesian grid. ie. integral(n f dz)/integral(n dz) where n is the number density and z is the direction of the sightline.

Arguments: f - (N,) array of the conserved quantity that you want the surface density of (e.g. particle masses) x - (N,3) array of particle positions h - (N,) array of particle smoothing lengths center - (2,) array containing the coordinates of the center of the map size - side-length of the map res - resolution of the grid

meshoid.grid_deposition.GridDensity(f, x, h, center, size, res=100, box_size=-1.0)

Estimates the density of the conserved quantity f possessed by each particle (e.g. mass, momentum, energy) on a 3D grid

Arguments: f - (N,) array of the conserved quantity that you want the surface density of (e.g. particle masses) x - (N,3) array of particle positions h - (N,) array of particle smoothing lengths center - (2,) array containing the coorindates of the center of the map size - side-length of the map res - resolution of the grid

meshoid.grid_deposition.GridSurfaceDensity(f, x, h, center, size, res=100, box_size=-1, parallel=True, conservative=False)

Computes the surface density of conserved quantity f colocated at positions x with smoothing lengths h. E.g. plugging in particle masses would return mass surface density. The result is on a Cartesian grid of sightlines, the result being the density of quantity f integrated along those sightlines.

Parameters:
  • (N (h -) –

  • masses) () array of the conserved quantity that you want the surface density of (e.g. particle) –

  • (N

  • positions (3) array of particle) –

  • (N

  • lengths () array of particle smoothing) –

  • (2 (center -) –

  • map (size - side-length of the) –

  • map

  • grid (res - resolution of the) –

  • parallel (parallel - whether to run in) –

  • cores (if numeric then how many) –

meshoid.grid_deposition.GridSurfaceDensity_conservative_core(f, x, h, center, size, res=100, box_size=-1)

Computes the surface density of conserved quantity f colocated at positions x with smoothing lengths h. E.g. plugging in particle masses would return mass surface density.

This method performs a kernel-weighted deposition so that the quantity deposited to the grid is conserved to machine precision.

Arguments: f - (N,) array of the conserved quantity that you want the surface density of (e.g. particle masses) x - (N,3) array of particle positions h - (N,) array of particle smoothing lengths center - (2,) array containing the coorindates of the center of the map size - side-length of the map res - resolution of the grid

meshoid.grid_deposition.GridSurfaceDensity_core(f, x, h, center, size, res=100, box_size=-1)

Computes the surface density of conserved quantity f colocated at positions x with smoothing lengths h. E.g. plugging in particle masses would return mass surface density. The result is on a Cartesian grid of sightlines, the result being the density of quantity f integrated along those sightlines.

Arguments: f - (N,) array of the conserved quantity that you want the surface density of (e.g. particle masses) x - (N,3) array of particle positions h - (N,) array of particle smoothing lengths center - (2,) array containing the coorindates of the center of the map size - side-length of the map res - resolution of the grid

meshoid.grid_deposition.Grid_PPZ_DataCube(f, x, h, center, size, z, h_z, res, box_size=-1)

A modified version of the GridSurfaceDensity script, it computes the PPZ datacube of conserved quantity f, where Z is an arbitrary data dimension (e.g. using line-of-sight velocity as Z gives the usual astro PPV datacubes, using position gives a PPP cube, which is just the density).

Arguments: f - (N,) array of the conserved quantity that you want the surface density of (e.g. particle masses) x - (N,3) array of particle positions h - (N,) array of particle smoothing lengths z - (N,) array of particle positions in the Z dimension (e.g. line of sight velocity values) h_z - (N,) array of uncertainties (“smoothing lengths”) in the Z dimension (e.g. thermal velocity) center - (3,) array containing the coordinates of the center of the PPZ map size - (2) side-length of the map, first value for the PP map, second value is for Z (e.g. max velocity - min velocity) res - (2) resolution of the PPX grid, first value is for the PP map, second is for Z (e.g. [128, 16] means a 128x128x16 PPX cube)

meshoid.grid_deposition.Grid_PPZ_DataCube_Multigrid(f, x, h, center, size, z, h_z, res, box_size=-1, N_grid_kernel=8)

Faster, multigrid version of Grid_PPZ_DataCube. Since the third dimension is separate from the spatial ones, we only do the multigrid approach on the spatial grid. See Grid_PPZ_DataCube for desription of inputs

meshoid.grid_deposition.WeightedGridInterp3D(f, wt, x, h, center, size, res=100, box_size=-1)

Peforms a weighted grid interpolation of quantity f onto a 3D grid

Arguments: f - (N,) array of the function defined on the point set that you want to interpolate to the grid wt - (N,) array of weights x - (N,3) array of particle positions h - (N,) array of particle smoothing lengths center - (3,) array containing the coorindates of the center of the map size - side-length of the map res - resolution of the grid

meshoid.grid_deposition.grid_dx_from_coordinate(x: float, index: int, grid_length: float, grid_center: float, N_grid: int, box_size=None)

Returns the nearest image coordinate difference from a given coordinate x to the cell-centered grid-point residing at index

Parameters:
  • x (float) – The original coordinate, from which to compute coordinate difference

  • index (int) – Grid index, running from 0 to N_grid-1

  • grid_length (float) – Side-length of the grid

  • N_grid (int) – Number of cells per grid dimension

  • grid_center (float) – Coordinate of the grid center

  • box_size – Size of the periodic domain - if not None, we assume the domain coordinates run from [0,box_size)

Returns:

dx – The nearest-image Cartesian coordinate difference from x to grid cell i

Return type:

float

meshoid.grid_deposition.grid_index_to_coordinate(index: int, grid_length: float, grid_center: float, N_grid: int, box_size=None)

Convert the index of a grid to the cell-centered coordinate of that grid cell, in a given dimension

Parameters:
  • index (int) – Grid index, running from 0 to N_grid-1

  • grid_length (float) – Side-length of the grid

  • grid_center (float) – Coordinate of the grid center

  • N_grid (int) – Number of cells per grid dimension

  • box_size – Size of the periodic domain - if not None, we assume the domain coordinates run from [0,box_size)

Returns:

x – Coordinate of the center of the grid cell

Return type:

float

Functions for computing the matrices weights for least-squares derivative operators

meshoid.derivatives.get_num_derivs(dim: int, order: int) int

Returns number of unique derivatives to compute for a given matrix order and dimension

meshoid.derivatives.gradient_weights(pos, ngb, kernel_radius, indices, boxsize=None, weighted=True, order=1)

Computes the N_particles (dim x N_ngb) matrices that encode the least- squares gradient operators

meshoid.derivatives.kernel_dx_and_weights(i: int, pos: ndarray, ngb: ndarray, kernel_radius: float, boxsize=None, weighted: bool = True)

Computes coordinate differences and weights for the neighbors in the kernel of particle i.

meshoid.derivatives.nearest_image(dx_coord, boxsize)

Returns separation vector for nearest image, given the coordinate difference dx_coord and assuming coordinates run from 0 to boxsize

meshoid.derivatives.polyfit_leastsq_matrices(dx: ndarray, weights: ndarray, order: int = 1)

Return the left-hand side and right-hand side matrices for the system of equations for a weighted least-squares fit to a polynomial needed to compute the least-squares gradient matrices:

(A^T W A) X = (A^T W) Y

where A is the polynomial Vandermonde matrix, W are the weights, X Y are the differences in function values, and X are the unknown derivative components we wish to solve for.

Routines for radiative transfer calculations

meshoid.radiation.radtransfer.dust_abs_opacity(wavelength_um: ndarray = array([150, 250, 350, 500]), XH=0.71, Zd=1.0) ndarray

Returns the dust+PAH absorption opacity in cm^2/g at a set of wavelengths in micron

Parameters:
  • wavelength_um (array_like) – Shape (num_bands) array of wavelengths in micron

  • XH (float, optional) – Mass fraction of hydrogen (needed to convert from per H to per g)

  • Zd (float, optional) – Dust-to-gas ratio normalized to Solar neighborhood

Returns:

kappa – shape (num_bands,) array of opacities in cm^2/g

Return type:

array_like

meshoid.radiation.radtransfer.dust_emission_map(x_pc, m_msun, h_pc, Tdust, size_pc, res, wavelengths_um=array([150, 250, 350, 500]), center_pc=0) ndarray

Generates a map of dust emission in cgs units for specified wavelengths, neglecting scattering (OK for FIR/submm wavelengths)

Parameters:
  • x_pc (array_like) – Shape (N,3) array of coordinates in pc

  • m_msun (array_like) – Shape (N,) array of masses in msun

  • h_pc (array_like) – Shape (N,) array of kernel radii in pc

  • Tdust (array_like) – Shape (N,) array of dust temperatures in K

  • wavelengths_um (array_like) – Shape (num_bands,) array of wavelengths in micron

  • size_pc (float) – Size of the image in pc

  • res (int) – Image resolution

  • center_pc (array_like, optional) – Shape (3,) array providing the coordinate center of the image

Returns:

intensity – shape (res,res,num_bands) datacube of dust emission intensity in erg/s/cm^2/sr

Return type:

array_like

meshoid.radiation.radtransfer.radtransfer(j, m, kappa, x, h, gridres, L, center=0, i0=0)

Simple radiative transfer solver

Solves the radiative transfer equation with emission and absorption along a grid of sightlines, at multiple frequencies. Raytraces in parallel, with number of threads determined by numba set_num_threads

Parameters:
  • j (array_like) – shape (N, num_freqs) array of particle emissivities per unit mass (e.g. erg/s/g/Hz)

  • m (array_like) – shape (N,) array of particle masses

  • kappa (array_like) – shape (N,) array of particle opacities (dimensions: length^2 / mass)

  • x (array_like) – shape (N,3) array of particle positions

  • h (array_like) – shape (N,) array containing kernel radii of the particles

  • gridres (int) – image resolution

  • center (array_like) – shape (3,) array containing the center coordinates of the image

  • L (float) – size of the image window in length units

  • i0 (array_like, optional) – shape (num_freqs,) or (gridres,greidres,num_freqs) array of background intensities

Returns:

image – shape (res,res) array of integrated intensities, in your units of power / length^2 / sr / frequency (this is the quantity I_ν in RT)

Return type:

array_like

meshoid.radiation.radtransfer.radtransfer_singlethread(j, m, kappa, x, h, gridres, L, center=0, i0=0)

Simple radiative transfer solver

Solves the radiative transfer equation with emission and absorption along a grid of sightlines, at multiple frequencies

Parameters:
  • j (array_like) – shape (N, num_freqs) array of particle emissivities per unit mass (e.g. erg/s/g/Hz)

  • m (array_like) – shape (N,) array of particle masses

  • kappa (array_like) – shape (N,) array of particle opacities (dimensions: length^2 / mass)

  • x (array_like) – shape (N,3) array of particle positions

  • h (array_like) – shape (N,) array containing kernel radii of the particles

  • gridres (int) – image resolution

  • center (array_like) – shape (3,) array containing the center coordinates of the image

  • L (float) – size of the image window in length units

  • i0 (array_like, optional) – shape (num_freqs,) or (gridres,greidres,num_freqs) array of background intensities

  • return_taumap (boolean, optional) – Whether to return the optical depth map alongside the intensity map (default False)

Returns:

intensity, tau – shape (res,res,num_bands) array of integrated intensity, in your units of power / length^2 / sr / frequency (this is the quantity I_ν in RT), and of optical depth in the different bands along the lines of sight

Return type:

array_like, array_like

meshoid.radiation.radtransfer.thermal_emissivity(kappa, T, wavelengths_um=array([150, 250, 350, 500]))

Returns the thermal emissivity j_ν = 4 pi kappa_ν B_ν(T) in erg/s/g for a specified list of wavelengths, temperatures, and opacities defined at those wavelengths

Parameters:
  • kappa (array_like) – shape (N,num_bands) array of opacities

  • T (array_like) – shape (N,) array of temperatures

  • wavelengths (array_like) – shape (num_bands,) array of wavelengths

Returns:

j – shape (N,num_bands) array of thermal emissivities

Return type:

array_like

meshoid.radiation.modified_blackbody.blackbody_residual_jacobian(freqs, sed_error, logtau, beta, logT)

Jacobian for least-squares fit to modified blackbody

Parameters:
  • freqs (array_like) – Photon frequencies in Hz

  • sed_error (array_like) – SED errors

  • logtau (float) – log10 of optical depth at 500um

  • beta (float) – Dust spectral index beta

  • logT (float) – log10 of temperature in K

Returns:

jac – Shape (N,3) matrix of the gradient of the modified blackbody function with respect to logtau, beta, and logT

Return type:

np.ndarray

meshoid.radiation.modified_blackbody.modified_blackbody_fit_gaussnewton(sed, sed_error, wavelengths, p0=(1.0, 1.5, 30.0))

Fits a single SED to a modified blackbody using Gauss-Newton method

Parameters:
  • sed (array_like) – Shape (num_bands,) array of dust emission intensities

  • sed – Shape (num_bands,) array of dust emission intensity errors

  • wavelengths – Shape (num_bands,) array of wavelengths at which the SEDs are computed

  • p0 (tuple, optional) – Shape (3,) tuple of initial guesses for tau(500um), beta, T

Returns:

params

Shape (3,) array of best-fitting modified blackbody parameters:
  • tau(500um): optical depth at 500 micron

  • beta: spectral index

  • T: temperature

Return type:

np.ndarray

meshoid.radiation.modified_blackbody.modified_blackbody_fit_image(image, image_error, wavelengths)

Fit each pixel in a datacube to a modified blackbody

Parameters:
  • image (array_like) – Shape (N,N,num_bands) datacube of dust emission intensities

  • image_error (array_like) – Shape (N,N,num_bands) datacube of dust emission intensity errors

  • wavelengths – Shape (num_bands,) array of wavelengths at which the SEDs are computed

Returns:

params

Shape (N,N,3) array of best-fitting modified blackbody parameters:
  • tau: optical depth at 500 micron

  • beta: spectral index

  • T: temperature

Return type:

np.ndarray

meshoid.radiation.modified_blackbody.modified_planck_function(freqs, logtau, beta, T)

Modified blackbody function:

I = (1-exp(-tau(f))) B(T)

where

B(T) = 2 h f^3/c^2 / (exp(h f/(k_B T) -1 ))

and tau(f) = tau0 (f/f0)^beta

Parameters:
  • freqs (array_like) – Photon frequencies in Hz

  • logtau (float) – log10 of optical depth at 500um

  • beta (float) – Dust spectral index beta

  • T (float) – Temperature in K

Returns:

I – Intensity of modified blackbody in erg/s/cm^2/sr

Return type:

array_like