Emitter module#

lineshapes#

cherab.iter.emitter.lineshapes.gauss_ls(mass, wavelength, tion, wlcenters, intensities)#

Multi-gaussian 2D (space, wavelength) line shape profile

Parameters:
  • mass (float) – Isotope mass in eV

  • wavelength (ndarray) – Wavelength array (nm)

  • tion (ndarray) – Ion temperature profile (eV)

  • wlcenters (ndarray) – Spatial profile of the centers of single gaussian components (nm)

  • intensities (ndarray) – Spatial profile of the intensities of single gaussian components

Returns:

2D gaussian profile

Return type:

ndarray

class cherab.iter.emitter.lineshapes.AdasLineShape(path='/home/koyo/Documents/cherab/cherab_iter/cherab/iter/emitter/adas/data/adas603', mass=None, element='Be', line='527nm')#

Doppler-Zeeman spectral line shape with the fine structure taken from ADAS (ADAS603)

Parameters:
  • path (str) – Path to the profiles extracted from ADAS: Wavelength(Bmag), Intensity_parallel(Bmag), Intensity_perpendicular(Bmag)

  • mass (float | None) – Element/isotope mass in eV

  • element (str) – Element/isotope symbol, by default “Be”.

  • line (str) – Spectral line, by default “527nm”.

__call__(wavelength, tion, bmag, costheta)#

Returns spectral emission profile with respect to magnetic field’s direction.

Parameters:
  • wavelength (ndarray) – Wavelength array

  • tion (ndarray) – Temperature spatial profile

  • bmag (ndarray) – Magnetic field strength spatial profile

  • costheta (ndarray) – Cosine of the angle between the magnetic field direction and the line of sight

Returns:

returns from gauss_ls

Return type:

ndarray

component(wavelength, tion, bmag, parallel=True)#

Returns spectral emission profile parallel or perpendicular to magnetic field’s direction.

Parameters:
  • wavelength (ndarray) – Wavelength array

  • tion (ndarray) – Temperature spatial profile

  • bmag (ndarray) – Magnetic field strength spatial profile

  • parallel (bool) – True for the parallel component and False for the perpendicular component, by default True

Returns:

returns from gauss_ls

Return type:

numpy.ndarray

class cherab.iter.emitter.lineshapes.SimpleLineShape(mass=None, element='Be', linecenter=457.27)#

Simple Doppler-Zeeman triplet (in the case of Paschen-Back effect)

Parameters:
  • mass (float | None) – Element/isotope mass in eV

  • element (str) – Element/isotope symbol, by default “Be”.

  • linecenter (float) – Spectral line center in nm, by default 457.27.

__call__(wavelength, tion, bmag, costheta)#

Returns spectral emission profile with respect to magnetic field’s direction.

Parameters:
  • wavelength (ndarray) – Wavelength array

  • tion (ndarray) – Temperature spatial profile

  • bmag (ndarray) – Magnetic field strength spatial profile

  • costheta (ndarray) – Cosine of the angle between the magnetic field direction and the line of sight

Returns:

returns from gauss_ls

Return type:

ndarray

component(wavelength, tion, bmag, parallel=True)#

Returns spectral emission profile parallel or perpendicular to magnetic field’s direction.

Parameters:
  • wavelength (ndarray) – Wavelength array

  • tion (ndarray) – Temperature spatial profile

  • bmag (ndarray) – Magnetic field strength spatial profile

  • parallel (bool) – True for the parallel component and False for the perpendicular component

Returns:

returns from __call__

Return type:

ndarray

Ray Transfer Matrix Toroidal Emitter#

class cherab.iter.emitter.rtm_toroidal_emitter.RtmRegularRPhiZIntegrator(step=0.001)#

Bases: VolumeIntegrator

Calculates ray transfer matrix (or geometry matrix) for the light sources defined on a regular grid in cylindrical coordinate system: \((R, \phi, Z)\). The value for each light source is stored in respective spectral bin. Custom geometry may be applied to the light sources via the material.map array (see RtmRegularRPhiZEmitter), which can map multiple voxels in \((R, \phi, Z)\) space to a single spectral bin (and thus to a single light source). It is assumed that emitter is periodic in phi direction with a period equal to material.phimax. The path length of the ray passing throuh the voxel is calculated approximately and the accuracy depends on the integration step.

Parameters:

step (float) – Integration step (in meters), by default 0.001.

integrate(spectrum, world, ray, primitive, material, start_point, end_point, world_to_primitive, primitive_to_world)#

Performs a customised integration of the emission through a volume emitter.

