radmc3dPy package

Submodules

radmc3dPy.analyze module

This module contains functions to read and write input/output data for RADMC-3D and to do some simple analysis/diagnostics of the model.

radmc3dPy.analyze.findContainerLeafID(cellCRD=None, cellHW=None, xi=None, yi=None, zi=None, childID=None, isLeaf=None, nChild=None, crd=None)

Function to find the tree index of a leaf cell containing a given point in space, i.e. if the following is true : xcell - dxcell <= xpoint < xcell + dxcell for each dimension. This function is to be used in multiprocessing.

Parameters:

cellCRD : ndarray

Array with dimensions [ncell, 3] containing the cell centre coordiantes of the tree

cellHW : ndarray

Array with dimensions [ncell, 3] containing the half width of cells in the tree

xi : ndarray

Array of cell interface indices in the base grid in the first dimension

yi : ndarray

Array of cell interface indices in the base grid in the second dimension

zi : ndarray

Array of cell interface indices in the base grid in the third dimension

childID : ndarray

Child index array

isLeaf : ndarray

Boolean array containing the node type for each cell (True - leaf, False - branch)

nChild : int

Number of children (8,4,2 for 3,2,1 active dimensions)

crd : ndarray

Array of length 3 containing the coordinates of the point whose container leaf is to be found

Returns:

Tree index of the container leaf if it is found. If the point is outside of the base grid -1 is returned.

radmc3dPy.analyze.findContainerLeafIDRec(x=None, y=None, z=None, dx=None, dy=None, dz=None, childID=None, isLeaf=None, nChild=None, crd=(), cellID=None)

Recursive function to find the leaf cell in the tree that contains a given point in space

Parameters:

x : ndarray

Tree cell center array in the first dimension

y : ndarray

Tree cell center array in the second dimension

z : ndarray

Tree cell center array in the tird dimension

dx : ndarray

Tree cell halfwidth array in the first dimension

dy : ndarray

Tree cell halfwidth array in the second dimension

dz : ndarray

Tree cell halfwidth array in the third dimension

childID : list

List of children indices. Each list element is an ndarray with nChild elements containing the child indices

isLeaf : ndarray

Boolean array for the cell type (True - leaf, False - branch)

nChild : int

Nr of children (i.e. 8, 4, or 2 for 3, 2, 1 active dimensions, respectively)

crd : ndarray

Three element list/tuple/array containing the point coordinates

cellID : int

Index of cell to be tested

radmc3dPy.analyze.gdensMinMax(x=None, y=None, z=None, dx=None, dy=None, dz=None, model=None, ppar=None, **kwargs)

Example function to be used as decision function for resolving cells in tree building. It calculates the gas density at a random sample of coordinates within a given cell than take the ratio of the max/min density. If it is larger than a certain threshold value it will return True (i.e. the cell should be resolved) if the density variation is less than the threshold it returns False (i.e. the cell should not be resolved)

Parameters:

x : ndarray

Cell centre coordinates of the cells in the first dimension

y : ndarray

Cell centre coordinates of the cells in the second dimension

z : ndarray

Cell centre coordinates of the cells in the third dimension

dx : ndarray

Half size of the cells in the first dimension

dy : ndarray

Half size of the cells in the second dimension

dz : ndarray

Half size of the cells in the third dimension

model : object

A radmc3dPy model (must contain a getGasDensity() function)

ppar : dictionary

All parameters of the problem (from the problem_params.inp file). It is not used here, but must be present for compatibility reasons.

**kwargs: dictionary

Parameters used to decide whether the cell should be resolved. It should the following keywords; ‘nsample’, which sets the number of random points the gas desity is sampled at within the cell and ‘threshold’ that sets the threshold value for max(gasdens)/min(gasdens) above which the cell should be resolved.

radmc3dPy.analyze.getDensVstruct(data=None, vmean_temp=False, ispec_tgas=0, gsize=None, idust=None, mstar=None, mu=None)

Calculates the vertical hydrostatic equilibrium

Parameters:

data : radmc3dData

An instance of the radmc3DData class containing the density structure of the model

vmean_temp : bool

If True (T(z) = T(-z) = 0.5*(T(z) + T(-z))) if False (T(z)!=T(-z))

idust : list

List of dust indices whose structure must be calculated

mstar : float

Stellar mass

ispec_tgas : int

Index of dust species whose temperature is taken to be the gas temperature

gsize : ndarray, optional

Dust grain sizes - If specified, the gas temperature is calculated as the average temperature of all dust grains in the grid cell weighted by the total surface area of dust grains with given size - NOTE: this approach assumes that all dust grains of a given size have the same bulk density

mu : float, optional

Mean molecular weight (default: 2.3)

Returns

——-

Returns an ndarray with the dust density

radmc3dPy.analyze.gmass(x=None, y=None, z=None, dx=None, dy=None, dz=None, model=None, ppar=None, **kwargs)

Example function to be used as decision function for resolving cells in tree building. It calculates the gas density at a random sample of coordinates within a given cell than take the ratio of the max/min density. If it is larger than a certain threshold value it will return True (i.e. the cell should be resolved) if the density variation is less than the threshold it returns False (i.e. the cell should not be resolved)

Parameters:

x : ndarray

Cell centre coordinates of the cells in the first dimension

y : ndarray

Cell centre coordinates of the cells in the second dimension

z : ndarray

Cell centre coordinates of the cells in the third dimension

dx : ndarray

Half size of the cells in the first dimension

dy : ndarray

Half size of the cells in the second dimension

dz : ndarray

Half size of the cells in the third dimension

model : object

A radmc3dPy model (must contain a getGasDensity() function)

ppar : dictionary

All parameters of the problem (from the problem_params.inp file). It is not used here, but must be present for compatibility reasons.

**kwargs: dictionary

Parameters used to decide whether the cell should be resolved. It should contain the following keywords; ‘nsample’, which sets the number of random points the gas desity is sampled at within the cell and ‘threshold’ that sets the threshold value for max(gasdens)/min(gasdens) above which the cell should be resolved.

radmc3dPy.analyze.interpolateOctree(data=None, x=None, y=None, z=None, var=None, nproc=1)

Nearest neighbour inteprolation on an octree

data : radmc3dData
Data container
x : ndarray
Coordiantes of the point to be interpolated on in the first dimension
y : ndarray
Coordiantes of the point to be interpolated on in the second dimension
z : ndarray
Coordiantes of the point to be interpolated on in the third dimension
var : list
Name of the variables to be interpolated, supported names are: ddens, dtemp, gdens, ndens, gtemp, gvel, vturb
nproc : int
Number of processes to be used (for parallel computing)
radmc3dPy.analyze.plotDustOpac(opac=None, var='kabs', idust=0, ax=None, xlabel=None, ylabel=None, fmt='', **kwargs)

Plots the dust opacity as a function of wavelength

Parameters:

opac : radmc3dDustOpac

Dust opacity container

var : {‘kabs’, ‘ksca’, ‘kext’, ‘g’}

Variable to be plotted

idust : int, optional

Dust index in opac to be plotted (default=0)

ax : Axis, optional

Matplotlib axis to plot on. If not set the current axis will be used.

xlabel : str, optional

Label of the x-axis

ylabel : str, optional

Label of the y-axis

fmt : str, optional

Format of the plotted line. The same as the third non-keyword argument of matplotlib.pyplot.plot()

Keyword Arguments:

Any further keyword argument that will be passed to matplotlib.pyplot.plot()

Returns:

The returned list by matplotlib.pyplot.plot()

radmc3dPy.analyze.plotScatmat(opac=None, var='z11', idust=0, iwav=None, wav=None, xvar='ang', iang=None, ang=None, ax=None, xlabel=None, ylabel=None, title=None, fmt='', **kwargs)

Plots the scattering matrix elements either as a function of scattering angle at a specific wavelength (default) or as a function of wavelength at a specific scattering angle

Parameters:

opac : radmc3dDustOpac

Dust opacity container

var : {‘kabs’, ‘ksca’, ‘kext’, ‘g’}

Variable to be plotted

idust : int, optional

Dust index in opac to be plotted (default=0)

iwav : int, optional

Wavelength index (used only if xvar=’ang’) to be plotted.

wav : float, optional

Wavelength at which the plot should be made (used only if xvar=’ang’). In practice, instead of interpolating to this wavelength, the nearest wavelength in the wavelength grid will be used.

xvar : {‘ang’, ‘wav’}

Variable for plotting the scattering matrix elements against (default=’ang’)

iang : int, optional

Scattering angle index (used only if xvar=’wav’)

ang : int, optional

Scattering angle in degrees at which the plot should be made (used only if xvar=’wav’). In practice, instead of interpolating to this scattering angle, the nearest angle in the angular grid will be used.

ax : Axis, optional

Matplotlib axis to plot on. If not set the current axis will be used.

xlabel : str, optional

Label of the x-axis

ylabel : str, optional

Label of the y-axis

title : str, optional

Title of the plot. If not set either the wavelength (for xvar=’ang’) or the scattering angle (for xvar=’wav’) the plot was made at.

fmt : str, optional

Format of the plotted line. The same as the third non-keyword argument of matplotlib.pyplot.plot()

Keyword Arguments:

Any further keyword argument that will be passed to matplotlib.pyplot.plot()

Returns:

The returned list by matplotlib.pyplot.plot()

radmc3dPy.analyze.plotSlice2D(data=None, var='ddens', plane='xy', crd3=0.0, icrd3=None, ispec=-1, xlim=(), ylim=(), log=False, linunit='cm', angunit='rad', nx=100, ny=100, showgrid=False, gridcolor='k', gridalpha=1.0, nproc=1, contours=False, clev=None, clmin=None, clmax=None, ncl=None, cllog=False, clcol='k', cllabel=False, cllabel_fontsize=10, cllabel_fmt='%.1f', clalpha=1.0, ax=None, lattitude=True, **kwargs)

Function to plot an axis-aligned 2D slice of the variables in the model. Any additional keyword argument above the listed ones will be passed on to matplotlib.pylab.pcolormesh(). For an octree grid the variables are interpolated onto a regular grid using nearest neighbour interpolation before plotting. The size and resolution of the regular image grid can be set at input.

Parameters:

data : radmc3dData

Instance of radmc3dData containing the field variable to be displayed

var : {‘ddens’, ‘dtemp’, ‘gdens’, ‘ndens’, ‘gtemp’, ‘vturb’, ‘vx’, ‘vy’, ‘vz’, ‘taux’, ‘tauy’}

Variable to be displayed

plane : {‘xy’, ‘xz’, ‘yz’, ‘yx, ‘zx’, ‘yz’}

Plane to be displayed

crd3 : float

Coordinate of the third dimension (i.e. when plotting a slice in the x-y plane, crd3 is the z-coordinate)

icrd3 : int

Index of the third coordinate in the grid (only for regular grid!)

ispec : int

Dust species index. If negative dust densities will be summed up and the total cumulative density will be displayed

xlim : tuple

Coordinate boundaries in the first dimension of the plot (also the coordinate boundary of the regular grid data on AMR grids are interpolated to)

ylim : tuple

Coordinate boundaries in the second dimension of the plot (also the coordinate boundary of the regular grid data on AMR grids are interpolated to)

log : bool

If True the contour/image will be displayed on a logarithmic stretch

linunit : {‘cm’, ‘au’, ‘pc’, ‘rs’}

Unit selection for linear image coordinate axes.

nx : int

Number of horizontal pixels in the interpolated image if the data is defined in an Octree

ny : int

Number of vertical pixels in the interpolated image if the data is defined in an Octree

showgrid : bool

If True the spatial grid will be overlayed

gridcolor : str

Color of the spatial grid overlay

gridalpha : float

Opacity of the lines in the spatial grid overlay (0.0 - fully transparent, 1.0 - fully opaque)

angunit : {‘rad’, ‘deg’}

Unit selection for angular image coordinate axes (only if spherical coordinate system is used).

nproc : int

Number of parallel processes to be used for interpolation.

contours : bool

If True contour lines are plotted, if False a colorscale plot will be created

clev : ndarray

A numpy ndarray containing the levels to be displayed with contour lines. If clev is set then clmin, clmax and ncl are omitted

clmin : float

Min. contour level (for setting auto-contours between clmin and clmax at ncl values)

clmax : float

Max. contour level (for setting auto-contours between clmin and clmax at ncl values)

ncl : float

Number of contour levels (for setting auto-contours between clmin and clmax at ncl values)

cllog : bool

If clmin, clmax and ncl are used to generate the contour levels, then if cllog is True the contours will be log-scaled

clcol : str

Color-code for the contour lines for single color contours

cllabel : bool

If True the contour line values will be displayed, if False only the contour lines will be displayed (default = False)

cllabel_fontsize: int

Size of the font used to displaye the contour line values

cllabel_fmt : str

Format of the contour line labels (default “%.1f”)

clalpha : float

Transparency of the contour lines (1.0 fully opaque, 0.0 fully transparent)

lattitude : bool

If the coordinate sytem used in RADMC-3D is spherical, then the 2nd coordiante is the co-lattitude. If lattitude is set to True then the 2nd coordinate in the RADMC-3D grid will be transformet to true lattitude (i.e. pi/2.-colattitude). If set to false the original co-lattitude will be used.

ax : matplotlib.axes.Axes

Matplotlib axis to plot to

Keyword Arguments :

All other keyword arugments will be passed to pcolormesh() or countour()

radmc3dPy.analyze.plotSpectrum(a, ev=False, kev=False, micron=False, jy=False, lsun=False, lnu=False, nulnu=False, fnu=False, nufnu=False, dpc=1.0, oplot=False, xlg=False, ylg=False, obs=False, mol=None, ilin=None)

Plot the spectrum / SED

Parameters:

a : ndarray

A 2D array of size [Nfreq,2] returned by readSpectrum(). [:,0] - wavelength in micrometer, or for line data the velocity in km/s [:,1] - flux density in erg/s/cm/cm/Hz

ev : bool

True –> energy in electronvolt (default=Hz)

kev : bool

True –> energy in kiloelectronvolt (default=Hz)

micron : bool

True –> wavelength in micron (default=Hz)

jy : bool

True –> Flux in Jansky

lnu : bool

True –> L_nu (default L_nu)

fnu : bool

True –> F_nu in units of erg/s/cm^2/Hz(default L_nu)

nufnu : bool

True –> nu*F_nu in units of erg/s/cm^2 (default L_nu)

nulnu : bool

True –> nu*L_nu (default F_nu)

lsun : bool

True –> nu*L_nu in units of solar luminosity

dpc : bool

Distance of observer in units of parsec (Default: 1 pc)

oplot : bool

True –> Plot without refreshing subplot

xlg : bool

True –> logarithmic x-axis

ylg : bool

True –> logarithmic y-axis

obs : bool

True –> Treat the spectrum as an observation

(i.e. do not scale with dpc^(-2))

mol : radmc3dMolecule

(optional) Molecule data (see radmc3dMolecule class)

This is required if you want to plot a line spectrum with on the x-axis the radial velocity in km/s

ilin : bool

(if set) the index of the line (of mol; starting,

as in RADMC-3D, with the index 1) which shall act as the 0 km/s wavelength reference. If ilin is set the x axis will be in km/s (overriding other settings)

