magritte.core module

core.__version__ = '0.5.8'

Current version of Magritte.

Core module of Magritte: a modern software library for 3D radiative transfer.

class magritte.core.Boundary

Bases: pybind11_builtins.pybind11_object

Class containing the model boundary.

get_boundary_condition(self: magritte.core.Boundary, arg0: int) magritte.core.BoundaryCondition

Getter for the boundary condition.

read(self: magritte.core.Boundary, arg0: magritte.core.Io) None

Read object from file.

set_boundary_condition(self: magritte.core.Boundary, arg0: int, arg1: magritte.core.BoundaryCondition) magritte.core.BoundaryCondition

Setter for the boundary condition.

write(self: magritte.core.Boundary, arg0: magritte.core.Io) None

Write object to file.

property boundary2point

Array with point index for each boundary point.

property boundary_temperature

Array with radiative temperature for each boundary point (only relevant for thermal boundary conditions).

property point2boundary

Array with boundary index for each point.

class magritte.core.Chemistry

Bases: pybind11_builtins.pybind11_object

Class containing the chemistry.

read(self: magritte.core.Chemistry, arg0: magritte.core.Io) None

Read object from file.

write(self: magritte.core.Chemistry, arg0: magritte.core.Io) None

Write object to file.

property species

Species object.

class magritte.core.CollisionPartner

Bases: pybind11_builtins.pybind11_object

Class containing collision partner data.

read(self: magritte.core.CollisionPartner, arg0: magritte.core.Io, arg1: int, arg2: int) None

Read object from file.

write(self: magritte.core.CollisionPartner, arg0: magritte.core.Io, arg1: int, arg2: int) None

Write object to file.

property Cd

Array with collisional de-excitation rates.

property Ce

Array with collisional excitation rates.

property icol

Array with upper levels of the collisional transitions.

property jcol

Array with lower levels of the collisional transitions.

property ncol

Number of collisional transitions.

property ntmp

Number of temperatures at which data is given.

property num_col_partner

Number of the species corresponding to the collision partner.

property orth_or_para_H2

In case the collision partner is H2, this indicates whether it is ortho (o) or para (p).

property tmp

Array with temperatures corresponding to the collisional data.

class magritte.core.Frequencies

Bases: pybind11_builtins.pybind11_object

Class containing the (local) frequency bins.

read(self: magritte.core.Frequencies, arg0: magritte.core.Io) None

Read object from file.

write(self: magritte.core.Frequencies, arg0: magritte.core.Io) None

Write object to file.

property nu

Array of frequency bins for each point in the model.

class magritte.core.Geometry

Bases: pybind11_builtins.pybind11_object

Class containing the model geometry.

read(self: magritte.core.Geometry, arg0: magritte.core.Io) None

Read object from file.

write(self: magritte.core.Geometry, arg0: magritte.core.Io) None

Write object to file.

property boundary

Boundary object

property lengths

Array containing the lengths of the rays for each direction and point.

property points

Points object.

property rays

Rays object.

class magritte.core.Image

Bases: pybind11_builtins.pybind11_object

Image class, 2D point cloud of intensities for each frequency bin.

property I

Intensity of the points in the image (for each frequency bin).

property ImX

X-coordinates of the points in the image plane.

property ImY

Y-coordinates of the points in the image plane.

property freqs

Frequency bins of the image.

property imagePointPosition

Position of image points (model points or projection surface).

property imageType

Type of image (intensity or optical depth).

property image_direction_x

coordinate of image x direction in a 3D general geometry.

property image_direction_y

coordinate of image y direction in a 3D general geometry.

property image_direction_z

direction in which the image is taken.

property nfreqs

Number of frequency bins in the image.

property ray_nr

Number of the ray along which the image is taken.

class magritte.core.ImagePointPosition

Bases: pybind11_builtins.pybind11_object

Members:

AllModelPoints

ProjectionSurface

property name
class magritte.core.ImageType

Bases: pybind11_builtins.pybind11_object

Members:

Intensity

OpticalDepth

property name
class magritte.core.Io

Bases: pybind11_builtins.pybind11_object

Abstract input/output base class.

class magritte.core.IoPython

Bases: magritte.core.Io

Input/output class for io, implemented in Python.

property implementation

Type of input/output. Either “hdf5” or “text”.

property io_file

Name of the folder or file to read from or write to.

class magritte.core.IoText

Bases: magritte.core.Io