Parameters:
  • spectrum (Spectrum) – pectrum measured so far along ray path. Add your emission to this spectrum, don’t override it.

  • world (World) – The world scene-graph.

  • ray (Ray) – The ray being traced.

  • primitive (Primitive) – The geometric primitive to which this material belongs (i.e. a cylinder or a mesh).

  • material (InhomogeneousVolumeEmitter) – The material whose emission needs to be integrated.

  • start_point (Point3D) – The start point for integration in world space.

  • end_point (Point3D) – The end point for integration in world space.

  • world_to_primitive (AffineMatrix3D) – Affine matrix defining the coordinate transform from world space to the primitive’s local space.

  • primitive_to_world (AffineMatrix3D) – Affine matrix defining the coordinate transform from the primitive’s local space to world space.

Returns:

integrated emission

Return type:

Spectrum

class cherab.iter.emitter.rtm_toroidal_emitter.RtmRegularXYZIntegrator(step=0.001)#

Bases: VolumeIntegrator

Calculates ray transfer matrix (or geometry matrix) for the light sources defined on a regular grid in Cartesian coordinate system. The value for each light source is stored in respective spectral bin. Custom geometry may be applied to the light sources via the material.map array (see RtmRegularXYZEmitter), which can map multiple voxels in \((X, Y, Z)\) space to a single spectral bin (and thus to a single light source). The path length of the ray passing throuh the voxel is calculated approximately and the accuracy depends on the integration step.

Parameters:

step (float) – Integration step (in meters), by default 0.001.

integrate(spectrum, world, ray, primitive, material, start_point, end_point, world_to_primitive, primitive_to_world)#

Performs a customised integration of the emission through a volume emitter.

Parameters:
  • spectrum (Spectrum) – pectrum measured so far along ray path. Add your emission to this spectrum, don’t override it.

  • world (World) – The world scene-graph.

  • ray (Ray) – The ray being traced.

  • primitive (Primitive) – The geometric primitive to which this material belongs (i.e. a cylinder or a mesh).

  • material (InhomogeneousVolumeEmitter) – The material whose emission needs to be integrated.

  • start_point (Point3D) – The start point for integration in world space.

  • end_point (Point3D) – The end point for integration in world space.

  • world_to_primitive (AffineMatrix3D) – Affine matrix defining the coordinate transform from world space to the primitive’s local space.

  • primitive_to_world (AffineMatrix3D) – Affine matrix defining the coordinate transform from the primitive’s local space to world space.

Returns:

integrated emission

Return type:

Spectrum

class cherab.iter.emitter.rtm_toroidal_emitter.RtmRegularRPhiZEmitter(nr, nz, dr, dz, rmin, nphi=0, dphi=0.0, phimax=0.0, voxel_map=None, mask=None, integrator=None)#

Bases: InhomogeneousVolumeEmitter

Unity emitter defined on a regular 2D (RZ plane) or 3D (R, phi, Z) grid, which can be used to calculate ray transfer matrix (geometry matrix)

Parameters:
  • nr (int) – Number of grid points in R direction

  • nz (int) – Number of grid points in Z direction

  • dr (float) – Grid step in R direction (in meters)

  • dz (float) – Grid step in Z direction (in meters)

  • rmin (float) – Lower bound of grid in R direction (in meters).

  • nphi (int) – Number of grid points in Phi direction. Default: 0 (2D grid).

  • dphi (float) – Grid step in phi direction (in degree) Used only if nphi > 0.

  • phimax (float) – Upper bound of RPhiZ grid in phi directon (in degree) and a period in toroidal direction. (emiss[:, -1, :] is defined between phimax - dphi and phimax). It is assumed that the grid starts from 0 in phi direction. Just rotate the emitter if not. Used only if nphi > 0.

  • voxel_map (ndarray | None) – 2D (nphi = 0) or 3D (nphi > 0) numpy.ndarray of integers. Array that contain indeces of light sources. It maps voxels in \((R, \phi, Z)\) space to respective light sources. The voxels with identical indeces in voxel_map array form a single light source. If voxel_map[ir, iphi, iz] == -1, the voxel (ir, iphi, iz) will not be mapped to any light source. This parameters allows to apply custom geometry (pixelated though) to the light sources.

  • mask (ndarray | None) – 2D (nphi = 0) or 3D (nphi > 0) numpy.ndarray of booleans. Ignored if voxel_map is provided. The ray tranfer matrix will be calculated only for those voxels for which mask == True.

  • integrator (VolumeIntegrator | None) – integrator object, by default RtmRegularRPhiZIntegrator.

map_from_mask(mask=None)#

masking map array

Parameters:
  • mask (Optional[ndarray]) – mask array, by default None

  • optional – mask array, by default None