radmc3dPy.analyze.readData(ddens=False, dtemp=False, gdens=False, gtemp=False, gvel=False, ispec=None, vturb=False, grid=None, binary=True, old=False, octree=False)

Reads the physical variables of the model (e.g. density, velocity, temperature).

Parameters:

ddens : bool

If True dust density will be read (all dust species and grain sizes)

dtemp : bool

If True dust temperature will be read (all dust species and grain sizes)

gdens : bool

If True gas density will be read (NOTE: the gas density will be number density in 1/cm^3)

gtemp : bool

If True gas temperature will be read (all dust species and grain sizes)

gvel : bool

If True the velocity field will be read

ispec : str

Name of the molecule in the ‘molecule_ispec.inp’ filename

vturb : bool

If True the microturbulent velocity field will be read

grid : radmc3dGrid

An instance of radmc3dGrid containing the spatial and frequency grid of the model. If the grid is passed to the function it will not be read again from file. This can be useful for octree models to save time.

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

binary: bool

Set it to True for C-style binary and False for formatted ASCII files

octree: bool

True for models with octree AMR and False for models with regular grid

Returns:

Returns an instance of the radmc3dData class

radmc3dPy.analyze.readGrid(sgrid=True, wgrid=True, sgrid_fname=None, wgrid_fname=None, old=False)

Reads the spatial and frequency grid. This function is an interface to radmc3dGrid.readGrid().

Parameters:

sgrid : bool

If True the spatial grid will be read

wgrid : bool

If True the wavelength grid will be read

sgrid_fname : str

File containing the spatial grid (default: amr_grid.inp)

wgrid_fname : str

File containing the wavelength grid (default: wavelength_micron.inp)

old : bool

If True the format of the old 2D version of radmc3d (radmc) will be used

Returns

——-

Returns an instance of the radmc3dGrid (for regular grid) or radmc3dOctree (for octree AMR) class

radmc3dPy.analyze.readMol(mol=None, fname=None)

Wrapper around the radmc3dMolecule.read() method

Parameters:

mol : str

molecule name (e.g. ‘co’) if the file name is in the form of ‘molecule_<mol>.inp’

fname : str

full file name

radmc3dPy.analyze.readOpac(ext=None, idust=None, scatmat=None, old=False)

Reads the dust opacity files. This function is an interface to radmc3dDustOpac.readOpac()

Parameters:

ext : list

Each element of the list is be a string, the file name extension (file names should look like ‘dustkappa_ext.inp’)

idust : list

Each element of the list is an integer, the index of the dust species in the master opacity file (dustopac.inp’)

scatmat: list

If specified, its elements should be booleans indicating whether the opacity file contains also the full scattering matrix (True) or only dust opacities (False)

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

Returns:

Returns an instance of the radmc3dDustOpac class

radmc3dPy.analyze.readParams()

Reads the problem_params.inp file. This function is an interface to radmc3dPar.readPar().

Returns:Returns an instance of the radmc3dPar class
radmc3dPy.analyze.readSpectrum(fname='', old=False)

Reads the spectrum / SED

Parameters:

fname : str, optional

Name of the file to be read

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

Returns:

Returns an ndarray with [Nwavelength, 2] dimensions

[Nwavelength,0] is the wavelength / velocity and [Nwavelength,1] is the flux density

radmc3dPy.analyze.readStars(fname='')

Reads the data (mass, radius, temperature, spectrum) of discrete stellar sources

Parameters:

fname : str

Name of the file to be read (if omitted the default value is stars.inp)

Returns:

An instance of radmc3dRadSources containing the stellar data

radmc3dPy.analyze.writeDefaultParfile(model='', fname='')

Writes a parameter file (problem_params.inp) with default parameters for a given model.

Parameters:

model : str

Name of the model whose parameter should be written to the file

fname : str, optional

Name of the parameter file to be written (if omitted problem_params.inp will be used)

radmc3dPy.conf module

radmc3dPy.crd_trans module

This module contains functions for coordinate transformations (e.g. rotation). For help on the syntax or functionality of each function see the help of the individual functions

radmc3dPy.crd_trans.csrot(crd=None, ang=None, xang=0.0, yang=0.0, zang=0.0, deg=False)

Performs coordinate system rotation.

Parameters:

crd : numpy ndarray

Three element vector containing the coordinates of a given point in a cartesian system

ang : list, ndarray

Three element list/ndarray describing the rotation angles around the x, y and z axes, respectively

xang: float

Rotation around the x-axis

yang: float

Rotation around the y-axis

zang: float

Rotation around the z-axis

deg : float, optional

If True angles should be given in degree instead of radians (as by default)

Returns:

list

Returns a three element list with the rotated coordinates

Notes

Rotation matrices

Around the x-axis:

\[\begin{split}\left(\begin{matrix} 1 & 0 & 0 \\ 0 & cos(\alpha) & -sin(\alpha)\\ 0 & sin(\alpha) & cos(\alpha) \end{matrix}\right)\end{split}\]

Around the y-axis:

\[\begin{split}\left(\begin{matrix} cos(\beta) & 0 & -sin(\beta) \\ 0 & 1 & 0\\ sin(\beta)& 0 & cos(\beta) \end{matrix}\right)\end{split}\]

Around the z-axis

\[\begin{split}\left(\begin{matrix} cos(\gamma) & -sin\gamma) & 0 \\ sin(\gamma) & cos(\gamma) & 0 \\ 0 & 0 & 1 \end{matrix}\right)\end{split}\]
radmc3dPy.crd_trans.ctransSph2Cart(crd=None, reverse=False)

Transform coordinates between spherical to cartesian systems

Parameters:

crd : ndarray

Three element array containing the input coordinates [x,y,z] or [r,theta,phi] by default the coordinates assumed to be in the cartesian system

reverse : bool

If True calculates the inverse transformation (cartesian -> spherical). In this case crd should be [r,theta,phi]

Returns:

Returns a three element array containig the output coordinates [r,theta,phi] or [x,y,z]

radmc3dPy.crd_trans.vrot(crd=None, v=None, ang=None)

Rotates a vector in spherical coordinate system. First transforms the vector to cartesian coordinate system, then does the rotation then makes the inverse transformation

Parameters:

crd : ndarray

Three element array containing the coordinates of a given point in the cartesian system

v : ndarray

Three element array, angles of rotation around the x,y,z axes

ang : ndarray

Three element arrray containing the angles to rotate around the x, y, z, axes, respectively

radmc3dPy.crd_trans.vtransSph2Cart(crd=None, v=None, reverse=False)

Transform velocities between spherical to cartesian systems

Parameters:

crd : ndarray

Three element array containing the input coordinates [x,y,z] or [r,theta,phi] by default the coordinates assumed to be in the cartesian system

v : ndarray

Three element array containing the input velocities in the same coordinate system as crd

reverse : bool

If True it calculates the inverse trasnformation (cartesian -> spherical)

Returns:

Returns a three element array containg the output velocities [vr,vphi,vtheta] or [vx,vy,vz]

radmc3dPy.data module

This module contains a class for handling variable data in radmc-3d

class radmc3dPy.data.radmc3dData(grid=None)

Bases: object

RADMC-3D data class.
Reading and writing dust density/temperature, gas density/temperature/velocity, generating a legacy vtk file for visualization.

Attributes

grid (radmc3dGrid, radmc3dOctree) Instance of the radmc3dGrid class, contains the spatial and frequency grids
rhodust (ndarray) Dust density in g/cm^3
dusttemp (ndarray) Dust temperature in K
rhogas (ndarray) Gas density in g/cm^3
ndens_mol (ndarray) Number density of the molecule [molecule/cm^3]
ndens_cp (ndarray) Number density of the collisional partner [molecule/cm^3]
gasvel (ndarray) Gas velocity in cm/s
gastemp (ndarray) Gas temperature in K
vturb (ndarray) Mictroturbulence in cm/s
taux (ndarray) Optical depth along the x (cartesian) / r (cylindrical) / r (spherical) dimension
tauy (ndarray) Optical depth along the y (cartesian) / theta (cylindrical) / theta (spherical) dimension
tauz (ndarray) Optical depth along the z (cartesian) / z (cylindrical) / phi (spherical) dimension
sigmadust (ndarray) Dust surface density in g/cm^2
sigmagas (ndarray) Gas surface density in molecule/cm^2 (or g/cm^2 depending on the dimension of rhogas)

Methods

getDustMass([idust]) Calculates the dust mass in radmc3dData.rhodust
getGasMass([mweight, rhogas]) Calculates the gas mass in radmc3dData.ndens_mol or radmc3dData.rhogas
getSigmaDust([idust]) Calculates the dust surface density.
getSigmaGas() Calculates the gas surface density.
getTau([idust, axis, wav, kappa, old]) Calculates the optical depth along any given combination of the axes.
getTauOneDust([idust, axis, kappa]) Calculates the optical depth of a single dust species along any given combination of the axes.
readDustDens([fname, binary, old, octree]) Reads the dust density.
readDustTemp([fname, binary, octree, old]) Reads the dust temperature.
readGasDens([fname, ispec, binary, octree]) Reads the gas density.
readGasTemp([fname, binary, octree]) Reads the gas temperature.
readGasVel([fname, binary, octree]) Reads the gas velocity.
readVTurb([fname, binary, octree]) Reads the turbulent velocity field.
writeDustDens([fname, binary, old, octree]) Writes the dust density.
writeDustTemp([fname, binary, octree]) Writes the dust density.
writeGasDens([fname, ispec, binary, octree]) Writes the gas density.
writeGasTemp([fname, binary, octree]) Writes the gas temperature.
writeGasVel([fname, binary, octree]) Writes the gas velocity.
writeVTK([vtk_fname, ddens, dtemp, idust, ...]) Writes physical variables to a legacy vtk file.
writeVTurb([fname, binary, octree]) Writes the microturbulence file.
getDustMass(idust=-1)

Calculates the dust mass in radmc3dData.rhodust

Parameters:

idust : int

Dust index whose dust should be calculated. If it is set to -1 (default) the total dust mass is calculated summing up all dust species

Returns:

A single float being the dust mass in gramm

getGasMass(mweight=2.3, rhogas=False)

Calculates the gas mass in radmc3dData.ndens_mol or radmc3dData.rhogas

Parameters:

mweight : float

Molecular weight [atomic mass unit / molecule, i.e. same unit as mean molecular weight]

rhogas : bool, optional

If True the gas mass will be calculated from radmc3dData.rhogas, while if set to False the gas mass will be calculated from radmc3dData.ndens_mol. The mweight parameter is only required for the latter.

Returns:

A single float being the gas mass in gramm

getSigmaDust(idust=-1)

Calculates the dust surface density.

Parameters:

idust : int, optional

Index of the dust species for which the surface density should be calculated if omitted the calculated surface density will be the sum over all dust species

getSigmaGas()

Calculates the gas surface density. This method uses radmc3dData.rhogas to calculate the surface density, thus the unit of surface density depends on the unit of radmc3dData.rhogas (g/cm^2 or molecule/cm^2)

getTau(idust=None, axis='xy', wav=0.0, kappa=None, old=False)

Calculates the optical depth along any given combination of the axes.

Parameters:

idust : list

List of dust component indices whose optical depth should be calculated If multiple indices are set the total optical depth is calculated summing over all dust species in idust

axis : str

Name of the axis/axes along which the optical depth should be calculated (e.g. ‘x’ for the first dimension or ‘xyz’ for all three dimensions)

wav : float

Wavelength at which the optical depth should be calculated

kappa : list, tuple

If set it should be a list of mass extinction coefficients at the desired wavelength The number of elements in the list should be equal to that in the idust keyword

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

getTauOneDust(idust=0, axis='', kappa=0.0)

Calculates the optical depth of a single dust species along any given combination of the axes.

Parameters:

idust : int

Index of the dust species whose optical depth should be calculated

axis : str

Name of the axis/axes along which the optical depth should be calculated (e.g. ‘x’ for the first dimension or ‘xyz’ for all three dimensions)

kappa : float

Mass extinction coefficients of the dust species at the desired wavelength

Returns:

Returns a dictionary with the following keys

taux : ndarray

optical depth along the first dimension

tauy : ndarray

optical depth along the second dimension

(tauz is not yet implemented)

readDustDens(fname=None, binary=True, old=False, octree=False)

Reads the dust density.

Parameters:

fname : str, optional

Name of the file that contains the dust density. If omitted ‘dust_density.inp’ is used (or if binary=True the ‘dust_density.binp’ is used).

binary : bool, optional

If true the data will be read in binary format, otherwise the file format is ascii

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

octree : bool, optional

If the data is defined on an octree-like AMR

readDustTemp(fname=None, binary=True, octree=False, old=False)

Reads the dust temperature.

Parameters:

fname : str, optional

Name of the file that contains the dust temperature.

binary : bool, optional

If true the data will be read in binary format, otherwise the file format is ascii

octree : bool, optional

If the data is defined on an octree-like AMR

old : bool, optional

If True dust temperature will be written in the old RADMC format

readGasDens(fname=None, ispec='', binary=True, octree=False)

Reads the gas density.

Parameters:

ispec : str

File name extension of the ‘numberdens_ispec.inp’ (or if binary=True ‘numberdens_ispec.binp’) file.

binary : bool

If true the data will be read in binary format, otherwise the file format is ascii

octree : bool, optional

If the data is defined on an octree-like AMR

readGasTemp(fname=None, binary=True, octree=False)

Reads the gas temperature.

Parameters:

fname : str,optional

Name of the file that contains the gas temperature. If omitted ‘gas_temperature.inp’ (or if binary=True ‘gas_tempearture.binp’) is used.

binary : bool

If true the data will be read in binary format, otherwise the file format is ascii

octree : bool, optional

If the data is defined on an octree-like AMR

readGasVel(fname=None, binary=True, octree=False)

Reads the gas velocity.

Parameters:

fname : str, optional

Name of the file that contains the gas velocity If omitted ‘gas_velocity.inp’ (if binary=True ‘gas_velocity.binp’)is used.

binary : bool

If true the data will be read in binary format, otherwise the file format is ascii

octree : bool, optional

If the data is defined on an octree-like AMR

readVTurb(fname=None, binary=True, octree=False)

Reads the turbulent velocity field.

Parameters:

fname : str, optional

Name of the file that contains the turbulent velocity field If omitted ‘microturbulence.inp’ (if binary=True ‘microturbulence.binp’) is used.

binary : bool

If true the data will be read in binary format, otherwise the file format is ascii

octree : bool, optional

If the data is defined on an octree-like AMR

writeDustDens(fname='', binary=True, old=False, octree=False)

Writes the dust density.

Parameters:

fname : str, optional

Name of the file into which the dust density should be written. If omitted ‘dust_density.inp’ is used.

binary : bool

If true the data will be written in binary format, otherwise the file format is ascii

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

octree : bool, optional

If the data is defined on an octree-like AMR

writeDustTemp(fname='', binary=True, octree=False)

Writes the dust density.

Parameters:

fname : str, optional

Name of the file into which the dust density should be written. If omitted ‘dust_density.inp’ is used.

binary : bool

If true the data will be written in binary format, otherwise the file format is ascii

octree : bool, optional

If the data is defined on an octree-like AMR

writeGasDens(fname=None, ispec='', binary=True, octree=False)