Intput/output class for text-based io, implemented in C++.

property io_file

Name of the folder or file to read from or write to.

class magritte.core.LineProducingSpecies

Bases: pybind11_builtins.pybind11_object

Class containing a line producing species.

index(self: magritte.core.LineProducingSpecies, arg0: int, arg1: int) int
read(self: magritte.core.LineProducingSpecies, arg0: magritte.core.Io, arg1: int) None

Read object from file.

write(self: magritte.core.LineProducingSpecies, arg0: magritte.core.Io, arg1: int) None

Write object to file.

property J

Isotropic radiation field.

property J2_0

Anisotropic radiation field tensor element 0

property J2_1_Im

Anisotropic radiation field tensor element 1 (imaginary part)

property J2_1_Re

Anisotropic radiation field tensor element 1 (real part)

property J2_2_Im

Anisotropic radiation field tensor element 2 (imaginary part)

property J2_2_Re

Anisotropic radiation field tensor element 2 (real part)

property Jeff

Array with effective mean intensities in the lines.

property linedata

Linedata object.

property population

Array with level populations for each point.

property population_tot

Array with the sum of all level populations at each point. (Should be equal to the abundance of the species.)

property quadrature

Quadrature object.

class magritte.core.Linedata

Bases: pybind11_builtins.pybind11_object

Class containing line data.

read(self: magritte.core.Linedata, arg0: magritte.core.Io, arg1: int) None

Read object from file.

write(self: magritte.core.Linedata, arg0: magritte.core.Io, arg1: int) None

Write object to file.

property A

Array with Einstein A coefficients.

property Ba

Array with Einstein B (absorption) coeficients.

property Bs

Array with Einstsin B (stimulated emission) coefficients.

property colpar

Vector of collision partner objects.

property energy

Array with the energy for each level.

property frequency

Array with frequencies of the line transitions.

property inverse_mass

Inverse mass (1/mass) of the species.

property irad

Array with upper levels of the radiative transitions.

property jrad

Array with lower levels of the radiative transitions.

property ncolpar

Number of collision partners.

property nlev

Number of energy levels.

property nrad

Number of radiative transitions.

property num

Number of the species corresponding the line data.

property sym

Chemical symbol of the species.

property weight

Array with the statistical weight for each level.

class magritte.core.Lines

Bases: pybind11_builtins.pybind11_object

Class containing the lines and their data.

read(self: magritte.core.Lines, arg0: magritte.core.Io) None

Read object from file.

resize_LineProducingSpecies(self: magritte.core.Lines, arg0: int) None

Resize the vector LineProducingSpecies.

set_emissivity_and_opacity(self: magritte.core.Lines) None

Set the emissivity and opacity arrays.

write(self: magritte.core.Lines, arg0: magritte.core.Io) None

Write object to file.

property emissivity

Array with emissivities for each point and frequency bin.

property inverse_width

Array with the inverse widths for each line at each point.

property line

Array with the line (centre) frequencies.

property lineProducingSpecies

Vector of LineProducingSpecies.

property opacity

Array with opacities for each point and frequency bin.

class magritte.core.Model

Bases: pybind11_builtins.pybind11_object

Class containing the Magritte model.

compute_Jeff(self: magritte.core.Model) int

Compute the effective mean intensity in the line.

compute_LTE_level_populations(self: magritte.core.Model) int

Compute the level populations for the model assuming local thermodynamic equilibrium (LTE).

compute_image(self: magritte.core.Model, arg0: int) int

Compute an image for the model along the given ray.

compute_image_for_point(self: magritte.core.Model, arg0: int, arg1: int) int

Compute image (single pixel) for a single point.

compute_image_new(*args, **kwargs)

Overloaded function.

  1. compute_image_new(self: magritte.core.Model, arg0: int) -> int

Compute an image of the model along the given ray direction, using the new imager.

  1. compute_image_new(self: magritte.core.Model, arg0: int, arg1: int, arg2: int) -> int

Compute an image of the model along the given ray direction, using the new imager, specifying the image resolution.

  1. compute_image_new(self: magritte.core.Model, arg0: float, arg1: float, arg2: float, arg3: int, arg4: int) -> int

Compute an image of the model along the given ray direction, using the new imager, specifying the ray direction and image resolution.

compute_image_optical_depth(self: magritte.core.Model, arg0: int) int

Compute an image of the optical depth for the model along the given ray.