Return type:

masked voxel_map

class cherab.iter.emitter.rtm_toroidal_emitter.RtmRegularXYZEmitter(nx, ny, nz, dx, dy, dz, voxel_map=None, mask=None, integrator=None)#

Bases: InhomogeneousVolumeEmitter

Unity emitter defined on a regular 3D \((X, Y, Z)\) grid, which can be used to calculate ray transfer matrices

Parameters:
  • nx (int) – Number of grid points in X direction

  • ny (int) – Number of grid points in Y direction

  • nz (int) – Number of grid points in Z direction

  • dx (float) – Grid step in X direction (in meters)

  • dy (float) – Grid step in Y direction (in meters)

  • dz (float) – Grid step in Z direction (in meters)

  • voxel_map (ndarray | None) – Array that contain indeces of light sources. It maps voxels in \((X, Y, Z)\) space to respective light sources. The voxels with identical indeces in voxel_map array form a single light source. If voxel_map[ix, iy, iz] == -1, the voxel (ix, iy, iz) will not be mapped to any light source. This parameters allows to apply custom geometry (pixelated though) to the light sources.

  • mask (ndarray | None) – Mask array. Ignored if voxel_map is provided. The ray tranfer matrix will be calculated only for those voxels for which mask == True.

  • integrator (VolumeIntegrator) – integrator object, by default RtmRegularXYZIntegrator

map_from_mask(mask=None)#

masking map array

Parameters:
  • mask (Optional[ndarray]) – mask array, by default None

  • optional – mask array, by default None

Return type:

masked voxel_map

Toroidal Emitter#

class cherab.iter.emitter.toroidal_emitter.RegularRPhiZIntegrator(step=0.005)#

Bases: VolumeIntegrator

Integrates the emission along ray’s trajectory. Use with RegularRPhiZEmitter and DiscreteSpectralRPhiZEmitter material.

Parameters:

step (float) – Integration step (in meters), by default 0.005.

integrate(spectrum, world, ray, primitive, material, start_point, end_point, world_to_primitive, primitive_to_world)#

Performs a customised integration of the emission through a volume emitter.

Parameters:
  • spectrum (Spectrum) – pectrum measured so far along ray path. Add your emission to this spectrum, don’t override it.

  • world (World) – The world scene-graph.

  • ray (Ray) – The ray being traced.

  • primitive (Primitive) – The geometric primitive to which this material belongs (i.e. a cylinder or a mesh).

  • material (InhomogeneousVolumeEmitter) – The material whose emission needs to be integrated.

  • start_point (Point3D) – The start point for integration in world space.

  • end_point (Point3D) – The end point for integration in world space.

  • world_to_primitive (AffineMatrix3D) – Affine matrix defining the coordinate transform from world space to the primitive’s local space.

  • primitive_to_world (AffineMatrix3D) – Affine matrix defining the coordinate transform from the primitive’s local space to world space.

Returns:

integrated emission

Return type:

Spectrum

class cherab.iter.emitter.toroidal_emitter.ZeemanRegularRPhiZIntegrator(step=0.005, threshold=1e-05)#

Bases: VolumeIntegrator

Integrates spectral emission along ray’s trajectory in RegularRPhiZEmitter material with respect to magnetic field direction

Parameters:
  • step (float) – Integration step (in meters), by default 0.005.

  • threshold (float) – Zero threshold for the emission, by default 1.e-5.

integrate(spectrum, world, ray, primitive, material, start_point, end_point, world_to_primitive, primitive_to_world)#

Performs a customised integration of the emission through a volume emitter.

Parameters:
  • spectrum (Spectrum) – pectrum measured so far along ray path. Add your emission to this spectrum, don’t override it.

  • world (World) – The world scene-graph.

  • ray (Ray) – The ray being traced.

  • primitive (Primitive) – The geometric primitive to which this material belongs (i.e. a cylinder or a mesh).

  • material (InhomogeneousVolumeEmitter) – The material whose emission needs to be integrated.

  • start_point (Point3D) – The start point for integration in world space.

  • end_point (Point3D) – The end point for integration in world space.

  • world_to_primitive (AffineMatrix3D) – Affine matrix defining the coordinate transform from world space to the primitive’s local space.

  • primitive_to_world (AffineMatrix3D) – Affine matrix defining the coordinate transform from the primitive’s local space to world space.

Returns:

integrated emission

Return type:

Spectrum

class cherab.iter.emitter.toroidal_emitter.RegularRPhiZEmitter(emiss, rmin, dr, dz, phimax=0.0, dphi=0.0, spectralfunc=None, wavelength=None, ti=None, bmag=None, bn=None, integrator=None)#