Writes the gas density.

Parameters:

fname : str, optional

Name of the file into which the data will be written. If omitted “numberdens_xxx.inp” and “numberdens_xxx.binp” will be used for ascii and binary format, respectively (xxx is the name of the molecule).

ispec : str

File name extension of the ‘numberdens_ispec.inp’ (if binary=True ‘numberdens_ispec.binp’) file into which the gas density should be written

binary : bool

If true the data will be written in binary format, otherwise the file format is ascii

octree : bool, optional

If the data is defined on an octree-like AMR

writeGasTemp(fname='', binary=True, octree=False)

Writes the gas temperature.

Parameters:

fname : str, optional

Name of the file into which the gas temperature should be written. If omitted ‘gas_temperature.inp’ (if binary=True ‘gas_tempearture.binp’) is used.

binary : bool

If true the data will be written in binary format, otherwise the file format is ascii

octree : bool, optional

If the data is defined on an octree-like AMR

writeGasVel(fname='', binary=True, octree=False)

Writes the gas velocity.

Parameters:

fname : str, optional

Name of the file into which the gas temperature should be written. If omitted ‘gas_velocity.inp’ (if binary=True ‘gas_velocity.binp’) is used.

binary : bool

If true the data will be written in binary format, otherwise the file format is ascii

octree : bool, optional

If the data is defined on an octree-like AMR

writeVTK(vtk_fname='', ddens=False, dtemp=False, idust=None, gdens=False, gvel=False, gtemp=False)

Writes physical variables to a legacy vtk file.

Parameters:

vtk_fname : str

Name of the file to be written, if not specified ‘radmc3d_data.vtk’ will be used

ddens : bool

If set to True the dust density will be written to the vtk file

dtemp : bool

If set to True the dust temperature will be written to the vtk file

idust : list

List of indices that specifies which dust component should be written if not set then the first dust species (zero index) will be used

gdens : bool

If set to True the gas density will be written to the vtk file

gtemp : bool

If set to True the gas temperature will be written to the vtk file

gvel : bool

If set to True the gas velocity will be written to the vtk file

writeVTurb(fname='', binary=True, octree=False)

Writes the microturbulence file.

Parameters:

fname : str, optional

Name of the file into which the turubulent velocity field should be written. If omitted ‘microturbulence.inp’ (if binary=True ‘microturbuulence.binp’) is used.

binary : bool

If true the data will be written in binary format, otherwise the file format is ascii

octree : bool, optional

If the data is defined on an octree-like AMR

radmc3dPy.dustopac module

This module contains classes for handling dust opacities

radmc3dPy.dustopac.computeDustOpacMie(fname='', matdens=None, agraincm=None, lamcm=None, theta=None, logawidth=None, wfact=3.0, na=20, chopforward=0.0, errtol=0.01, verbose=False, extrapolate=False, return_type=1)

Compute dust opacity with Mie theory based on the optical constants in the optconst_file. Optionally also the scattering phase function in terms of the Mueller matrix elements can be computed. To smear out the resonances that appear due to the perfect sphere shape, you can optionally smear out the grain size distribution a bit with setting the width of a Gaussian grain size distribution.

Parameters:

fname : str

File name of the optical constants file. This file should contain three columns: first the wavelength in micron, then the n-coefficient and then the k-coefficient. See Jena optical constants database: http://www.astro.uni-jena.de/Laboratory/Database/databases.html

matdens : float

Material density in g/cm^3

agraincm : float

Grain radius in cm

lamcm : ndarray

Wavelength grid in cm

theta : ndarray, optional

Angular grid (a numpy array) between 0 and 180 which are the scattering angle sampling points at which the scattering phase function is computed.

logawidth : float, optional

If set, the size agrain will instead be a sample of sizes around agrain. This helps to smooth out the strong wiggles in the phase function and opacity of spheres at an exact size. Since in Nature it rarely happens that grains all have exactly the same size, this is quite natural. The value of logawidth sets the width of the Gauss in ln(agrain), so for logawidth<<1 this give a real width of logawidth*agraincm.

wfact : float

Grid width of na sampling points in units of logawidth. The Gauss distribution of grain sizes is cut off at agrain * exp(wfact*logawidth) and agrain * exp(-wfact*logawidth). Default = 3

na : int

Number of size sampling points (if logawidth set, default=20)

chopforward : float

If >0 this gives the angle (in degrees from forward) within which the scattering phase function should be kept constant, essentially removing the strongly peaked forward scattering. This is useful for large grains (large ratio 2*pi*agraincm/lamcm) where the forward scattering peak is extremely strong, yet extremely narrow. If we are not interested in very forward-peaked scattering (e.g. only relevant when modeling e.g. the halo around the moon on a cold winter night), this will remove this component and allow a lower angular grid resolution for the theta grid.

errtol : float

Tolerance of the relative difference between kscat and the integral over the zscat Z11 element over angle. If this tolerance is exceeded, a warning is given.

verbose : bool

If set to True, the code will give some feedback so that one knows what it is doing if it becomes slow.

extrapolate : bool

If set to True, then if the wavelength grid lamcm goes out of the range of the wavelength grid of the optical constants file, then it will make a suitable extrapolation: keeping the optical constants constant for lamcm < minimum, and extrapolating log-log for lamcm > maximum.

return_type : {0, 1}

If 0 a dictionary is returned (original return type) if 1 an instance of radmc3dDustOpac will be returned

Returns:

A dictionary with the following keys:

  • kabs : ndarray

    Absorption opacity kappa_abs_nu (a numpy array) in units of cm^2/gram

  • ksca : ndarray

    Scattering opacity kappa_abs_nu (a numpy array) in units of cm^2/gram

  • gsca : ndarray

    The <cos(theta)> g-factor of scattering

  • theta : ndarray (optional, only if theta is given at input)

    The theta grid itself (just a copy of what was given)

  • zscat : ndarray (optional, only if theta is given at input)

    The components of the scattering Mueller matrix Z_ij for each wavelength and each scattering angel. The normalization of Z is such that kscat can be reproduced (as can be checked) by the integral: 2*pi*int_{-1}^{+1}Z11(mu)dmu=kappa_scat. For symmetry reasons only 6 elements of the Z matrix are returned: Z11, Z12, Z22, Z33, Z34, Z44. Note that Z21 = Z12 and Z43 = -Z34. The scattering matrix is normalized such that if a plane wave with Stokes flux

    Fin = (Fin_I,Fin_Q,Fin_U,Fin_V)

    hits a dust grain (which has mass mgrain), then the scattered flux

    Fout = (Fout_I,Fout_Q,Fout_U,Fout_V)

    at distance r from the grain at angle theta is given by

    Fout(theta) = (mgrain/r^2) * Zscat . Fin

    where . is the matrix-vector multiplication. Note that the Stokes components must be such that the horizontal axis in the “image” is pointing in the scattering plane. This means that radiation with Fin_Q < 0 is scattered well, because it is vertically polarized (along the scattering angle axis), while radiation with Fin_Q > 0 is scatterd less well because it is horizontally polarized (along the scattering plane).

  • kscat_from_z11 : ndarray (optional, only if theta is given at input)

    The kscat computed from the (above mentioned) integral of Z11 over all angles. This should be nearly identical to kscat if the angular grid is sufficiently fine. If there are strong differences, this is an indication that the angular gridding (the theta grid) is not fine enough. But you should have then automatically gotten a warning message as well (see errtol).

  • wavmic : ndarray (optional, only if extrapolate is set to True)

    The original wavelength grid from the optical constants file, with possibly an added extrapolated

  • ncoef : ndarray (optional, only if extrapolate is set to True)

    The optical constant n at that grid

  • kcoef : ndarray (optional, only if extrapolate is set to True)

    The optical constant k at that grid

  • agr : ndarray (optional, only if logawidth is not None)

    Grain sizes

  • wgt : ndarray (optional, only if logawidth is not None)

    The averaging weights of these grain (not the masses!) The sum of wgt.sum() must be 1.

  • zscat_nochop : ndarray (optional, only if chopforward > 0)

    The zscat before the forward scattering was chopped off

  • kscat_nochop : ndarray (optional, only if chopforward > 0)

    The kscat originally from the bhmie code

class radmc3dPy.dustopac.radmc3dDustOpac

Bases: object

Class to handle dust opacities.

Attributes

wav (list) Each element of the list contains an ndarray with the wavelength grid
freq (list) Each element of the list contains an ndarray with the frequency grid
nwav (list) Each element of the list contains an integer with the number of wavelengths
kabs (list) Each element of the list contains an ndarray with the absorption coefficient per unit mass
ksca (list) Each element of the list contains an ndarray with the scattering coefficient per unit mass
phase_g (list) Each element of the list contains an ndarray with the hase function
ext (list) Each element of the list contains a string wht the file name extension of the duskappa_ext.Kappa file
therm (list) Each element of the list contains a bool, if it is set to False the dust grains are quantum-heated (default: True)
idust (lisintt) Each element of the list contains an integer with the index of the dust species in the dust density distribution array
scatmat (list) Each element is a boolean indicating whether the dust opacity table includes (True) the full scattering matrix or not (False)
nang (list) Each element is a string, containing the number of scattering angles in the scattering matrix if its given
scatang (list) Each element is a numpy ndarray containing the scattering angles in the scattering matrix if its given
z11 (list) Each element is a numpy ndarray containing the (1,1) element of the scattering angles in the scattering matrix if its given
z12 (list) Each element is a numpy ndarray containing the (1,2) element of the scattering angles in the scattering matrix if its given
z22 (list) Each element is a numpy ndarray containing the (2,2) element of the scattering angles in the scattering matrix if its given
z33 (list) Each element is a numpy ndarray containing the (3,3) element of the scattering angles in the scattering matrix if its given
z34 (list) Each element is a numpy ndarray containing the (3,4) element of the scattering angles in the scattering matrix if its given
z44 (list) Each element is a numpy ndarray containing the (4,4) element of the scattering angles in the scattering matrix if its given

Methods

makeOpac([ppar, wav, old, code, theta, ...]) Createst the dust opacities using a Mie code distributed with RADMC-3D.
makeopacRadmc2D([ext]) Creates dust opacities (dustopac_*.inp files) for the previous 2D version of radmc
mixOpac([ppar, mixnames, mixspecs, mixabun, ...]) Mixes dust opacities.
readMasterOpac() Reads the master opacity file ‘dustopac.inp’.
readOpac([ext, idust, scatmat, old]) Reads the dust opacity files.
runMakedust([freq, gmin, gmax, ngs, ...]) Interface function to the F77 code makedust to calculate mass absorption coefficients.
writeMasterOpac([ext, therm, ...]) Writes the master opacity file ‘dustopac.inp’.
writeOpac([fname, ext, idust, scatmat]) Writes dust opacities to file
makeOpac(ppar=None, wav=None, old=False, code='python', theta=None, logawidth=None, wfact=3.0, na=20, chopforward=0.0, errtol=0.01, verbose=False, extrapolate=False)

Createst the dust opacities using a Mie code distributed with RADMC-3D.

Parameters:

ppar : dictionary

Parameters of the simulations

wav : ndarray, optional

Wavelength grid on which the mass absorption coefficients should be calculated

code : {‘python’, ‘fortran’}

Version of the mie scattering code BHMIE to be used. ‘fortran’ - use the original fortran77 code of Bruce Drain (should be downloaded separately, compiled and its path added to the PATH environment variable), ‘python’ a python version of BHMIE by Kees Dullemond (radmc3dPy.miescat).

theta : ndarray, optional

Angular grid (a numpy array) between 0 and 180 which are the scattering angle sampling points at which the scattering phase function is computed.

logawidth : float, optional

If set, the size agrain will instead be a sample of sizes around agrain. This helps to smooth out the strong wiggles in the phase function and opacity of spheres at an exact size. Since in Nature it rarely happens that grains all have exactly the same size, this is quite natural. The value of logawidth sets the width of the Gauss in ln(agrain), so for logawidth<<1 this give a real width of logawidth*agraincm.

wfact : float

Grid width of na sampling points in units of logawidth. The Gauss distribution of grain sizes is cut off at agrain * exp(wfact*logawidth) and agrain * exp(-wfact*logawidth). Default = 3

na : int

Number of size sampling points (if logawidth set, default=20)

chopforward : float

If >0 this gives the angle (in degrees from forward) within which the scattering phase function should be kept constant, essentially removing the strongly peaked forward scattering. This is useful for large grains (large ratio 2*pi*agraincm/lamcm) where the forward scattering peak is extremely strong, yet extremely narrow. If we are not interested in very forward-peaked scattering (e.g. only relevant when modeling e.g. the halo around the moon on a cold winter night), this will remove this component and allow a lower angular grid resolution for the theta grid.

errtol : float

Tolerance of the relative difference between kscat and the integral over the zscat Z11 element over angle. If this tolerance is exceeded, a warning is given.

verbose : bool

If set to True, the code will give some feedback so that one knows what it is doing if it becomes slow.

extrapolate : bool

If set to True, then if the wavelength grid lamcm goes out of the range of the wavelength grid of the optical constants file, then it will make a suitable extrapolation: keeping the optical constants constant for lamcm < minimum, and extrapolating log-log for lamcm > maximum.

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

makeopacRadmc2D(ext=None)

Creates dust opacities (dustopac_*.inp files) for the previous 2D version of radmc It takes the input dust opacity files and interpolates them onto the used frequency grid

Parameters:

ext : list

List of dustkappa file name extensions, i.e. the input file name has to be named as dustkappa_ext[i].inp

static mixOpac(ppar=None, mixnames=None, mixspecs=None, mixabun=None, writefile=True)

Mixes dust opacities.

Parameters:

ppar : dictionary, optional

All parameters of the actual model setup.

mixnames : list, optional

Names of the files into which the mixed dust opacities will be written (not needed if writefile=False)

mixspecs : list, optional

Names of the files from which the dust opacities are read (not needed if readfile=False)

mixabun : list, optional

Abundances of different dust species

writefile : bool

If False the mixed opacities will not be written out to files given in mixnames.

NOTE, either ppar or mixname, mixspecs, and mixabun should be set.

static readMasterOpac()

Reads the master opacity file ‘dustopac.inp’. It reads the dustkappa filename extensions (dustkappa_ext.inp) corresponding to dust species indices

Returns:

Returns a dictionary with the following keys:

*ext : list of dustkappa file name extensions

*therm : a list of integers specifying whether the dust grain is thermal or quantum heated (0 - thermal, 1 - quantum heated)

readOpac(ext=None, idust=None, scatmat=None, old=False)

Reads the dust opacity files.

Parameters:

ext : list

File name extension (file names should look like ‘dustkappa_ext.inp’)

idust: list

Indices of the dust species in the master opacity file (dustopac.inp’) - starts at 0

scatmat: list

If specified, its elements should be booleans indicating whether the opacity file contains also the full scattering matrix (True) or only dust opacities (False)

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

static runMakedust(freq=None, gmin=None, gmax=None, ngs=None, lnk_fname=None, gdens=None)

Interface function to the F77 code makedust to calculate mass absorption coefficients.

Parameters:

freq : ndarray

Contains the frequency grid on which the opacities should be calculated

gmin : float

Minimum grain size

gmax : float

Maximum grain size