compute_image_optical_depth_new(*args, **kwargs)

Overloaded function.

  1. compute_image_optical_depth_new(self: magritte.core.Model, arg0: int) -> int

Compute an image of the optical depth for the model along the given ray direction, using the new imager.

  1. compute_image_optical_depth_new(self: magritte.core.Model, arg0: int, arg1: int, arg2: int) -> int

Compute an image of the optical depth for the model along the given ray direction, using the new imager, specifying the image resolution.

  1. compute_image_optical_depth_new(self: magritte.core.Model, arg0: float, arg1: float, arg2: float, arg3: int, arg4: int) -> int

Compute an image of the optical depth for the model along the given ray direction, using the new imager, specifying the ray direction and image resolution.

compute_inverse_line_widths(self: magritte.core.Model) int

Compute the inverse line widths for the model. (Needs to be recomputed whenever the temperature of turbulence changes.)

compute_level_populations(self: magritte.core.Model, arg0: bool, arg1: int) int

Compute the level populations for the model assuming statistical equilibrium until convergence, optionally using Ng-acceleration, and for the given maximum number of iterations.

compute_level_populations_from_stateq(self: magritte.core.Model) int

Compute the level populations for the model assuming statistical equilibrium.

compute_level_populations_shortchar(self: magritte.core.Model, arg0: bool, arg1: int) int

Compute the level populations using the short-characteristics solver for the model assuming statistical equilibrium until convergence, optionally using Ng-acceleration, and for the given maximum number of iterations.

compute_level_populations_sparse(self: magritte.core.Model, arg0: bool, arg1: int) int

Compute the level populations for the model assuming statistical equilibrium until convergence, optionally using Ng-acceleration, and for the given maximum number of iterations. (Memory sparse option.)

compute_radiation_field_feautrier_order_2(self: magritte.core.Model) int

Compute the radiation field for the modle using the 2nd-order Feautrier solver.

compute_radiation_field_feautrier_order_2_anis(self: magritte.core.Model) int

Compute the radiation field for the modle using the 2nd-order Feautrier solver, anisotropic case.

compute_radiation_field_feautrier_order_2_sparse(self: magritte.core.Model) int

Compute the radiation field for the modle using the 2nd-order Feautrier solver.

compute_radiation_field_shortchar_order_0(self: magritte.core.Model) int

Compute the radiation field for the modle using the 0th-order short-characteristics methods.

compute_spectral_discretisation(*args, **kwargs)

Overloaded function.

  1. compute_spectral_discretisation(self: magritte.core.Model) -> int

Compute the spectral discretisation for the model tailored for line (Gauss-Hermite) quadrature.

  1. compute_spectral_discretisation(self: magritte.core.Model, arg0: float) -> int

Compute the spectral discretisation for the model tailored for images with the given spectral width.

  1. compute_spectral_discretisation(self: magritte.core.Model, arg0: float, arg1: float) -> int

Compute the spectral discretisation for the model tailored for images with the given min and max frequency.

  1. compute_spectral_discretisation(self: magritte.core.Model, arg0: float, arg1: float, arg2: int) -> int

Compute the spectral discretisation for the model tailored for images with the given min and max frequency. Can also specify the amount of frequency bins to use (instead of defaulting to parameters.nfreqs).

read(*args, **kwargs)

Overloaded function.

  1. read(self: magritte.core.Model) -> None

Read model file (using the hdf5 implementation of IoPython and using the model name specified in parameters.)

  1. read(self: magritte.core.Model, arg0: magritte.core.Io) -> None

Read model file using the given Io object.

  1. read(self: magritte.core.Model, arg0: str) -> None

Read model file (assuming HDF5 file format).

set_boundary_condition(self: magritte.core.Model) int

Set boundary condition (internally).

set_column(self: magritte.core.Model) int

Set column (internally).

set_eta_and_chi(self: magritte.core.Model, arg0: int) int

Set latest emissivity and opacity for the model in the eta and chi variables respectively.

write(*args, **kwargs)

Overloaded function.

  1. write(self: magritte.core.Model) -> None

Write model file (using the hdf5 implementation of IoPython and using the model name specified in parameters.)

  1. write(self: magritte.core.Model, arg0: magritte.core.Io) -> None

Write model file using the given Io object.

  1. write(self: magritte.core.Model, arg0: str) -> None

Write model file (assuming HDF5 file format).

property error_max

Max relative error in the level populaitons.