Bases: InhomogeneousVolumeEmitter

Emitter defined on a regular 3D grid in \((R, \phi, Z)\) coordinates (or 2D grid in \((R, Z)\) coordinates). The emitter is periodic in toroidal direction.

Parameters:
  • emiss (ndarray) – Emissivity:math:(R, phi, Z) or Emissivity:math:(R, Z) defined on a regular 3D or 2D gird.

  • rmin (float) – Lower bound of RPhiZ grid in R directon (in meters) (emiss[0, :, :] is defined between :obj:.`rmin` and rmin + dr).

  • dr (float) – Grid step in R direction (in meters)

  • dz (float) – Grid step in Z direction (in meters)

  • phimax (float) – Upper bound of RPhiZ grid in phi directon (in degree) and a period in toroidal direction. (emiss[:, -1, :] is defined between phimax - dphi and phimax). It is assumed that the grid starts from 0 in phi direction. Just rotate the material if not. Used only in the case of 3D emitter (emiss.ndim > 2).

  • dphi (float) – Grid step in phi direction (in degree) Used only in the case of 3D emitter (emiss.ndim > 2).

  • spectralfunc (Callable[[ndarray, float, float, int], ndarray] | None) – A callable that returns spectral emission profile parallel and perpendicular to magnetic field’s direction spectralfunc(wavelength, ti, bmag, along), where wavelength is array-like, ti is 1D temperature profile, bmag is 1D magnetic field strength profile, and along is 1 for the parallel component and 0 for the perpendicular component

  • ti (ndarray | None) – Ion (or atom) temperature (in eV) defined on a regular 2D gird (RZ-plane). Required when used with ZeemanRegularRPhiZIntegrator.

  • bmag (ndarray | None) – Magnetic field strngth (in T) defined on a regular 2D gird (RZ-plane). Required when used with ZeemanRegularRPhiZIntegrator.

  • bn (ndarray | None) – Magnetic field direction (normalized) defined on a regular 2D gird (RZ-plane). bn[:, :, 0] - R-component, bn[:, :, 1] - phi-component, bn[:, :, 2] - Z-component. bn[:, :, 0]**2 + bn[:, 1]**2 + bn[:, :, 2]**2 = 1 Required when used with ZeemanRegularRPhiZIntegrator.

  • integrator (VolumeIntegrator | None) – Volume integrator: RegularRPhiZIntegrator (defualt) or ZeemanRegularRPhiZIntegrator.

  • wavelength (ndarray | None) –

class cherab.iter.emitter.toroidal_emitter.DiscreteSpectralRPhiZEmitter(emiss, lines_wavelength, min_wavelength, max_wavelength, bins, rmin, dr, dz, phimax=0, dphi=0, integrator=None)#

Bases: InhomogeneousVolumeEmitter

Discrete spectral emitter defined on a regular 3D grid in \((R, \phi, Z)\) coordinates (or 2D grid in \((R, Z)\) coordinates). The emitter is periodic in toroidal direction. This emitter allows to calculate the radiance from multiple spectral lines in a single run of camera.observe() function.

Note

Spectral line shapes are not calculated, only the total intensity of the spectral lines.

Parameters:
  • emiss (ndarray) – Emissivity(R, phi, Z, lines_wavelength) or Emissivity(R, Z, lines_wavelength) defined on a regular 3D or 2D gird for the set of spectral lines.

  • lines_wavelength (ndarray) – Wavelengths of the spectral lines for which the emissivity is provided.

  • min_wavelength (float) – Lower wavelength bound for sampled spectral range. Must be equal to camera.min_wavelength.

  • max_wavelength (float) – Upper wavelength bound for sampled spectral range. Must be equal to camera.max_wavelength.

  • bins (int) – The number of spectral samples over the wavelength range. Must be equal to camera.spectral_bins.

  • rmin (float) – Lower bound of RPhiZ grid in R directon (in meters) (emiss[0, :, :] is defined between rmin and rmin + dr).

  • dr (float) – Grid step in R direction (in meters)

  • dz (float) – Grid step in Z direction (in meters)

  • phimax (float) – Upper bound of RPhiZ grid in phi directon (in degree) and a period in toroidal direction. (emiss[:, -1, :] is defined between phimax - dphi and phimax). It is assumed that the grid starts from 0 in phi direction. Just rotate the material if not. Used only in the case of 3D emitter (emiss.ndim > 2).

  • dphi (float) – Grid step in phi direction (in degree) Used only in the case of 3D emitter (emiss.ndim > 2).

  • integrator (VolumeIntegrator | None) – Volume integrator, by default RegularRPhiZIntegrator.