ngs : int

Number of grain sizes

gdens : float

Density of the dust grain in g/cm^3

lnk_fname : str

Name of the file in which the optical constants are stored

Returns:

Returns an ndarray with [nfreq,ngs] dimensions containing the resulting opacities

static writeMasterOpac(ext=None, therm=None, scattering_mode_max=1, old=False)

Writes the master opacity file ‘dustopac.inp’.

Parameters:

ext : list

List of dustkappa file name extensions

therm : list

List of integers specifying whether the dust grain is thermal or quantum heated (0-thermal, 1-quantum)

scattering_mode_max : int

Scattering mode code in radmc3d : 0 - no scattering, 1 - isotropic scattering, 2 - anisotropic scattering with Henyei-Greenstein phase function, 5 - anisotropic scattering using the full scattering matrix and stokes vectors.

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

writeOpac(fname=None, ext=None, idust=None, scatmat=False)

Writes dust opacities to file

Parameters:

fname : str

Name of the file to write the dust opacties into

ext : str

If fname is not specified, the output file name will be generated as dustkappa_EXT.inp or dustkapscatmat_EXT.inp depending on the file format

idust : int

Dust species index whose opacities should be written to file

scatmat : bool

If True the full scattering matrix will be written to file on top of the opacities (i.e. the file name should be dustkapscatmat_EXT.inp). If False only the dust opacities and the asymmetry parameter (if present) will be written to file (dustkappa_EXT.inp type files)

radmc3dPy.image module

This module contains classes/functions to create and read images with RADMC-3D and to calculate interferometric visibilities and write fits files For help on the syntax or functionality of each function see the help of the individual functions

radmc3dPy.image.cmask(im=None, rad=0.0, au=False, arcsec=False, dpc=None)
Simulates a coronographic mask.
Sets the image values to zero within circle of a given radius around the image center.
Parameters:

im : radmc3dImage

A radmc3dImage class containing the image

rad : float

The raadius of the mask

au : bool

If true the radius is taken to have a unit of AU

arcsec : bool

If true the radius is taken to have a unit of arcsec (dpc should also be set)

dpc : float

Distance of the source (required if arcsec = True)

NOTE if arcsec=False and au=False rad is taken to have a unit of pixel

Returns:

Returns a radmc3dImage class containing the masked image

radmc3dPy.image.getPSF(nx=None, ny=None, psfType='gauss', pscale=None, fwhm=None, pa=None, tdiam_prim=8.2, tdiam_sec=0.94, wav=None)

Calculates a two dimensional Gaussian PSF.

Parameters:

nx : int

Image size in the first dimension

ny : int

Image size in the second dimension

psfType : {‘gauss’, ‘airy’}

Shape of the PSF. If psfType=’gauss’, fwhm and pa should also be given. If psfType=’airy’, the tdiam_prim, tdiam_sec and wav parameters should also be specified.

pscale : list

Pixelscale of the image, if set fwhm should be in the same unit, if not set unit of fwhm is pixels

fwhm : list, optional

Full width at half maximum of the psf along the two axis (should be set only if psfType=’gauss’)

pa : float, optional

Position angle of the gaussian if the gaussian is not symmetric (should be set only if psfType=’gauss’)

tdiam_prim : float, optional

Diameter of the primary aperture of the telescope in meter. (should be set only if psfType=’airy’)

tdiam_sec : float, optional

Diameter of the secondary mirror (central obscuration), if there is any, in meter. If no secondary mirror/obscuration is present, this parameter should be set to zero. (should be set only if psfType=’airy’)

wav : float, optional

Wavelength of observation in micrometer (should be set only if psfType=’airy’)

Returns:

Returns a dictionary with the following keys:

  • psf : ndarray
    The two dimensional psf
  • x : ndarray
    The x-axis of the psf
  • y : ndarray
    The y-axis of the psf
radmc3dPy.image.makeImage(npix=None, incl=None, wav=None, sizeau=None, phi=None, posang=None, pointau=None, fluxcons=True, nostar=False, noscat=False, widthkms=None, linenlam=None, vkms=None, iline=None, lambdarange=None, nlam=None, stokes=False, binary=False)

Calculates a rectangular image with RADMC-3D

Parameters:

npix : int

Number of pixels on the rectangular images

sizeau : float

Diameter of the image in au

incl : float

Inclination angle of the source

wav : float

Wavelength of the image in micron

phi : float, optional

Azimuthal rotation angle of the source in the model space

posang : float, optional

Position angle of the source in the image plane

pointau : Float, optional

Three elements list of the cartesian coordinates of the image center

widthkms : float, optional

Width of the frequency axis of the channel maps

linenlam : int, optional

Number of wavelengths to calculate images at

vkms : float, optional

A single velocity value at which a channel map should be calculated

iline : int, optional

Line transition index

lambdarange : list, optional

Two element list with the wavelenght boundaries between which multiwavelength images should be calculated

nlam : int, optional

Number of wavelengths to be calculated in lambdarange

fluxcons : bool, optional

This should not even be a keyword argument, it ensures flux conservation (adaptive subpixeling) in the rectangular images

nostar : bool, optional

If True the calculated images will not contain stellar emission

noscat : bool, optional

If True, scattered emission will be neglected in the source function, however,

extinction will contain scattering if kappa_scat is not zero.

stokes : bool, optional

If True, images in all four stokes parameters (IQUV) will be calculated, if False only the intensity will be calculated

binary : bool, optional

If True the output image will be written in a C-style binary format, if False the image format will be ASCII

radmc3dPy.image.plotImage(image=None, arcsec=False, au=False, log=False, dpc=None, maxlog=None, saturate=None, bunit='norm', ifreq=0, cmask_rad=None, interpolation='nearest', cmap=<matplotlib.colors.LinearSegmentedColormap object>, stokes='I', fig=None, ax=None, projection='polar', deg=True, rmax=None, rlog=True, **kwargs)

Plots a radmc3d image.

Parameters:

image : radmc3dImage

A radmc3dImage class returned by readimage

arcsec : bool

If True image axis will have the unit arcsec (NOTE: dpc keyword should also be set!)

au : bool

If True image axis will have the unit AU

log : bool

If True image scale will be logarithmic, otherwise linear

dpc : float

Distance to the source in parsec (This keywords should be set if arcsec=True, or bunit!=’norm’)

maxlog : float

Logarithm of the lowest pixel value to be plotted, lower pixel values will be clippde

saturate : float

Highest pixel values to be plotted in terms of the peak value, higher pixel values will be clipped

bunit : {‘norm’, ‘inu’, ‘snu’, ‘jy/beam’, ‘jy/pixel’}

Unit of the image, (‘norm’ - Inu/max(Inu), ‘inu’ - Inu, ‘snu’ - Jy/pixel, ‘jy/pixel’ - Jy/pixel, ‘jy/beam’ - Jy/beam), default is ‘norm’. The ‘snu’ keyword value is kept for backward compatibility as it is fully equivalent with the ‘jy/pixel’ keyword value.

ifreq : int

If the image file/array consists of multiple frequencies/wavelengths ifreq denotes the index of the frequency/wavelength in the image array to be plotted

cmask_rad : float

Simulates coronographyic mask : sets the image values to zero within this radius of the image center The unit is the same as the image axis (au, arcsec, cm) NOTE: this works only on the plot, the image array is not changed (for that used the cmask() function)

interpolation : str

interpolation keyword for imshow (e.g. ‘nearest’, ‘bilinear’, ‘bicubic’)

cmap : matplotlib color map

stokes : {‘I’, ‘Q’, ‘U’, ‘V’, ‘PI’, ‘P’}

What to plot for full stokes images, Stokes I/Q/U/V, PI - polarised intensity (PI = sqrt(Q^2 + U^2 + V^2)) P - polarisation fraction (i.e. sqrt(Q^2 + U^2 + V^2) / I) PIL - polarised intensity (PI = sqrt(Q^2 + U^2)) PL - fraction of linear polarisation (i.e. sqrt(Q^2 + U^2) / I)

fig : matplotlig.figure.Figure, optional

An instance of a matplotlib Figure. If not provided a new Figure will be generated. If provided plotImage will add a single Axes to the Figure, using add_subplots() with the appropriate projection. If the desired plot is to be made for a multi-panel plot, the appropriate Axes instance can be passed to the ax keyword. This keyword is only used for circular images.

ax : matplotlib.axes.Axes, optional

An instance of a matplotlib Axes to draw the plot on. Note, that the projection of the axes should be the same as the projection keyword passed to plotImage. This keyword is only used for circular images.

projection : {‘polar’, ‘cartesian’}

Projection of the plot. For cartesian plots a rectangular plot will be drawn, with the horizontal axis being the azimuth angle, and the vertical axis the radius.

deg : bool

If True the unit of the azimuthal coordinates will degree, if False it will be radian. Used only for circular images and for cartesian projection.

rmax : float

Maximum value of the radial coordinate for polar projection. Used only for circular images.

rlog : bool

If True the radial coordiante axis will be set to logarithmic for cartesian projection. Used only for circular images.

radmc3dPy.image.plotPolDir(image=None, arcsec=False, au=False, dpc=None, ifreq=0, cmask_rad=None, color='w', nx=20, ny=20)

Function to plot the polarisation direction for full stokes images

Parameters:

image : radmc3dImage

A radmc3dImage class returned by readimage

arcsec : bool

If True image axis will have the unit arcsec (NOTE: dpc keyword should also be set!)

au : bool

If True image axis will have the unit AU

dpc : float

Distance to the source in parsec (This keywords should be set if arcsec=True, or bunit!=’norm’)

ifreq : int

If the image file/array consists of multiple frequencies/wavelengths ifreq denotes the index of the frequency/wavelength in the image array to be plotted

cmask_rad : float

Simulates coronographyic mask : sets the image values to zero within this radius of the image center The unit is the same as the image axis (au, arcsec, cm) NOTE: this works only on the plot, the image array is not changed (for that used the cmask() function)

color : str

Color for the polarisation direction plot

nx : int

Number of grid points along the horizontal axis at which the direction should be displayed

ny : int

Number of grid points along the vertical axis at which the direction should be displayed

class radmc3dPy.image.radmc3dCircimage

Bases: object

RADMC-3D circular image class

Attributes

image (ndarray) The image as calculated by radmc3d (the values are intensities in erg/s/cm^2/Hz/ster)
rc (ndarray) Radial cell center coordinate of the image [cm]
ri (ndarray) Radial cell interface coordinate of the image [cm]
phic (ndarray) Azimuthal cell center coordinate of the image [rad]
phii (ndarray) Azimuthal cell interface coordinate of the image [rad]
nr (int) Number of pixels in the radial direction
nphi (int) Number of pixels in the azimuthal direction
nfreq (int) Number of frequencies in the image cube
freq (ndarray) Frequency grid in the image cube
nwav (int) Number of wavelengths in the image cube (same as nfreq)
wav (ndarray) Wavelength grid in the image cube
filename (str) Name of the file the image data was read from
stokes (bool) If True the image data contain the full stokes vector (I,Q,U,V)

Methods

getPixelSize() Calculates the pixel size
readImage([filename, old]) Reads a circular image
getPixelSize()

Calculates the pixel size

Returns:The pixel size in cm^2
readImage(filename='circimage.out', old=False)

Reads a circular image

Parameters:

filename : str

Name of the file to be read.

old : bool

If True the image format of the old 2D code (radmc) will be used. If False (default) the RADMC-3D format is used.

class radmc3dPy.image.radmc3dImage

Bases: object

RADMC-3D image class

Attributes

image (ndarray) The image as calculated by radmc3d (the values are intensities in erg/s/cm^2/Hz/ster)
imageJyppix (ndarray) The image with pixel units of Jy/pixel
x (ndarray) x coordinate of the image [cm]
y (ndarray) y coordinate of the image [cm]
nx (int) Number of pixels in the horizontal direction
ny (int) Number of pixels in the vertical direction
sizepix_x (float) Pixel size in the horizontal direction [cm]
sizepix_y (float) Pixel size in the vertical direction [cm]
nfreq (int) Number of frequencies in the image cube
freq (ndarray) Frequency grid in the image cube
nwav (int) Number of wavelengths in the image cube (same as nfreq)
wav (ndarray) Wavelength grid in the image cube

Methods

getClosurePhase([bl, pa, dpc]) Calculates clusure phases for a given model image for any arbitrary baseline triplet.
getMomentMap([moment, nu0, wav0]) Calculates moment maps.
getVisibility([bl, pa, dpc]) Calculates visibilities for a given set of projected baselines and position angles with Discrete Fourier Transform.
imConv([dpc, psfType, fwhm, pa, tdiam_prim, ...]) Convolves a RADMC-3D image with a two dimensional Gaussian psf.
plotMomentMap([moment, nu0, wav0, dpc, au, ...]) Plots moment maps
readImage([fname, binary, old]) Reads an image calculated by RADMC-3D
writeFits([fname, dpc, coord, bandwidthmhz, ...]) Writes out a RADMC-3D image data in fits format.
getClosurePhase(bl=None, pa=None, dpc=None)

Calculates clusure phases for a given model image for any arbitrary baseline triplet.

Parameters:

bl : list/ndarray

A list or ndrray containing the length of projected baselines in meter.

pa : list/ndarray

A list or Numpy array containing the position angles of projected baselines in degree.

dpc : distance of the source in parsec

Returns:

Returns a dictionary with the following keys:

  • bl : projected baseline in meter
  • pa : position angle of the projected baseline in degree
  • nbl : number of baselines
  • u : spatial frequency along the x axis of the image
  • v : spatial frequency along the v axis of the image
  • vis : complex visibility at points (u,v)
  • amp : correlation amplitude
  • phase : Fourier phase
  • cp : closure phase
  • wav : wavelength
  • nwav : number of wavelengths

Notes

bl and pa should either be an array with dimension [N,3] or if they are lists each element of the list should be a list of length 3, since closure phases are calculated only for closed triangles

getMomentMap(moment=0, nu0=None, wav0=None)

Calculates moment maps.

Parameters:

moment : int

Moment of the channel maps to be calculated

nu0 : float

Rest frequency of the line in Hz

wav0 : float

Rest wavelength of the line in micron

Returns:

Ndarray with the same dimension as the individual channel maps

getVisibility(bl=None, pa=None, dpc=None)

Calculates visibilities for a given set of projected baselines and position angles with Discrete Fourier Transform.

Parameters:

bl : list/ndarray

A list or ndrray containing the length of projected baselines in meter.

pa : list/ndarray

A list or Numpy array containing the position angles of projected baselines in degree.

dpc : distance of the source in parsec

Returns:

Returns a dictionary with the following keys:

  • bl : projected baseline in meter
  • pa : position angle of the projected baseline in degree
  • nbl : number of baselines
  • u : spatial frequency along the x axis of the image
  • v : spatial frequency along the v axis of the image
  • vis : complex visibility at points (u,v)
  • amp : correlation amplitude
  • phase : Fourier phase
  • wav : wavelength
  • nwav : number of wavelengths
imConv(dpc=1.0, psfType='gauss', fwhm=None, pa=None, tdiam_prim=8.2, tdiam_sec=0.94)

Convolves a RADMC-3D image with a two dimensional Gaussian psf. The output images will have the same brightness units as the input images.