property error_mean

Mean relative error in the level populations.

class magritte.core.Parameters

Bases: pybind11_builtins.pybind11_object

Class containing the model parameters.

adaptive_ray_tracing(self: magritte.core.Parameters) bool

Whether or not to use adaptive ray tracing.

dimension(self: magritte.core.Parameters) int

Spatial dimension of the model.

hnrays(self: magritte.core.Parameters) int

Half the number of rays.

model_name(self: magritte.core.Parameters) str

Model name.

nboundary(self: magritte.core.Parameters) int

Number of boundary points.

nfreqs(self: magritte.core.Parameters) int

Number of frequency bins.

nlines(self: magritte.core.Parameters) int

Number of lines.

nlspecs(self: magritte.core.Parameters) int

Number of line producing species.

npoints(self: magritte.core.Parameters) int

Number of points.

nquads(self: magritte.core.Parameters) int

Number of quadrature points.

nrays(self: magritte.core.Parameters) int

Number of rays.

nspecs(self: magritte.core.Parameters) int

Number of species.

read(self: magritte.core.Parameters, arg0: magritte.core.Io) None

Rread object from file.

set_adaptive_ray_tracing(self: magritte.core.Parameters, arg0: bool) None

Set whether or not to use adaptive ray tracing.

set_dimension(self: magritte.core.Parameters, arg0: int) None

Set spatial dimension of the model.

set_hnrays(self: magritte.core.Parameters, arg0: int) None

Set the half of the number of rays.

set_model_name(self: magritte.core.Parameters, arg0: str) None

Set model name.

set_nboundary(self: magritte.core.Parameters, arg0: int) None

Set number of boundary points.

set_nfreqs(self: magritte.core.Parameters, arg0: int) None

Set number of frequency bins.

set_nlines(self: magritte.core.Parameters, arg0: int) None

Set number of lines.

set_nlspecs(self: magritte.core.Parameters, arg0: int) None

Set number of line producing species.

set_npoints(self: magritte.core.Parameters, arg0: int) None

Set number of points.

set_nquads(self: magritte.core.Parameters, arg0: int) None

Set number of quadrature points.

set_nrays(self: magritte.core.Parameters, arg0: int) None

Set number of rays.

set_nspecs(self: magritte.core.Parameters, arg0: int) None

Set number of species.

set_spherical_symmetry(self: magritte.core.Parameters, arg0: bool) None

Set whether or not to use spherical symmetry

set_use_scattering(self: magritte.core.Parameters, arg0: bool) None

Set whether or not to use scattering.

spherical_symmetry(self: magritte.core.Parameters) bool

Whether or not to use spherical symmetry.

use_scattering(self: magritte.core.Parameters) bool

Whether or not to use scattering.

version(self: magritte.core.Parameters) str

Magritte version number.

write(self: magritte.core.Parameters, arg0: magritte.core.Io) None

Write object to file.

property Ng_acceleration_mem_limit

Memory limit for ng acceleration.

property Ng_acceleration_remove_N_its

Number of iterations to ignore when using ng-acceleration

property adaptive_Ng_acceleration_min_order

Minimal order of ng acceleration when using adaptive Ng acceleration. Has to be larger than 1.

property adaptive_Ng_acceleration_use_max_criterion

Whether of not to use the maximum relative change as criterion for the adaptive ng-acceleration

property convergence_fraction

Fraction of levels that should obey the convergence criterion.

property max_distance_opacity_contribution

The distance (scaled to line width) at which we ignore the opacity/emissivity contribution.

property max_width_fraction

Max tolerated Doppler shift as fraction of line width (default=0.5).

property min_dtau

Minimum optical depth increment that will be assumed in the solver.

property min_opacity

Minimum opacity that will be assumed in the solver.

property min_rel_pop_for_convergence

Minimum relative level population to be considered in the convergence criterion.

property n_off_diag

Bandwidth of the ALO (0=diagonal, 1=tri-diaginal, 2=penta-diagonal, …)

property one_line_approximation

Whether or not to use one line approximation.

property pop_prec

Required precision for ALI.

property population_inversion_fraction

Threshold factor for population inversion required for LTE to be used; set this higher than 1

property store_intensities

Whether or not to store intensities.

property sum_opacity_emissivity_over_all_lines

Whether or not to sum the opacity, emissivity over all lines (including the zero contributions).

property use_adaptive_Ng_acceleration

