oimComponent

Components defined in Fourier or image planes

Functions:

getFourierComponents()

A function to get the list of all available components deriving from the oimComponentFourier class

Classes:

oimComponent(**kwargs)

The OImComponent class is the parent abstract class for all types of components that can be added to a OImModel.

oimComponentFourier(**kwargs)

Class for all component analytically defined in the Fourier plan.

oimComponentImage(**kwargs)

Base class for components define in 2D : x,y (regular grid) in the image plan

oimComponentRadialProfile(**kwargs)

Base class for components define by their radial profile

oimComponentFitsImage([fitsImage, useinternalPA])

Component load load images or chromatic-cubes from fits files

oimodeler.oimComponent.getFourierComponents()

A function to get the list of all available components deriving from the oimComponentFourier class

Returns:

res – list of all available components deriving from the oimComponentFourier class.

Return type:

list

class oimodeler.oimComponent.oimComponent(**kwargs)

The OImComponent class is the parent abstract class for all types of components that can be added to a OImModel.

It has a similar interface than the oimModel and allow to compute images (or image cubes fore wavelength dependent models or time dependent models) and complex coherentFluxes for a vector of u,v,wl, and t coordinates

Variables

name:

The name of the component

shortname:

Short name for the component

description:

Detailed description of the component

params:

The dictionary of the component parameters

_firstInit = True
name = 'Generic component'
shortname = 'Gen comp'
description = 'This is the class from which all components derived'
_paramstr()
property _wl: ndarray

Gets the wavelengths.

property _t: ndarray

Gets the times.

_eval(checkParam=True, **kwargs) None
getComplexCoherentFlux(u, v, wl=None, t=None) ndarray

Compute and return the complex coherent flux for an array of u,v (and optionally wavelength and time ) coordinates

Parameters:
  • u (list or numpy array) – spatial coordinate u (in cycles/rad)

  • v (list or numpy array) – spatial coordinate vu (in cycles/rad) .

  • wl (list or numpy array, optional) – wavelength(s) in meter. The default is None.

  • t (list or numpy array, optional) – time in s (mjd). The default is None.

Returns:

The complex coherent flux.

Return type:

A numpy array of the same size as u & v

getImage(dim, pixSize, wl=None, t=None) ndarray

Compute and return an image or and image cube (if wavelength and time are given).

The returned image as the x,y dimension dim in pixel with an angular pixel size pixSize in rad. Image is returned as a numpy array unless the keyword fits is set to True. In that case the image is returned as an astropy.io.fits hdu

Parameters:
  • dim (integer) – image x & y dimension in pixels

  • pixSize (float) – pixel angular size in rad

  • wl (integer, list or numpy array, optional) – wavelength(s) in meter. The default is None

  • t (integer, list or numpy array, optional) – time in s (mjd). The default is None

  • fits (bool, optional) – if True returns result as a fits hdu. The default is False

Returns:

  • a numpy 2D array (or 3 or 4D array if wl, t or both are given) or an

  • astropy.io.fits hdu. image hdu if fits=True. – The image of the component with given size in pixels and rad

_ftTranslateFactor(ucoord, vcoord, wl, t) ndarray
_directTranslate(x, y, wl, t)
getNonRegularImage(xx, yy, wl=None, t=None)

Compute and return a non-regular image function at the xx, yy and optional wl and t coordinates)

serialize() Dict[str, Any]

Serializes the oimComponent.

classmethod deserialize(ser: Dict[str, Any]) oimComponent

Deserializes into an oimComponent.

pickle(f: str | Path | BufferedWriter, **kwargs) None

Save the pickled representation of the object into an already open file or opens a file from a string or pathlib.Path.

classmethod unpickle(f: Path | TextIOWrapper, **kwargs) object

Read the pickled representation from an open file or reads a string or pathlib.Path into a file to return the reconstituted object.

class oimodeler.oimComponent.oimComponentFourier(**kwargs)

Class for all component analytically defined in the Fourier plan. Inherit from the oimComponent.

Implements translation in direct and Fourier space, getImage from the Fourier definition of the object, ellipticity (i.e. flatening) Children classes should only implement the _visFunction and _imageFunction functions.

elliptic = False
flat = False
extincted = False
getComplexCoherentFlux(ucoord, vcoord, wl=None, t=None)

Compute and return the complex coherent flux for an array of u,v (and optionally wavelength and time ) coordinates

Parameters:
  • u (list or numpy array) – spatial coordinate u (in cycles/rad)

  • v (list or numpy array) – spatial coordinate vu (in cycles/rad) .

  • wl (list or numpy array, optional) – wavelength(s) in meter. The default is None.

  • t (list or numpy array, optional) – time in s (mjd). The default is None.

Returns:

The complex coherent flux.

Return type:

A numpy array of the same size as u & v

_visFunction(ucoord, vcoord, rho, wl, t)
getImage(dim, pixSize, wl=None, t=None)

Compute and return an image or and image cube (if wavelength and time are given).