Parameters:

dpc : float

Distance of the source in pc.

psfType : {‘gauss’, ‘airy’}

Shape of the PSF. If psfType=’gauss’, fwhm and pa should also be given. If psfType=’airy’, the tdiam_prim, tdiam_sec and wav parameters should also be specified.

fwhm : list, optional

A list of two numbers; the FWHM of the two dimensional psf along the two principal axes. The unit is assumed to be arcsec. (should only be set if psfType=’gauss’)

pa : float, optional

Position angle of the psf ellipse (counts from North counterclockwise, should only be set if psfType=’gauss’)

tdiam_prim : float, optional

Diameter of the primary aperture of the telescope in meter. (should be set only if psfType=’airy’)

tdiam_sec : float, optional

Diameter of the secondary mirror (central obscuration), if there is any, in meter. If no secondary mirror/obscuration is present, this parameter should be set to zero. (should only be set if psfType=’airy’)

Returns:

Returns a radmc3dImage

plotMomentMap(moment=0, nu0=None, wav0=None, dpc=1.0, au=False, arcsec=False, cmap=None, vclip=None)

Plots moment maps

Parameters:

moment : int

Moment of the channel maps to be calculated

nu0 : float

Rest frequency of the line in Hz

wav0 : float

Rest wavelength of the line in micron

dpc : float

Distance of the source in pc

au : bool

If True displays the image with AU as the spatial axis unit

arcsec : bool

If True displays the image with arcsec as the spatial axis unit (dpc should also be set!)

cmap : matplotlib colormap

Color map to be used to display the moment map

vclip : list/ndarray

Two element list / Numpy array containin the lower and upper limits for the values in the moment

map to be displayed

readImage(fname=None, binary=False, old=False)

Reads an image calculated by RADMC-3D

Parameters:

fname : str, optional

File name of the radmc3d output image (if omitted ‘image.out’ is used)

old : bool

If set to True it reads old radmc-2d style image

binary : bool, optional

False - the image format is formatted ASCII if True - C-compliant binary (omitted if old=True)

writeFits(fname='', dpc=1.0, coord='03h10m05s -10d05m30s', bandwidthmhz=2000.0, casa=False, nu0=0.0, stokes='I', fitsheadkeys=[], ifreq=None)

Writes out a RADMC-3D image data in fits format.

Parameters:

fname : str

File name of the radmc3d output image (if omitted ‘image.fits’ is used)

dpc : float

Distance of the source in pc

coord : str

Image center coordinates

bandwidthmhz : float

Bandwidth of the image in MHz (equivalent of the CDELT keyword in the fits header)

casa : bool

If set to True a CASA compatible four dimensional image cube will be written

nu0 : float

Rest frequency of the line (for channel maps)

stokes : {‘I’, ‘Q’, ‘U’, ‘V’, ‘PI’}

Stokes parameter to be written if the image contains Stokes IQUV (possible choices: ‘I’, ‘Q’, ‘U’, ‘V’, ‘PI’ -Latter being the polarized intensity)

fitsheadkeys : dictionary

Dictionary containing all (extra) keywords to be added to the fits header. If the keyword is already in the fits header (e.g. CDELT1) it will be updated/changed to the value in fitsheadkeys, if the keyword is not present the keyword is added to the fits header.

ifreq : int

Frequency index of the image array to write. If set only this frequency of a multi-frequency array will be written to file.

radmc3dPy.image.readImage(fname=None, binary=False, old=False)
Reads an image calculated by RADMC-3D.
This function is an interface to radmc3dImage.readImage().
Parameters:

fname : str, optional

File name of the radmc3d output image (if omitted ‘image.out’ is used)

old : bool

If set to True it reads old radmc-2d style image

binary : bool, optional

False - the image format is formatted ASCII if True - C-compliant binary (omitted if old=True)

radmc3dPy.image.readcircimage(filename='circimage.out', old=False)

A convenience function to read circular images

Parameters:

filename : str

Name of the file to be read.

old : bool

If True the image format of the old 2D code (radmc) will be used. If False (default) the RADMC-3D format is used.

radmc3dPy.miescat module

This module contains functions for Mie-scattering and to write dust opacity files

radmc3dPy.miescat.bhmie(x=None, refrel=None, theta=None)

The famous Bohren and Huffman Mie scattering code. This version was ported to Python from the f77 code from Bruce Draine, which can be downloaded from:

The code originates from the book by Bohren & Huffman (1983) on “Absorption and Scattering of Light by Small Particles”. This python version was created by Cornelis Dullemond, February 2017.

Parameters:

x : ndarray

Size parameter (2*pi*radius_grain/lambda)

refrel : complex

Complex index of refraction (example: 1.5 + 0.01*1j)

theta : ndarray

Scattering angles between 0 and 180 deg

Returns:

A list with the following elements:

  • [0] S1 : ndarray
    Complex phase function S1 (E perpendicular to scattering plane) as a function of theta
  • [1] S2 : ndarray
    Complex phase function S1 (E parallel to scattering plane) as a function of theta
  • [2] Qext : ndarray
    Efficiency factor for extinction (C_ext/pi*a**2)
  • [3] Qsca : ndarray
    Efficiency factor for scatterin(C_sca/pi*a**2)
  • [4] Qback : ndarray
    Backscattering efficiency ((dC_sca/domega)/pi*a**2 ) Note, this is (1/4*pi) smaller than the “radar backscattering efficiency” - see Bohren & Huffman 1983 pp. 120-123]
  • [5] gsca : <cos(theta)> for scattering
radmc3dPy.miescat.compute_opac_mie(fname='', matdens=None, agraincm=None, lamcm=None, theta=None, logawidth=None, wfact=3.0, na=20, chopforward=0.0, errtol=0.01, verbose=False, extrapolate=False)

Compute dust opacity with Mie theory based on the optical constants in the optconst_file. Optionally also the scattering phase function in terms of the Mueller matrix elements can be computed. To smear out the resonances that appear due to the perfect sphere shape, you can optionally smear out the grain size distribution a bit with setting the width of a Gaussian grain size distribution.

Parameters:

fname : str

File name of the optical constants file. This file should contain three columns: first the wavelength in micron, then the n-coefficient and then the k-coefficient. See Jena optical constants database: http://www.astro.uni-jena.de/Laboratory/Database/databases.html

matdens : float

Material density in g/cm^3

agraincm : float

Grain radius in cm

lamcm : ndarray

Wavelength grid in cm

theta : ndarray, optional

Angular grid (a numpy array) between 0 and 180 which are the scattering angle sampling points at which the scattering phase function is computed.

logawidth : float, optional

If set, the size agrain will instead be a sample of sizes around agrain. This helps to smooth out the strong wiggles in the phase function and opacity of spheres at an exact size. Since in Nature it rarely happens that grains all have exactly the same size, this is quite natural. The value of logawidth sets the width of the Gauss in ln(agrain), so for logawidth<<1 this give a real width of logawidth*agraincm.

wfact : float

Grid width of na sampling points in units of logawidth. The Gauss distribution of grain sizes is cut off at agrain * exp(wfact*logawidth) and agrain * exp(-wfact*logawidth). Default = 3

na : int

Number of size sampling points (if logawidth set, default=20)

chopforward : float

If >0 this gives the angle (in degrees from forward) within which the scattering phase function should be kept constant, essentially removing the strongly peaked forward scattering. This is useful for large grains (large ratio 2*pi*agraincm/lamcm) where the forward scattering peak is extremely strong, yet extremely narrow. If we are not interested in very forward-peaked scattering (e.g. only relevant when modeling e.g. the halo around the moon on a cold winter night), this will remove this component and allow a lower angular grid resolution for the theta grid.

errtol : float

Tolerance of the relative difference between kscat and the integral over the zscat Z11 element over angle. If this tolerance is exceeded, a warning is given.

verbose : bool

If set to True, the code will give some feedback so that one knows what it is doing if it becomes slow.

extrapolate : bool

If set to True, then if the wavelength grid lamcm goes out of the range of the wavelength grid of the optical constants file, then it will make a suitable extrapolation: keeping the optical constants constant for lamcm < minimum, and extrapolating log-log for lamcm > maximum.

Returns:

A dictionary with the following keys:

  • kabs : ndarray

    Absorption opacity kappa_abs_nu (a numpy array) in units of cm^2/gram

  • ksca : ndarray

    Scattering opacity kappa_abs_nu (a numpy array) in units of cm^2/gram

  • gsca : ndarray

    The <cos(theta)> g-factor of scattering

  • theta : ndarray (optional, only if theta is given at input)

    The theta grid itself (just a copy of what was given)

  • zscat : ndarray (optional, only if theta is given at input)

    The components of the scattering Mueller matrix Z_ij for each wavelength and each scattering angel. The normalization of Z is such that kscat can be reproduced (as can be checked) by the integral: 2*pi*int_{-1}^{+1}Z11(mu)dmu=kappa_scat. For symmetry reasons only 6 elements of the Z matrix are returned: Z11, Z12, Z22, Z33, Z34, Z44. Note that Z21 = Z12 and Z43 = -Z34. The scattering matrix is normalized such that if a plane wave with Stokes flux

    Fin = (Fin_I,Fin_Q,Fin_U,Fin_V)

    hits a dust grain (which has mass mgrain), then the scattered flux

    Fout = (Fout_I,Fout_Q,Fout_U,Fout_V)

    at distance r from the grain at angle theta is given by

    Fout(theta) = (mgrain/r^2) * Zscat . Fin

    where . is the matrix-vector multiplication. Note that the Stokes components must be such that the horizontal axis in the “image” is pointing in the scattering plane. This means that radiation with Fin_Q < 0 is scattered well, because it is vertically polarized (along the scattering angle axis), while radiation with Fin_Q > 0 is scatterd less well because it is horizontally polarized (along the scattering plane).

  • kscat_from_z11 : ndarray (optional, only if theta is given at input)

    The kscat computed from the (above mentioned) integral of Z11 over all angles. This should be nearly identical to kscat if the angular grid is sufficiently fine. If there are strong differences, this is an indication that the angular gridding (the theta grid) is not fine enough. But you should have then automatically gotten a warning message as well (see errtol).

  • wavmic : ndarray (optional, only if extrapolate is set to True)

    The original wavelength grid from the optical constants file, with possibly an added extrapolated

  • ncoef : ndarray (optional, only if extrapolate is set to True)

    The optical constant n at that grid

  • kcoef : ndarray (optional, only if extrapolate is set to True)

    The optical constant k at that grid

  • agr : ndarray (optional, only if logawidth is not None)

    Grain sizes

  • wgt : ndarray (optional, only if logawidth is not None)

    The averaging weights of these grain (not the masses!) The sum of wgt.sum() must be 1.

  • zscat_nochop : ndarray (optional, only if chopforward > 0)

    The zscat before the forward scattering was chopped off

  • kscat_nochop : ndarray (optional, only if chopforward > 0)

    The kscat originally from the bhmie code

radmc3dPy.miescat.write_radmc3d_kappa_file(package=None, name=None)
The RADMC-3D radiative transfer package
http://www.ita.uni-heidelberg.de/~dullemond/software/radmc-3d/

can perform dust continuum radiative transfer for diagnostic purposes. It is designed for astronomical applications. The code needs the opacities in a particular form. This subroutine writes the opacities out in that form. It will write it to the file dustkappa_<name>.inp. This is the simpler version of the opacity files, containing only kabs, kscat, gscat as a function of wavelength.

radmc3dPy.miescat.write_radmc3d_scatmat_file(package=None, name=None, comment=False)
The RADMC-3D radiative transfer package
http://www.ita.uni-heidelberg.de/~dullemond/software/radmc-3d/

can perform dust continuum radiative transfer for diagnostic purposes. It is designed for astronomical applications. The code needs the opacities in a particular form. This subroutine writes the opacities out in that form. It will write it to the file dustkapscatmat_<name>.inp.

radmc3dPy.molecule module

This module contains classes to handle molecular data

class radmc3dPy.molecule.radmc3dMolecule

Bases: object

RADMC-3D molecule class Based on the Leiden LAMDA database, but is in principle generic

NOTE: For now only the levels and lines are included, not the
collision rates.

Attributes

name (str) The name as listed in the molecule file
molweight (float) Molecular weight in units of proton mass
nlev (int) Nr of energy levels
nlin (int) Nr of lines
energycminv (float) Energy[ilev] of level ilev in 1/cm
energy (float) Energy[ilev] of level ilev in erg
wgt (float) Statistical weight[ilev] of level ilev
jrot (float) Quantum rotational J[ilev] of level ilev
iup (int) ilev of upper level of line ilin (starting with 0)
ilow (int) ilev of lower level of line ilin (starting with 0)
aud (float) Einstein A up low of line ilin in 1/second
freq (float) Frequency of line ilin in Hz
lam (float) Wavelength of line ilin in micron
temp (ndarray) Temperature grid of the partition function
pfunc (ndarray) Partition function

Methods

getPartitionFunction([temp, tmin, tmax, ...]) Calculates the partition function at a grid of temperatures
read([mol, fname]) Read the molecule_<mol>.inp file
getPartitionFunction(temp=None, tmin=None, tmax=None, ntemp=None, tlog=True)

Calculates the partition function at a grid of temperatures

Parameters:

temp : list,ndarray

Temperature(s) in Kelvin to calculate the partition function at

tmin : float, optional

Minimum temperature in the grid (if temp is None)

tmax : float, optional

Maximum temperature in the grid (if temp is None)

ntemp : int, optional

Number of temperature in the grid (if temp is None)

tlog : bool

If True the generated temperature grid will be logarithmic. If False the grid will be linear.

Returns:

The temperature grid and partition function will be put in the temp and pfunc data attributes of the class

read(mol=None, fname=None)

Read the molecule_<mol>.inp file

The file format is the format of the Leiden LAMDA molecular database

Parameters:

mol : str

molecule name (e.g. ‘co’) if the file name is in the form of ‘molecule_<mol>.inp’

fname : str

full file name

radmc3dPy.natconst module

This module contains natural constants in CGS units

Translated from RADMC’s IDL function problem_natconst.pro

List of natural constants:

Name Description
gg Gravitational constant
mp Mass of proton [g]
me Mass of electron [g]
kk Bolzmann’s constant [erg/K]
hh Planck’s constant [erg.s]
ee Unit charge
cc Light speed [cm/s]
st Thmpson cross-section [cm^2]
ss Stefan-Boltzmann const [erg/cm^2/K^4/s]
aa 4 * ss / cc
muh2 Mean molecular weight (H2 + He + metals)
ev Electronvolt [erg]
kev Kilo electronvolt [erg]
micr Micron [cm]
km Kilometer [cm]
angs Angstrom [cm]
ls Solar luminosity [erg/s]
rs Solar radius [cm]
ms Solar mass [g]
ts Solar effective temperature [K]
au Astronomical unit [cm]
pc Parsec [cm]
mea Mass of Earth [g]
rea Equatorila radius of Earth [cm]
mmo Mass of Moon [g]
rmo Radius of Moon [cm]
dmo Distance earth-moon (center-to-center) [cm]
mju Mass of Jupiter
rju Equatorial radius of Jupiter [cm]
dju Distance Jupiter-Sun [cm]
year Year [s]
hour Hour [s]
day Day [s]