Whether to use adaptive Ng acceleration.

class magritte.core.Points

Bases: pybind11_builtins.pybind11_object

Class containing the spatial points.

read(self: magritte.core.Points, arg0: magritte.core.Io) None

Read object from file.

write(self: magritte.core.Points, arg0: magritte.core.Io) None

Write object to file.

property cum_n_neighbors

Cumulative number of neighbours of each point.

property n_neighbors

Number of neighbours of each point.

property neighbors

Linearised array of neighbours of each point.

property position

Array with position vectors of the points.

property velocity

Array with velocity vectors of the points (as a fraction of the speed of light).

class magritte.core.Quadrature

Bases: pybind11_builtins.pybind11_object

Class containing the data for Gauss-Hermite quadrature.

read(self: magritte.core.Quadrature, arg0: magritte.core.Io, arg1: int) None

Read object from file.

write(self: magritte.core.Quadrature, arg0: magritte.core.Io, arg1: int) None

Write object to file.

property roots

Array containing the roots for the Gauss-Hermite quadrature.

property weights

Array containing the weights for the Gauss-Hermite quadrature.

class magritte.core.Radiation

Bases: pybind11_builtins.pybind11_object

Class containing the radiation field.

read(self: magritte.core.Radiation, arg0: magritte.core.Io) None

Read object from file.

write(self: magritte.core.Radiation, arg0: magritte.core.Io) None

Write object to file.

property I

Array containing the intensity for each ray, point, and frequency bin.

property J

Array containing the mean intensity for each point and frequency bin.

property frequencies

Frequencies object.

property u

Array containing the mean intensity up and down a ray for each ray pair, point, and frequency bin.

property v

Array containing the flux up and down a ray for each ray pair, point, and frequency bin.

class magritte.core.Rays

Bases: pybind11_builtins.pybind11_object

Class containing the (light) rays (directional discretisation).

read(self: magritte.core.Rays, arg0: magritte.core.Io) None

Read object from file.

write(self: magritte.core.Rays, arg0: magritte.core.Io) None

Write object to file.

property antipod

Array with the number of the antipodal ray for each ray.

property direction

Array with direction vector of each ray.

property weight

Array with the weights that each ray contributes in integrals over directions.

class magritte.core.Species

Bases: pybind11_builtins.pybind11_object

Class containing the chemical species.

read(self: magritte.core.Species, arg0: magritte.core.Io) None

Read object from file.

write(self: magritte.core.Species, arg0: magritte.core.Io) None

Write object to file.

property abundance

Array with the abundances at each point.

property symbol

Symbol of the species.

class magritte.core.Temperature

Bases: pybind11_builtins.pybind11_object

Class containing the temperature.

read(self: magritte.core.Temperature, arg0: magritte.core.Io) None

Read object from file.

write(self: magritte.core.Temperature, arg0: magritte.core.Io) None

Write object to file.

property gas

Kinetic temperature of the gas.

class magritte.core.Thermodynamics

Bases: pybind11_builtins.pybind11_object

Class containing the thermodynamics.

read(self: magritte.core.Thermodynamics, arg0: magritte.core.Io) None

Read object from file.

write(self: magritte.core.Thermodynamics, arg0: magritte.core.Io) None

Write object to file.

property temperature

Temperature object.

property turbulence

Turbulence object.

class magritte.core.Turbulence

Bases: pybind11_builtins.pybind11_object

Class containing the (micro) turbulence.

read(self: magritte.core.Turbulence, arg0: magritte.core.Io) None

Read object from file.

write(self: magritte.core.Turbulence, arg0: magritte.core.Io) None

Write object to file.

property vturb2

Square of the micro turbulence as a fraction of the speed of light.

magritte.core.pcmp_comm_rank() int

Get the rank of the current process (using MPI).

magritte.core.pcmp_comm_size() int

Get the size of the current communicator (using MPI).

magritte.core.pcmp_length(arg0: int) int

Get the size of an array with given length that is given to the current process (using MPI).

magritte.core.pcmp_start(arg0: int) int

Get the first index of an array with given length that will be given to this process (using MPI).

magritte.core.pcmp_stop(arg0: int) int

Get the last index of an array with given length that will be given to this process (using MPI).

magritte.core.pcmt_n_threads_avail() int

Get the number of available threads (using OpenMP).

magritte.core.pcmt_set_n_threads_avail(arg0: int) None

Set the number of available threads (using OpenMP).