The returned image as the x,y dimension dim in pixel with an angular pixel size pixSize in rad. Image is returned as a numpy array unless the keyword fits is set to True. In that case the image is returned as an astropy.io.fits hdu

Parameters:
  • dim (integer) – image x & y dimension in pixels

  • pixSize (float) – pixel angular size in rad

  • wl (integer, list or numpy array, optional) – wavelength(s) in meter. The default is None

  • t (integer, list or numpy array, optional) – time in s (mjd). The default is None

  • fits (bool, optional) – if True returns result as a fits hdu. The default is False

Returns:

  • a numpy 2D array (or 3 or 4D array if wl, t or both are given) or an

  • astropy.io.fits hdu. image hdu if fits=True. – The image of the component with given size in pixels and rad

getNonRegularImage(xx, yy, wl=None, t=None)

Compute and return a non-regular image function at the xx, yy and optional wl and t coordinates)

_imageFunction(xx, yy, wl, t)
class oimodeler.oimComponent.oimComponentImage(**kwargs)

Base class for components define in 2D : x,y (regular grid) in the image plan

elliptic = False
flat = False
extincted = False
getComplexCoherentFlux(ucoord, vcoord, wl=None, t=None)

Compute and return the complex coherent flux for an array of u,v (and optionally wavelength and time ) coordinates

Parameters:
  • u (list or numpy array) – spatial coordinate u (in cycles/rad)

  • v (list or numpy array) – spatial coordinate vu (in cycles/rad) .

  • wl (list or numpy array, optional) – wavelength(s) in meter. The default is None.

  • t (list or numpy array, optional) – time in s (mjd). The default is None.

Returns:

The complex coherent flux.

Return type:

A numpy array of the same size as u & v

getImage(dim, pixSize, wl=None, t=None)

Compute and return an image or and image cube (if wavelength and time are given).

The returned image as the x,y dimension dim in pixel with an angular pixel size pixSize in rad. Image is returned as a numpy array unless the keyword fits is set to True. In that case the image is returned as an astropy.io.fits hdu

Parameters:
  • dim (integer) – image x & y dimension in pixels

  • pixSize (float) – pixel angular size in rad

  • wl (integer, list or numpy array, optional) – wavelength(s) in meter. The default is None

  • t (integer, list or numpy array, optional) – time in s (mjd). The default is None

  • fits (bool, optional) – if True returns result as a fits hdu. The default is False

Returns:

  • a numpy 2D array (or 3 or 4D array if wl, t or both are given) or an

  • astropy.io.fits hdu. image hdu if fits=True. – The image of the component with given size in pixels and rad

getInternalImage(wl, t)
_internalImage()
_imageFunction(xx, yy, wl, t)
_getInternalGrid(simple=True, flatten=False, wl=None, t=None)
getPixelSize(mas=False)
class oimodeler.oimComponent.oimComponentRadialProfile(**kwargs)

Base class for components define by their radial profile

asymmetric = False
modulation_order = 1
elliptic = False
flat = False
extincted = False
_getInternalGrid(simple=True, flatten=False, wl=None, t=None)
_internalRadialProfile()
_radialProfileFunction(r=None, wl=None, t=None)
getInternalRadialProfile(wl, t)
static hankel(Ir, r, wlin, tin, sfreq, wl, t, precision=None)
getImage(dim, pixSize, wl=None, t=None)

Compute and return an image or and image cube (if wavelength and time are given).

The returned image as the x,y dimension dim in pixel with an angular pixel size pixSize in rad. Image is returned as a numpy array unless the keyword fits is set to True. In that case the image is returned as an astropy.io.fits hdu

Parameters:
  • dim (integer) – image x & y dimension in pixels

  • pixSize (float) – pixel angular size in rad

  • wl (integer, list or numpy array, optional) – wavelength(s) in meter. The default is None

  • t (integer, list or numpy array, optional) – time in s (mjd). The default is None

  • fits (bool, optional) – if True returns result as a fits hdu. The default is False

Returns:

  • a numpy 2D array (or 3 or 4D array if wl, t or both are given) or an

  • astropy.io.fits hdu. image hdu if fits=True. – The image of the component with given size in pixels and rad

getComplexCoherentFlux(ucoord, vcoord, wl=None, t=None)

Compute and return the complex coherent flux for an array of u,v (and optionally wavelength and time ) coordinates

Parameters:
  • u (list or numpy array) – spatial coordinate u (in cycles/rad)

  • v (list or numpy array) – spatial coordinate vu (in cycles/rad) .

  • wl (list or numpy array, optional) – wavelength(s) in meter. The default is None.

  • t (list or numpy array, optional) – time in s (mjd). The default is None.

Returns:

The complex coherent flux.

Return type:

A numpy array of the same size as u & v

class oimodeler.oimComponent.oimComponentFitsImage(fitsImage=None, useinternalPA=False, **kwargs)

Component load load images or chromatic-cubes from fits files

elliptic = False
name = 'Fits Image Component'
shortname = 'Fits_Comp'
loadImage(fitsImage, useinternalPA=False)
_internalImage()
getPixelSize(mas=False)