radmc3dPy.octree module

This module contains a class for handling Octree mesh

class radmc3dPy.octree.radmc3dOctree

Bases: object

Octree-like object with switchable resolution in each dimension

Attributes

xi (ndarray) Base grid cell interface grid in the first dimension
yi (ndarray) Base grid cell interface grid in the second dimension
zi (ndarray) Base grid cell interface grid in the third dimension
xc (ndarray) Base grid cell center grid in the first dimension
yc (ndarray) Base grid cell center grid in the second dimension
zc (ndarray) Base grid cell center grid in the third dimension
x (ndarray) Tree cell center array in the first dimension
y (ndarray) Tree cell center array in the second dimension
z (ndarray) Tree cell center array in the third dimension
dx (ndarray) Tree cell halfwidth array in the first dimension
dy (ndarray) Tree cell halfwidth array in the second dimension
dz (ndarray) Tree cell halfwidth array in the third dimension
leafID (ndarray) Leaf index array, mapping between a full tree and an array containing only the leaves
isLeaf (ndarray) Boolean array for the cell type (True - leaf, False - branch)
level (ndarray) Level array (base grid is level 0)
parentID (ndarray) Array containing the index of the parent cell (currently unused, only needed if we go up in the tree)
childID (list) List of children indices. Each list element is an ndarray with nChild elements containing the child indices
act_dim (list) A three element array to indicate which dimension is active, i.e. which dimensions are the cells resolved (0 - inactive, 1 - active)
nCell (int) Nr of cells (both branch and leaf) in the tree
nxRoot (int) Nr of cells in the base grid in the first dimension
nyRoot (int) Nr of cells in the base grid in the second dimension
nzRoot (int) Nr of cells in the base grid in the third dimension
nLeaf (int) Nr of leaf cells (i.e. true, unresolved grid cells)
nBranch (int) Nr of branches (i.e. resolved cells)
nChild (int) Nr of children (i.e. 8, 4, or 2 for 3, 2, 1 active dimensions, respectively)
levelMax (int) Highest actual level in the tree (Base grid has a level of 0 and the level increases)
levelMaxLimit (int) Highest allowed level in the tree (only used in tree building)
crd_sys ({‘car’, ‘sph’}) Coordinate system type cartesian or spherical

Methods

convArrLeaf2Tree([var]) Converts a leaf array to full tree size.
convArrTree2Leaf([var]) Converts a data array to leaf size.
generateLeafID() Function to generate the cell index mapping from arrays containing the full tree and those containing
getCellVolume([fullTree]) Calculates the grid cell volume
getContainerLeafID([crd]) Finds the tree index of a leaf that contains a given coordinate
makeSpatialGrid([ppar, levelMaxLimit, ...]) Function to create an octree-like AMR grid
makeWavelengthGrid([wbound, nw, ppar]) Creates the wavelength/frequency grid.
putNode([crd, cellsize, level, parentID, cellID]) Function to put the data of a single node into the tree.
readGrid() Reads the spatial and wavelength grids from files
readSpatialGrid([fname]) Reads the spatial grid from amr_grid.inp
readWavelengthGrid([fname]) Function to read the wavelength/frequency grid
resolveNodes([rsIDs]) Resolve multiple nodes simultaneously and add the children of the resolved node to the tree arrays extending
selfCheck() Performs a self-check of the tree allocation and report it to the screen
setModel([model]) Sets the model to be used for tree building
writeSpatialGrid([fname]) Writes the wavelength grid to a file (e.g.
writeWavelengthGrid([fname, old]) Wriites the wavelength grid to a file (e.g.
convArrLeaf2Tree(var=None)

Converts a leaf array to full tree size.

Parameters:

var : ndarray

A one or two dimensional ndarray with the first dimension is the size of the full tree

Returns:

A one or two dimensional ndarray with size of of the full tree in the first dimension

convArrTree2Leaf(var=None)

Converts a data array to leaf size. The input is a scalar or vector variable defined at all nodes and the returned variable will only represent values at leaf nodes thereby reduced in length compared to the input.

Parameters:

var : ndarray

A one or two dimensional ndarray with the first dimension is the size of the full tree

Returns:

A one or two dimensional ndarray with size of nLeaf in the first dimension

generateLeafID()

Function to generate the cell index mapping from arrays containing the full tree and those containing only the leaves

getCellVolume(fullTree=False)

Calculates the grid cell volume

Parameters:

fullTree : bool, optional

If True the cell volumes of the full tree (including both branches and leaves) will be calculated, while if set to False (default) the volume of only the leaf cells will be calculated

Returns:

An linear ndarray containing the cell volumes

getContainerLeafID(crd=())

Finds the tree index of a leaf that contains a given coordinate

Parameters:

crd : tuple

List/tuple/ndarray containing the tree dimensional coordinates of the point

makeSpatialGrid(ppar=None, levelMaxLimit=None, dfunc=None, model='', **kwargs)

Function to create an octree-like AMR grid

Parameters:

ppar : dictionary

Dictionary containing all input parameters of the model (from the problem_params.inp file)

model : str

Name of the model to be used in the tree building

dfunc : function

A user defined function that decides whether an AMR grid cell should be refined

levelMaxLimit : int, optional

Highest allowable level of the tree. This keyword is optional. If not specified at input as a separate keyword, levelMaxLimit should be present in the problem_params.inp file.

makeWavelengthGrid(wbound=None, nw=None, ppar=None)

Creates the wavelength/frequency grid.

Parameters:

wbound : list

Contains the wavelength boundaries of the wavelength grid (should contain at least two elements)

nw : list

Contains len(wbound)-1 elements containing the number of wavelengths between the bounds set by wbound

ppar : dictionary, optional

Contains all input parameters with the parameter names as keys

putNode(crd=(), cellsize=(), level=None, parentID=-1, cellID=None)

Function to put the data of a single node into the tree. This funcion assumes that all the arrays have already been allocated for the tree so input cell indices must refer to already existing array elements.

Parameters:

crd : tuple

Cell center coordinates of the node

cellsize : tuple

Full size of the cell in each dimension

level : int

Level of the cell in the tree

parentID : int

Tree index of the parent cell

cellID : int

Tree index of the cell to be added

readGrid()

Reads the spatial and wavelength grids from files

readSpatialGrid(fname='')

Reads the spatial grid from amr_grid.inp

Parameters:

fname : str, optional

File name from which the spatial grid should be read. If omitted ‘amr_grid.inp’ will be used.

readWavelengthGrid(fname='wavelength_micron.inp')

Function to read the wavelength/frequency grid

Parameters:

fname : str, optional

Name of the file to read the wavelength grid from (if not specified wavelenth_micron.inp will be used)

resolveNodes(rsIDs=None)

Resolve multiple nodes simultaneously and add the children of the resolved node to the tree arrays extending the tree array

Parameters:

rsIDs : list

List/tuple/array of indices of the resolvable cell in the tree array

selfCheck()

Performs a self-check of the tree allocation and report it to the screen

setModel(model='')

Sets the model to be used for tree building

Parameters:

model : str

Name of the model

writeSpatialGrid(fname='')

Writes the wavelength grid to a file (e.g. amr_grid.inp).

Parameters:

fname : str, optional

File name into which the spatial grid should be written. If omitted ‘amr_grid.inp’ will be used.

writeWavelengthGrid(fname='', old=False)

Wriites the wavelength grid to a file (e.g. wavelength_micron.inp).

Parameters:

fname : str, optional

File name into which the wavelength grid should be written. If omitted ‘wavelength_micron.inp’ will be used

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

radmc3dPy.params module

This module contains classes for handling of parameters of the model setup and for file I/O

class radmc3dPy.params.radmc3dPar

Bases: object

Parameter class for a RADMC-3D model.

Attributes

ppar (dictionary) Contains parameter values with parameter names as keys
pdesc (dictionary) Contains parameter description (comments in the parameter file) with parameter names as keys
pblock (dictionary) Contains the block names in the parameter file and parameter names as values
pvalstr: dictionary Contains parameter values as strings with parameter names as keys

Methods

loadDefaults([model, ppar, reset]) Sets all parameters to default values.
printPar() Prints the parameters of the current model.
readPar([fname]) Reads a parameter file.
setPar([parlist]) Sets a parameter in the radmc3DPar class.
writeParfile([fname]) Writes a parameter file.
loadDefaults(model='', ppar=None, reset=True)

Sets all parameters to default values.

Parameters:

model : str

Model name whose paraemters should also be loaded

ppar : dictionary

Contains parameter values as string and parameter names as keys Default values will be re-set to the values in this dictionary

reset : bool

If True the all class attributes will be re-initialized before the default values would be loaded. I.e. it will remove all entries from the dictionary that does not conain default values either in this function or in the optional ppar keyword argument

printPar()

Prints the parameters of the current model.

readPar(fname='')

Reads a parameter file. The parameters in the files should follow python syntax

Parameters:

fname : str, optional

File name to be read (if omitted problem_params.inp is used)

Returns:

Returns a dictionary with the parameter names as keys

setPar(parlist=None)

Sets a parameter in the radmc3DPar class. If the paramter is already defined its value will be modified

Parameters:

parlist : list

If the parameter is already defined parlist should be a two element list 1) parameter name, 2) parameter expression/value as a string

If the parameter is not yet defined parlist should be a four element list 1) parameter name, 2) parameter expression/value as a string 3) Parameter description (= comment field in the parameter file)

writeParfile(fname='')

Writes a parameter file.

Parameters:

fname : str, optional

File name to be read (if omitted problem_params.inp is used)

radmc3dPy.radsources module

This module contains classes for radiation sources

class radmc3dPy.radsources.radmc3dRadSources(ppar=None, grid=None)

Bases: object

Class of the radiation sources. Currently discrete stars and continuous starlike source, the latter only in spherical coordinates.

Attributes

wav (ndarray) Wavelength for the stellar spectrum
freq (ndarray) Frequency for the stellar spectrum
nwav (int) Number of wavelenghts in the stellar spectrum
nfreq (int) Number of frequencies in the stellar spectrum
mstar (list) List of stellar masses
tstar (list) List of stellar effective temperatures
rstar (list) List of stellar radii
lstar (list) List of stellar luminosities
nstar (int) Number of stars
pstar (list) Each element of the list contains a three element list, the cartesian coordinates of the stars
fnustar (ndarray) Stellar spectrum (flux@1pc in erg/s/cm/cm/Hz)
csdens (ndarray) Stellar density for continuous starlike source
csntemplate (int) Number of stellar templates
cstemp (ndarray) Stellar template
cstemptype (int) Stellar template type 1 - Blackbody given by the effective temperature 2 - Frequency dependent spectrum
cststar (ndarray) Stellar effective temperature
csmstar (ndarray) Stellar mass
csrstar (ndarray) Stellar radius
tacc (ndarray) Effective temperature of a viscous accretion disk as a function of cylindrical radius
accrate (float) Accretion rate of the viscous accretion disk [g/s]
fnuaccdisk (ndarray) Spatially integrated frequency-dependent flux density of the accretion disk @ 1pc distance
tspot (float) Temperature of the hot spot / boundary layer on the stellar surface
starsurffrac (float) Fraction of the stellar surface covered by the hot spot / boundary layer
fnustpot (ndarray) Frequency-dependent flux density of the hot spot / boundary layer @ 1pc distance

Methods

findPeakStarspec() Calculates the peak wavelength of the stellar spectrum.
getAccdiskSpectra([ppar, grid, incl]) Calculates the emergent spectra of an optically thick accretion disk at face-on orientation (incl=0deg).
getAccdiskStellarDensity([grid]) Calculates the stellar density for continuous starlike sources for modeling a viscous accretion disk.
getAccdiskStellarTemplates([ppar, grid]) Calculates the stellar template for continuous starlike sources for modeling a viscous accretion disk.
getAccdiskTemperature([ppar, grid]) Calculates the effective temperature of a viscous accretion disk.
getSpotSpectrum([ppar, grid]) Calculates the spectrum of a hot spot / boundary layer on the stellar surface
getStarSpectrum([tstar, rstar, lstar, ...]) Calculates a blackbody stellar spectrum.
getTotalLuminosities([readInput]) Calcultes the frequency integrated luminosities of all radiation sources.
readStarsinp([fname]) Reads the data of discrete stellar sources from the stars.inp file.
readStellarsrcDensity([fname, binary]) Reads the stellar density of a continuous starlike source.
readStellarsrcTemplates([fname]) Reads the stellar template of a continuous starlike source.
writeStarsinp([ppar, wav, freq, old]) Writes the input file for discrete stellar sources (stars.inp).
writeStellarsrcDensity([fname, binary]) Writes the stellar density of a continuous starlike source.
writeStellarsrcTemplates([fname]) Writes the stellar template of a continuous starlike source.
findPeakStarspec()

Calculates the peak wavelength of the stellar spectrum.

Returns:

The peak wavelength of the stellar spectrum in nu*Fnu for all

stars as a list

getAccdiskSpectra(ppar=None, grid=None, incl=0.0)

Calculates the emergent spectra of an optically thick accretion disk at face-on orientation (incl=0deg).

Parameters:

ppar : dictionary

Dictionary containing all input parameters keys should include * mstar : stellar mass * rstar : stellar radius * accrate : accretion rate

NOTE, that for the calculation of the effective disk temperature only the first star is used if more than one values are given in mstar and rstar.

incl : float, optional

Inclination angle in degrees at which the spectrum be calculated (default - 0deg)

grid : radmc3dGrid, optional

An instance of a radmc3dGrid class containing the spatial and wavelength grid

getAccdiskStellarDensity(grid=None)

Calculates the stellar density for continuous starlike sources for modeling a viscous accretion disk.

Parameters:

grid : radmc3dGrid, optional

An instance of a radmc3dGrid class containing the spatial and wavelength grid

getAccdiskStellarTemplates(ppar=None, grid=None)

Calculates the stellar template for continuous starlike sources for modeling a viscous accretion disk.

Parameters:

ppar : dictionary

Dictionary containing all input parameters keys should include: * mstar : stellar mass * rstar : stellar radius * accrate : accretion rate

NOTE, that for the calculation of the effective disk temperature only the first star is used if more than one values are given in mstar and rstar.

grid : radmc3dGrid, optional

An instance of a radmc3dGrid class containing the spatial and wavelength grid

getAccdiskTemperature(ppar=None, grid=None)

Calculates the effective temperature of a viscous accretion disk.

Parameters:

ppar : dictionary

Dictionary containing all input parameters keys should include * mstar : stellar mass * rstar : stellar radius * accrate : accretion rate

NOTE, that for the calculation of the effective disk temperature only the first star is used if more than one values are given in mstar and rstar.

grid : radmc3dGrid, optional

An instance of a radmc3dGrid class containing the spatial and wavelength grid

getSpotSpectrum(ppar=None, grid=None)

Calculates the spectrum of a hot spot / boundary layer on the stellar surface

Parameters:

ppar : dictionary

Dictionary containing all input parameters keys should include * mstar : stellar mass * rstar : stellar radius * accrate : accretion rate

NOTE, that for the calculation of the effective disk temperature only the first star is used if more than one values are given in mstar and rstar.

grid : radmc3dGrid, optional

An instance of a radmc3dGrid class containing the spatial and wavelength grid

getStarSpectrum(tstar=None, rstar=None, lstar=None, mstar=None, ppar=None, grid=None)

Calculates a blackbody stellar spectrum.

Parameters:

tstar : list

Effective temperature of the stars in [K]

rstar : list

Radius of the stars in [cm]

lstar : list

Bolometric luminosity of the star [erg/s] (either rstar or lstar should be given)

mstar : list

Stellar mass in [g] (only required if an atmosphere model is used to calculate logg)

ppar : dictionary

Dictionary containing all input parameters

grid : radmc3dGrid, optional

An instance of a radmc3dGrid class containing the spatial and wavelength grid

getTotalLuminosities(readInput=True)

Calcultes the frequency integrated luminosities of all radiation sources.

Parameters:

readInput : bool, optional

If true the input files of the radiation sources are read and the the total luminosities are calculated from them. If readInput is set to False, the luminosities are calculated by semi-analytic spectra.

Returns:

Returns a dictionary with the following keys

  • lnu_star : Luminosity of the discrete stars
  • lnu_accdisk : Luminosity of the accretion disk
  • lnu_spot : Luminosity of the hot spot / boundary layer on the stellar surface
readStarsinp(fname='')

Reads the data of discrete stellar sources from the stars.inp file.

Parameters:

fname : str, optional

File name of the file that should be read (if omitted stars.inp will be used)

readStellarsrcDensity(fname=None, binary=False)

Reads the stellar density of a continuous starlike source.

Parameters:

fname : str, optional

Name of the file from which the stellar templates will be read. If omitted the default ‘stellarsrc_templates.inp’ will be used.

binary : bool, optional

If True the file should contain a C-style binary stream, if False the file should be written as formatted ASCII

readStellarsrcTemplates(fname='stellarsrc_templates.inp')

Reads the stellar template of a continuous starlike source.

Parameters:

fname : str, optional

Name of the file from which the stellar templates will be read. If omitted the default ‘stellarsrc_templates.inp’ will be used.

writeStarsinp(ppar=None, wav=None, freq=None, old=False)

Writes the input file for discrete stellar sources (stars.inp).

Parameters:

ppar : dictionary

Dictionary containing all parameters of the model (only mandatory if accretion is switched on)

wav : ndarray, optional

Wavelength grid for the stellar spectrum

freq : ndarray, optional

Frequency grid for the stellar spectrum (either freq or wav should be set)

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

writeStellarsrcDensity(fname='', binary=False)

Writes the stellar density of a continuous starlike source.

Parameters:

fname : str, optional

Name of the file into which the stellar templates will be written. If omitted the default ‘stellarsrc_templates.inp’ will be used.

binary : bool, optional

If True the output will be written in a C-style binary stream, if False the output will be formatted ASCII

writeStellarsrcTemplates(fname='stellarsrc_templates.inp')

Writes the stellar template of a continuous starlike source.

Parameters:

fname : str, optional

Name of the file into which the stellar templates will be written. If omitted the default ‘stellarsrc_templates.inp’ will be used.

radmc3dPy.reggrid module

This module contains a class for handling regular wavelength and spatial grids

class radmc3dPy.reggrid.radmc3dGrid

Bases: object

Class for spatial and frequency grid used by RADMC-3D.

Attributes

act_dim (ndarray) A three element vector the i-th element is 1 if the i-th dimension is active, otherwize the i-th element is zero
crd_sys ({‘sph’, ‘cyl’, ‘car’}) coordinate system of the spatial grid
nx (int) Number of grid points in the x (cartesian) / r (cylindrical) / r (spherical) dimension
ny (int) Number of grid points in the y (cartesian) / theta (cylindrical) / theta (spherical) dimension
nz (int) Number of grid points in the z (cartesian) / z (cylindrical) / phi (spherical) dimension
nxi (int) Number of cell interfaces in the x (cartesian) / r (cylindrical) / r (spherical) dimension
nyi (int) Number of cell interfaces in the y (cartesian) / theta (cylindrical) / theta (spherical) dimension
nzi (int) Number of cell interfaces in the z (cartesian) / z (cylindrical) / phi (spherical) dimension
nwav (int) Number of wavelengths in the wavelength grid
nfreq (int) Number of frequencies in the grid (equal to nwav)
x (ndarray) Cell centered x (cartesian) / r (cylindrical) / r (spherical) grid points
y (ndarray) Cell centered y (cartesian) / theta (cylindrical) / theta (spherical) grid points
z (ndarray) Cell centered z (cartesian) / z (cylindrical) / phi (spherical) grid points
xi (ndarray) Cell interfaces in the x (cartesian) / r (cylindrical) / r (spherical) dimension
yi (ndarray) Cell interfaces in the y (cartesian) / theta (cylindrical) / theta (spherical) dimension
zi (ndarray) Cell interfaces in the z (cartesian) / z (cylindrical) / phi (spherical) dimension
wav (ndarray) Wavelengh grid
freq (ndarray) Frequency grid

Methods

getCellVolume() Calculates the volume of grid cells.
makeSpatialGrid([crd_sys, xbound, ybound, ...]) Calculates the spatial grid.
makeWavelengthGrid([wbound, nw, ppar]) Creates the wavelength/frequency grid.
readGrid([old]) Reads the spatial (amr_grid.inp) and frequency grid (wavelength_micron.inp).
readSpatialGrid([fname, old]) Reads the spatial grid
readWavelengthGrid([fname, old]) Reads the wavelength grid
writeSpatialGrid([fname, old]) Writes the wavelength grid to a file (e.g.
writeWavelengthGrid([fname, old]) Wriites the wavelength grid to a file (e.g.
getCellVolume()

Calculates the volume of grid cells.

makeSpatialGrid(crd_sys=None, xbound=None, ybound=None, zbound=None, nxi=None, nyi=None, nzi=None, ppar=None)

Calculates the spatial grid.

Parameters:

crd_sys : {‘sph’,’car’}

Coordinate system of the spatial grid

xbound : list

(with at least two elements) of boundaries for the grid along the first dimension

ybound : list

(with at least two elements) of boundaries for the grid along the second dimension

zbound : list

(with at least two elements) of boundaries for the grid along the third dimension

nxi : int

Number of grid points along the first dimension. List with len(xbound)-1 elements with nxi[i] being the number of grid points between xbound[i] and xbound[i+1]

nyi : int

Same as nxi but for the second dimension

nzi : int

Same as nxi but for the third dimension

ppar : Dictionary containing all input parameters of the model (from the problem_params.inp file)

if ppar is set all keyword arguments that are not set will be taken from this dictionary

makeWavelengthGrid(wbound=None, nw=None, ppar=None)

Creates the wavelength/frequency grid.

Parameters:

wbound : list

Contains the wavelength boundaries of the wavelength grid (should contain at least two elements)

nw : list

Contains len(wbound)-1 elements containing the number of wavelengths between the bounds set by wbound

ppar : dictionary, optional

Contains all input parameters with the parameter names as keys

readGrid(old=False)

Reads the spatial (amr_grid.inp) and frequency grid (wavelength_micron.inp).

Parameters:

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

readSpatialGrid(fname='', old=False)

Reads the spatial grid

Parameters:

fname : str, optional

File name from which the spatial grid should be read. If omitted ‘amr_grid.inp’ will be used.

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

readWavelengthGrid(fname=None, old=False)

Reads the wavelength grid

Parameters:

fname : str, optional

File name from which the spatial grid should be read. If omitted ‘wavelength_micron.inp’ will be used.

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

writeSpatialGrid(fname='', old=False)

Writes the wavelength grid to a file (e.g. amr_grid.inp).

Parameters:

fname : str, optional

File name into which the spatial grid should be written. If omitted ‘amr_grid.inp’ will be used.

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

writeWavelengthGrid(fname='', old=False)

Wriites the wavelength grid to a file (e.g. wavelength_micron.inp).

Parameters:

fname : str, optional

File name into which the wavelength grid should be written. If omitted ‘wavelength_micron.inp’ will be used

old : bool, optional

If set to True the file format of the previous, 2D version of radmc will be used

radmc3dPy.setup module

This module contains functions to set up a RADMC-3D model for dust and/or line simulations. For help on the syntax or functionality of each function see the help of the individual functions

radmc3dPy.setup.problemSetupDust(model=None, binary=True, writeDustTemp=False, old=False, dfunc=None, dfpar=None, **kwargs)

Function to set up a dust model for RADMC-3D

Parameters:

model : str

Name of the model that should be used to create the density structure. The file should be in a directory from where it can directly be imported (i.e. the directory should be in the PYTHON_PATH environment variable or it should be in the current working directory) and the file name should be ‘model_xxx.py’, where xxx stands for the string that should be specified in this variable

binary : bool, optional

If True input files will be written in binary format, if False input files are written as formatted ascii text.

writeDustTemp : bool, optional

If True a separate dust_temperature.inp/dust_tempearture.binp file will be written under the condition that the model contains a function getDustTemperature()

old : bool, optional

If set to True the input files for the old 2D version of radmc will be created

dfunc : function, optional

Decision function for octree-like amr tree building. It should take linear arrays of cell centre coordinates (x,y,z) and cell half-widhts (dx,dy,dz) in all three dimensions, a radmc3d model, a dictionary with all parameters from problem_params.inp and an other keyword argument (**kwargs). It should return a boolean ndarray of the same length as the input coordinates containing True if the cell should be resolved and False if not. An example for the implementation of such decision function can be found in radmc3dPy.analyze module (radmc3dPy.analyze.gdensMinMax()).

dfpar : dictionary

Dicionary of keyword arguments to be passed on to dfunc. These parameters will not be written to problem_params.inp. Parameters can also be passed to dfunc via normal keyword arguments gathered in **kwargs, however all keyword arguments in **kwargs will be written to problem_params.inp

**kwargs : Any varible name in problem_params.inp can be used as a keyword argument.

At first all variables are read from problem_params.in to a dictionary called ppar. Then if there is any keyword argument set in the call of problem_setup_dust the ppar dictionary is searched for this key. If found the value belonging to that key in the ppar dictionary is changed to the value of the keyword argument. If no such key is found then the dictionary is simply extended by the keyword argument. Finally the problem_params.inp file is updated with the new parameter values.

Notes

Files written by problemSetupDust() for RADMC-3D

  • dustopac.inp : Dust opacity master file.
  • wavelength_micron.inp : Wavelength grid.
  • amr_grid.inp : Spatial grid.
  • stars.inp : Input radiation field (discrete stellar sources).
  • stellarsrc_density.inp : Input radiation field (continuous stellar sources).
  • stellarsrc_templates.inp : Input radiation field (continuous stellar sources).
  • dust_density.inp : Dust density distribution.
  • radmc3d.inp : Parameters for RADMC-3D (e.g. Nr of photons to be used, scattering type, etc).
radmc3dPy.setup.problemSetupGas(model=None, fullsetup=False, binary=True, writeGasTemp=False, dfunc=None, dfpar=None, **kwargs)

Function to set up a gas model for RADMC-3D

Parameters:

model : str

Name of the model that should be used to create the density structure the file should be in a directory from where it can directly be imported (i.e. the directory should be in the PYTHON_PATH environment variable, or it should be the current working directory) and the file name should be ‘MODELNAME.py’, where MODELNAME stands for the string that should be specified in this variable

fullsetup : bool, optional

If False only the files related to the gas simulation is written out (i.e. no grid, stellar parameter file and radmc3d master command file is written) assuming that these files have already been created for a previous continuum simulation. If True the spatial and wavelength grid as well as the input radiation field and the radmc3d master command file will be (over)written.

binary : bool, optional

If True input files will be written in binary format, if False input files are written as formatted ascii text.

writeGasTemp : bool, optional

If True a separate gas_temperature.inp/gas_tempearture.binp file will be written under the condition that the model contains a function get_gas_temperature()

dfunc : function, optional

Decision function for octree-like amr tree building. It should take linear arrays of cell centre coordinates (x,y,z) and cell half-widhts (dx,dy,dz) in all three dimensions, a radmc3d model, a dictionary with all parameters from problem_params.inp and an other keyword argument (**kwargs). It should return a boolean ndarray of the same length as the input coordinates containing True if the cell should be resolved and False if not. An example for the implementation of such decision function can be found in radmc3dPy.analyze module (radmc3dPy.analyze.gdensMinMax()).

dfpar : dictionary

Dicionary of keyword arguments to be passed on to dfunc. These parameters will not be written to problem_params.inp. Parameters can also be passed to dfunc via normal keyword arguments gathered in **kwargs, however all keyword arguments in **kwargs will be written to problem_params.inp

**kwargs : Any varible name in problem_params.inp can be used as a keyword argument.

At first all variables are read from problem_params.in to a dictionary called ppar. Then if there is any keyword argument set in the call of problem_setup_gas the ppar dictionary is searched for such key. If found the value belonging to that key in the ppar dictionary is changed to the value of the keyword argument. If no such key is found then the dictionary is simply extended by the keyword argument. Finally the problem_params.inp file is updated with the new parameter values. Any additional keyword argument for the octree AMR mesh generation should also be passed here.

Notes

Files written by problemSetupGas()

  • lines.inp : Line mode master command file.
  • numberdens_xxx.inp : Number density of molecule/atomic species ‘xxx’
  • gas_velocity.inp : Gas velocity
  • microturbulence.inp : The standard deviation of the Gaussian line profile caused by turbulent
    broadening.
  • gas_temperature.inp : Gas temperature (which may be different from the dust temperature). If
    tgas_eq_tdust is set to zero in radmc3d.inp the gas temperature in this file will be used instead of the dust temperature.

If fullsetup is set to True the following additional files will be created

  • amr_grid.inp : Spatial grid.
  • wavelength_micron.inp : Wavelength grid.
  • stars.inp : Input radiation field.
  • radmc3d.inp : Parameters for RADMC-3D (e.g. Nr of photons to be used, scattering type, etc).
class radmc3dPy.setup.radmc3dModel(model=None, binary=None, old=False, dfunc=None, dfpar=None, parfile_update=True, **kwargs)

Bases: object

Parameters:

model : str

Name of the model that should be used to create the density structure. The file should be in a directory from where it can directly be imported (i.e. the directory should be in the PYTHON_PATH environment variable or it should be in the current working directory) and the file name should be ‘model_xxx.py’, where xxx stands for the string that should be specified in this variable

binary : bool, optional

If True input files will be written in binary format, if False input files are written as formatted ascii text.

old : bool, optional

If set to True the input files for the old 2D version of radmc will be created

dfunc : function, optional

Decision function for octree-like amr tree building. It should take linear arrays of cell centre coordinates (x,y,z) and cell half-widhts (dx,dy,dz) in all three dimensions, a radmc3d model, a dictionary with all parameters from problem_params.inp and an other keyword argument (**kwargs). It should return a boolean ndarray of the same length as the input coordinates containing True if the cell should be resolved and False if not. An example for the implementation of such decision function can be found in radmc3dPy.analyze module (radmc3dPy.analyze.gdensMinMax()).

dfpar : dictionary

Dicionary of keyword arguments to be passed on to dfunc. These parameters will not be written to problem_params.inp. Parameters can also be passed to dfunc via normal keyword arguments gathered in **kwargs, however all keyword arguments in **kwargs will be written to problem_params.inp

parfile_update : bool

If True the parameter file (problem_params.inp) will be updated / overwritten with any parameters that are possibly passed as keyword arguments. For False the problem_params.inp file will not be overwritten irrespectively of the parameter values the data members of the radmc3dModel instance would contain.

Keyword Arguments :

Any varible name in problem_params.inp can be used as a keyword argument. At first all variables are read from problem_params.in to a dictionary called ppar. Then if there is any keyword argument set in the call of problem_setup_dust the ppar dictionary is searched for this key. If found the value belonging to that key in the ppar dictionary is changed to the value of the keyword argument. If no such key is found then the dictionary is simply extended by the keyword argument. Finally the problem_params.inp file is updated with the new parameter values.

Attributes

binary (bool) If True, it sets the data file format to binary (overwriting also the rto_style parameter!), if False the file format will be formatted ascii
data (radmc3dData) Container for physical variables in the mode
grid (radmc3dGrid) Container for spatial and wavelength grid
model (str) Name of the model that has callable functions to generate physical variables for the model
old (bool) If True the model is meant to be for radmc, the pre-decessor code of radmc-3d.
opac (radmc3dOpac) Container for dust opacities
par (radmc3dPar) Container for the parameters of the model (i.e. the content of the problem_params.inp file)
radsources (radmc3dRadsources) Container for the radiation sources in the model

Methods

makeDustOpac() Generates dust opacities and writes them to file
makeGrid([sgrid, wgrid, writeToFile]) Generates a spatial and/or wavelength grid
makeRadSources([writeToFile])
Parameters:
makeVar([ddens, dtemp, gdens, gtemp, gvel, ...]) Generates variables and possibly writes them to file
readParams()
setupDust([sgrid, wgrid, radsources, ...]) Set up a model for RADMC-3D.
setupGas([sgrid, wgrid, radsources, gdens, ...]) Set up a model for RADMC-3D.
writeLinesInp() Writes the lines.inp master command file for line simulation in RADMC-3D
writeRadmc3dInp() Writes the radmc3d.inp master command file for RADMC-3D
makeDustOpac()

Generates dust opacities and writes them to file

makeGrid(sgrid=True, wgrid=True, writeToFile=False, **kwargs)

Generates a spatial and/or wavelength grid

Parameters:

wgrid : bool

Set to True to generate the wavelength grid

sgrid : bool

Set to True to generate the spatial grid

writeToFile : bool

If True the grid will be written to amr_grid.inp and/or wavelength_micron.inp

**kwargs : Any varible name in problem_params.inp can be used as a keyword argument.

At first all variables are read from problem_params.in to a dictionary called ppar. Then if there is any keyword argument set in the call of problem_setup_dust the ppar dictionary is searched for this key. If found the value belonging to that key in the ppar dictionary is changed to the value of the keyword argument. If no such key is found then the dictionary is simply extended by the keyword argument. Finally the problem_params.inp file is updated with the new parameter values.

makeRadSources(writeToFile=True)
Parameters:

writeToFile : bool

If True the radiation will be written to stars.inp and or stellarsrc_density.inp, stellarsrc_templates.inp

makeVar(ddens=False, dtemp=False, gdens=False, gtemp=False, gvel=False, vturb=False, writeToFile=False)

Generates variables and possibly writes them to file

Parameters:

ddens : bool, optional

If True generates the dust density

dtemp : bool, optional

If True generates the dust temperature (normally this would be calculated by RADMC-3D, but

in some cases, e.g. debugging, it might come handy)

gdens : bool, optional

If True generates the gas density

gtemp : bool, optional

If True generates the gas temperature

gvel : bool, optional

If True generates the gas velocity

vturb : bool, optional

If True generates the microturbulent velocity field

writeToFile : bool, optional

If True the generated variable(s) will be written to file

readParams()
setupDust(sgrid=True, wgrid=True, radsources=True, dustopac=True, ddens=True, dtemp=False)

Set up a model for RADMC-3D.

Parameters:

sgrid : bool

Should the spatial grid be generated? (default=True)

wgrid : bool

Should the wavelength grid be generated? (default=True)

radsources : bool

Should the radiation sources be generated? (default=True)

dustopac : bool

Should dust opacities be generated? (default=True)

ddens : bool

Should the dust density be generated? (default=True)

dtemp : bool

Should the dust temperature be generated (default=False)

Note, that by default all generated data will be written to files immediately. In case you wish to generate

the files but not write the data to files but only keep them instead in data attributes, please use the individual generator functions methods (e.g. makeGrid, makeRadSources, makeVar, etc)

setupGas(sgrid=False, wgrid=False, radsources=False, gdens=True, gtemp=True, gvel=True, vturb=True)

Set up a model for RADMC-3D.

Parameters:

sgrid : bool

Should the spatial grid be generated? (default=False)

wgrid : bool

Should the wavelength grid be generated? (default=False)

radsources : bool

Should the radiation sources be generated? (default=False)

gdens : bool

Should the gas density be generated? (default=True)

gtemp : bool

Should the gas temperature be generated (default=False)

gvel : bool

Should the gas veloctiy be generated (default=False)

vturb : bool

Should the turbulent velocity field be generated (default=False)

Note, that by default all generated data will be written to files immediately. In case you wish to generate

the files but not write the data to files but only keep them instead in data attributes, please use the individual generator functions methods (e.g. makeGrid, makeRadSources, makeVar, etc)

writeLinesInp()

Writes the lines.inp master command file for line simulation in RADMC-3D

writeRadmc3dInp()

Writes the radmc3d.inp master command file for RADMC-3D

radmc3dPy.setup.validateModel(model='', dustModel=False, gasModel=False, writeDustTemp=False, octree=False)

Function to validate a model. It checks three things: 1) whether or not the model can be imported, 2) whether the model has all the function to be used as dust and/or gas model, 3) if it has the right number of arguments. The function names tested are getDefaultParams, getDustDensity, getGasDensity, getGasAbundance, getVTurb, getVelocity, getDustTempearture (optional).

Parameters:

model : str

Name of the model to be tested

dustModel : bool

If True the existence of functions getDustDensity() and getDustTemperature() will be checked. The latter is only checked if writeDustTemp is set to True.

gasModel : bool

If True the existence of functions getGasDensity(), getGasAbundance(), getVTurb(), getVelocity() will be checked.

writeDustTemp: bool

If True the existence of the function getDustTemperature() will be checked.

octree : bool

If True the number of argument of the model functions will be checked. For regular grids only two arguments should be present for the grid instance and for the parameter dictionary (grid, ppar). For a model to be used with octree AMR three additional arguments for the three spatial coordiantes (x,y,z) should be present. The argument sequence should then be x, y, z, grid, ppar.

Returns:

A boolean True if the model is valid and False if it is not.

radmc3dPy.setup.writeLinesInp(ppar=None)

Writes the lines.inp master command file for line simulation in RADMC-3D

Parameters:

ppar : dictionary,

Contains all parameters of a RADMC-3D run

radmc3dPy.setup.writeRadmc3dInp(modpar=None)

Writes the radmc3d.inp master command file for RADMC-3D

Parameters:

modpar : radmc3dPar

Contains all parameters of a RADMC-3D run.

radmc3dPy.setup.writeRadmcInp(modpar=None, nphot=None)

Writes the radmc.inp master command file for the 2D version of radmc

Parameters:

modpar : radmc3dPar

Contains all parameters of a radmc run.

nphot : int

Number of photons used for the MC simulation

radmc3dPy.staratm module

This module contains functions to get various radiation sources in RADMC-3D. For help on the syntax or functionality of each function see the help of the individual functions

radmc3dPy.staratm.getAtmModel(teff=0.0, logg=None, mstar=None, lstar=None, rstar=None, iwav=None, model='kurucz', modeldir=None, wmax=7.0)

Interpolates the stellar model atmosphere on a pre-defined wavelength grid The model atmospheres are interpolated in logg and Teff and rebinned in wavelength to the input wavelength grid and scaled to have the same luminosity as specified

Parameters:

teff : float

Effective temperature of the star

logg : float

Logarithm of the surface gravity of the star

mstar : float

Mass of the star in gramm

lstar : float, optional

Luminosity of the star (either lstar or rstar should be specified)

rstar : float, optional

Radius of the star (either lstar or rstar should be specified)

iwav : ndarray

Wavelength grid for which the stellar atmosphere model should be calculated

model : {‘kurucz’, ‘nextgen’, ‘ames-dusty’}

Name of the model atmosphere family

modeldir : str

Path to the atmosphere models

wmax : float

Maximum wavelength until the model atmosphere is used on the interpolated grid Longwards of this wavelength the Rayleigh-Jeans approximation (lambda^(-2)) is used.

Returns:

Returns a dictionary with the following keys:

  • wav : ndarray
    Wavelength in micron (same as the input iwav keyword)
  • lnu : ndarray
    Monochromatic luminosity of the stellar atmosphere in erg/s/Hz
radmc3dPy.staratm.getSpectrumAmesDusty(teff=None, logg=None, mstar=None, lstar=None, rstar=None, modeldir=None, wav=None)

Interpolates the Ames-Dusty model atmospheres in logg and Teff

Parameters:

teff : float

Effective temperature of the star

logg : float

Logarithm of the surface gravity of the star

mstar : float, optional

Stellar mass in g (either logg or mstar and rstar should be specified)

lstar : float, optional

Luminosity of the star (either lstar or rstar should be specified)

rstar : float, optional

Radius of the star (either lstar or rstar should be specified)

modeldir : str

Path to the AmesDusty atmosphere model grid

wav : ndarray

Wavelength grid for which the stellar atmosphere model should be calculated

Returns:

Returns a dictionary with the following keys:

  • wav : ndarray
    Wavelength in micron (same as the input wav keyword)
  • lnu : ndarray
    Monochromatic luminosity of the stellar atmosphere in erg/s/Hz
  • bnu : ndarray
    Monochromatic luminosity of a blackbody stellar atmosphere with the same luminosity and effective temperature as the stellar model in erg/s/Hz
radmc3dPy.staratm.getSpectrumKurucz(teff=None, logg=None, mstar=None, lstar=None, rstar=None, modeldir=None, wav=None)

Interpolates the Kurucz model atmospheres in logg and Teff

Parameters:

teff : float

Effective temperature of the star

logg : float

Logarithm of the surface gravity of the star

mstar : float, optional

Stellar mass in g (either logg or mstar and rstar should be specified)

lstar : float, optional

Luminosity of the star (either lstar or rstar should be specified)

rstar : float, optional

Radius of the star (either lstar or rstar should be specified)

modeldir : str

Path to the Kurucz atmosphere model grid

wav : ndarray

Wavelength grid for which the stellar atmosphere model should be calculated

Returns:

Returns a dictionary with the following keys:

  • wav : ndarray
    Wavelength in micron (same as the input wav keyword)
  • lnu : ndarray
    Monochromatic luminosity of the stellar atmosphere in erg/s/Hz
  • lnucont : ndarray
    Monochromatic luminosity of the continuum stellar atmosphere in erg/s/Hz
radmc3dPy.staratm.getSpectrumNextGen(teff=None, logg=None, mstar=None, lstar=None, rstar=None, modeldir=None, wav=None)

Interpolates the NextGen model atmospheres in logg and Teff

Parameters:

teff : float

Effective temperature of the star

logg : float

Logarithm of the surface gravity of the star

mstar : float, optional

Stellar mass in g (either logg or mstar and rstar should be specified)

lstar : float, optional

Luminosity of the star (either lstar or rstar should be specified)

rstar : float, optional

Radius of the star (either lstar or rstar should be specified)

modeldir : str

Path to the NextGen atmosphere model grid

wav : ndarray

Wavelength grid for which the stellar atmosphere model should be calculated

Returns:

Returns a dictionary with the following keys:

  • wav : ndarray
    Wavelength in micron (same as the input wav keyword)
  • lnu : ndarray
    Monochromatic luminosity of the stellar atmosphere in erg/s/Hz
  • bnu : ndarray
    Monochromatic luminosity of a blackbody stellar atmosphere with the same luminosity and effective temperature as the stellar model in erg/s/Hz
radmc3dPy.staratm.readAmesDustySpectrum(fname='')

Reads the Ames-Dusty model atmosphere

Parameters:

fname : str

File name containing the Ames-Dusty model atmosphere (e.g. lte30-3.5-0.0.AMES-dusty.7)

Returns:

Returns a dictionary with the following keys:

  • wav : ndarray
    Wavelength in micron
  • nwav : int
    Number of wavelength points
  • inu : ndarray
    Specific intensity of the stellar model atmosphere in erg/s/cm/cm/Hz/ster
  • bnu : list
    Specific intensity of a blackbody stellar model atmosphere with the same luminosity and the same effective temperature in erg/s/cm/cm/Hz/ster
  • teff : float
    Effective temperature of the model
  • logg : float
    Logarithm of the surface gravity of the model
  • mph : float
    Metallicity of the atmosphere model
radmc3dPy.staratm.readKuruczGrid(fname='')

Reads the Kurucz model atmosphere grid. It reads a whole file from the Kurucz grid that contains a 2D grid of atmosphere models for Teff and logg.

Parameters:

fname : str

File name containing the Kurucz model atmosphere (e.g. fp00k2.pck)

Returns:

Returns a dictionary with the following keys:

  • wav : ndarray
    Wavelength in micron
  • nwav : int
    Number of wavelength points
  • inu : list
    Each element of the list contains an ndarray with the specific intensity of the stellar model atmosphere in erg/s/cm/cm/Hz/ster
  • inucont : list
    Each element of the list contains an ndarray with the specific intensity of the stellar model atmosphere continuum in erg/s/cm/cm/Hz/ster
  • teff : list
    Contains the Teff grid of the model grid
  • logg : list
    Contains the logg grid of the model grid
radmc3dPy.staratm.readNextGenSpectrum(fname='')

Reads the NextGen model atmosphere.

Parameters:

fname : str

File name containing the NextGen model atmosphere (e.g. nlte98-4.5-2.5.NextGen.spec)

Returns:

Returns a dictionary with the following keys:

  • wav : ndarray
    Wavelength in micron
  • nwav : int
    Number of wavelength points
  • inu : ndarray
    Specific intensity of the stellar model atmosphere in erg/s/cm/cm/Hz/ster
  • bnu : list
    Specific intensity of a blackbody stellar model atmosphere with the same luminosity and the same effective temperature in erg/s/cm/cm/Hz/ster
  • teff : float
    Effective temperature of the model
  • logg : float
    Logarithm of the surface gravity of the model
  • mph : float
    Metallicity of the atmosphere model
radmc3dPy.staratm.rebinSpectrum(wav=None, fnu=None, iwav=None)

Rebins the spectrum to a coarser wavelength grid

Parameters:

wav : ndarray

Wavelength grid of the spectrum to be rebinned

fnu : ndarray

Wavelength dependent spectrum (e.g. specific intensity, monochromatic luminosity etc) to be rebbinned

iwav : ndarray

Wavelength grid onto which the spectrum should be rebinned (it is assumed to be logarithmically spaced)

Returns:

Returns an ndarray containing the spectrum rebinned to the input wavelength grid

Module contents

RADMC-3D Python module (c) Attila Juhasz, 2011-2018