Doppler Maps

class starry.DopplerMap(ydeg=1, udeg=0, nc=1, nt=10, wav=None, wav0=None, wavc=None, oversample=None, interpolate=True, interp_order=1, interp_tol=1e-12, vsini_max=None, lazy=None, **kwargs)

Class representing a spectral-spatial stellar surface.

Instances of this class can be used to model the spatial and wavelength dependence of stellar surfaces and the spectral timeseries one would observe in the presence of rotation-induced Doppler shifts. This class is specifically suited to Doppler imaging studies. It is similar (but different in important ways; see below) to the traditional starry.Map class.

Important notes:

  • Currently, this class does not model occultations.

Parameters
  • ydeg (int, optional) – Degree of the spherical harmonic map. Default is 1.

  • udeg (int, optional) – Degree of the limb darkening filter. Default is 0.

  • nc (int, optional) – Number of spectral components. The surface is represented as the outer product of nc spherical harmonic coefficient vectors and nc associated spectra. Default is 1.

  • nt (int, optional) – Number of epochs. This is the number of spectra one wishes to model. Default is 10.

  • wav (ndarray, optional) – Wavelength grid on which the data is defined (and on which the model will be computed) in nanometers. Default is a grid of 200 bins uniformly spaced between 642.5 nm and 643.5 nm, the region spanning the classical FeI 6430 line used in Doppler imaging.

  • wav0 (ndarray, optional) – Wavelength grid on which the rest frame spectrum (the “template”) is defined. This can be the same as wav, but it is strongly recommended to pad it on either side to avoid edge effects. If wav0 is insufficiently padded, the values of the model near the edges of the wavelength grid will be incorrect, since they include contribution from Doppler-shifted spectral features beyond the edge of the grid. By default, if wav0 is not provided or set to None, starry will compute the required amount of padding automatically. See vsini_max below for more info.

  • wavc (scalar, optional) – A wavelength at which there are no spectral features in the observed spectrum. This is used for normalizing the model at each epoch. Default is the first element of wav.

  • oversample (float, optional) – The oversampling factor for the internal wavelength grid. The number of spectral bins used to perform the convolutions is equal to this value times the length of wav. If both, wav0 and wav are provided, default is the ratio of the lengths of the two arrays. Otherwise, default is 1.

  • interpolate (bool, optional) – The wavelength grid used internally is different from wav, since starry requires a grid that is uniformly spaced in the log of the wavelength. By default, interpolate is True, in which case starry interpolates to the user-defined wavelength grids (wav and wav0) when any method or property is accessed. This incurs a slight slowdown in the computation, so users can obtain higher efficiency by setting this to False and handling the interpolation on their own. If False, all spectra returned by the methods and properties of this class will instead be defined on the wav_ and wav0_ grids.

  • interp_order (int, optional) – Order of the spline interpolation between the internal (wav_ and wav0_) and user-facing (wav and wav0) wavelength grids. Default is 1 (linear).

  • interp_tol (float, optional) – For splines with interp_order greater than one, the interpolation can be made much faster by zeroing out elements of the interpolation operator whose absolute value is smaller than some small tolerance given by interp_tol. Default is 1e-12.

  • vsini_max (float, optional) – Maximum value of vsini for this map. This sets the size of the convolution kernel (and the amount of padding required in wav0) and should be adjusted to ensure that map.veq * sin(map.inc) is never larger than this quantity. Lower values of this quantity will result in faster evaluation times. Default is 100 km/s.

  • angle_unit (astropy.units.Unit, optional) – The unit used for angular quantities. Default deg.

  • velocity_unit (astropy.units.Unit, optional) – The unit used for velocity quantities. Default m/s.

  • inc (scalar, optional) – Inclination of the star in units of angle_unit. Default is 90.0.

  • obl (scalar, optional) – Obliquity of the star in units of angle_unit. Default is 0.0.

  • veq (scalar, optional) – Equatorial rotational velocity of the star in units of velocity_unit. Default is 0.0.

property N

Total number of map coefficients. Read-only

This is equal to \(N_\mathrm{y} + N_\mathrm{u}\).

property Nu

Number of limb darkening coefficients, including \(u_0\). Read-only

This is equal to \(u_\mathrm{deg} + 1\).

property Ny

Number of spherical harmonic coefficients. Read-only

This is equal to \((y_\mathrm{deg} + 1)^2\).

property angle_unit

An astropy.units unit defining the angle metric for this map.

baseline(theta=None, full=False)

Return the photometric baseline at each epoch.

Parameters
  • theta (vector, optional) – The angular phase(s) at which to compute the design matrix, in units of angle_unit. This must be a vector of size nt. Default is uniformly spaced values in the range [0, 2 * pi).

  • full (bool, optional) – If True, returns a matrix with the same shape as flux(). If False (default), returns a vector of length nc.

property deg

Total degree of the map. Read-only

This is equal to \(y_\mathrm{deg} + u_\mathrm{deg}\).

design_matrix(theta=None, fix_spectrum=False, fix_map=False)

Return the Doppler imaging design matrix.

This matrix dots into the spectral map to yield the model for the observed spectral timeseries (the flux).

Note that if this method is used to compute the spectral timeseries, the result should be reshaped into a matrix of shape (nt, nw) and optionally divided by the baseline() to match the return value of flux().

Parameters
  • theta (vector, optional) – The angular phase(s) at which to compute the design matrix, in units of angle_unit. This must be a vector of size nt. Default is uniformly spaced values in the range [0, 2 * pi).

  • fix_spectrum (bool, optional) – If True, returns the design matrix for a fixed spectrum; this can then be dotted into the spherical harmonic coefficient matrix to obtain the flux. See below for details. Default is False.

  • fix_map (bool, optional) – If True, returns the design matrix for a fixed map; this can then be dotted into the spectrum matrix to obtain the flux. See below for details. Default is False.

If fix_spectrum and fix_map are False (default), this method returns a sparse matrix of shape (nt * nw, Ny * nw0_). The flux may be computed from (assuming lazy is False)

If fix_spectrum is True, returns a dense matrix of shape (nt * nw, nc * Ny). The flux may be computed from

Finally, if fix_map is True, returns a dense matrix of shape (nt * nw, nc * nw0_). The flux may be computed from

Note

Instantiating this matrix is usually a bad idea, since it can be very slow and consume a lot of memory. Check out the dot() method instead for fast dot products with the design matrix.

dot(x, theta=None, transpose=False, fix_spectrum=False, fix_map=False)

Dot the Doppler design matrix into a given matrix or vector.

This method is useful for computing dot products between the design matrix and the spectral map (to compute the model for the spectrum) or between the design matrix and (say) a covariance matrix (when doing inference). This is in general much, much faster than instantiating the design_matrix() explicitly and dotting it in.

Parameters
  • x (vector or matrix) – The column vector or matrix into which the design matrix is dotted. This argument must have a specific shape that depends on the other arguments to this method; see below.

  • theta (vector, optional) – The angular phase(s) at which to compute the design matrix, in units of angle_unit. This must be a vector of size nt. Default is uniformly spaced values in the range [0, 2 * pi).

  • transpose (bool, optional) – If True, dots the transpose of the design matrix into x. Default is False.

  • fix_spectrum (bool, optional) – If True, performs the operation using the design matrix for a fixed spectrum. The current spectrum is “baked into” the design matrix, so the tensor x is a representation of the spherical harmonic decomposition of the surface. Default is False.

  • fix_map (bool, optional) – If True, performs the oprtation using the the design matrix for a fixed map. The current spherical harmonic map is “baked into” the design matrix, so the tensor x is a representation of the spectral decomposition of the surface. Default is False.

The shapes of x and of the matrix returned by this method depend on the method’s arguments. If transpose is False, this returns a dense matrix of shape (nt * nw, ...), where ... are any additional dimensions in x beyond the first. In this case, the shape of x should be as follows:

  • If fix_spectrum and fix_map are False (default), the input argument x must have shape (Ny * nw0_, ...).

  • If fix_spectrum is True, the input argument x must have shape (nc * Ny, ...).

  • If fix_map is True, the input argument x must have shape (nc * nw0_, ...).

If, instead, transpose is True, this returns a dense matrix of a shape that depends on the arguments:

  • If fix_spectrum and fix_map are False (default), the return value has shape (Ny * nw0_, ...).

  • If fix_spectrum is True, the return value has shape (nc * Ny, ...).

  • If fix_map is True, the return value has shape (nc * nw0_, ...).

When transpose is True, the input argument x must always have shape (nt * nw, ...).

Note that if this method is used to compute the spectral timeseries (with tranpose = False), the result should be reshaped into a matrix of shape (nt, nw) to match the return value of flux() and optionally divided by the baseline() to match the return value of flux().

flux(theta=None, normalize=True, method='dotconv')

Return the model for the full spectral timeseries.

Parameters
  • theta (vector, optional) – The angular phase(s) at which to compute the design matrix, in units of angle_unit. This must be a vector of size nt. Default is uniformly spaced values in the range [0, 2 * pi).

  • normalize (bool, optional) – Whether or not to normalize the flux. If True (default), normalizes the flux so that the continuum level is unity at all epochs. If False, preserves the natural changes in the continuum as the total flux of the star changes during its rotation. Note that this is not usually an observable effect, since spectrographs aren’t designed to preserve this information!

  • method (str, optional) – The strategy for computing the flux. Must be one of dotconv, convdot, conv, or design. Default is dotconv, which is the fastest method in most cases. All three of dotconv, convdot, and conv compute the flux via fast two-dimensional convolutions (via conv2d), while design computes the flux by instantiating the design matrix and dotting it in. This last method is usually extremely slow and memory intensive; its use is not recommended in general.

This method returns a matrix of shape (nt, nw) corresponding to the model for the observed spectrum (evaluated on the wavelength grid wav) at each of nt epochs.

property inc

The inclination of the rotation axis in units of angle_unit.

intensity(n=0, **kwargs)

Return the intensity of a surface map component. See the documentation for starry.Map().intensity() for details.

Parameters

n (int, optional) – The index of the surface map component. Default is 0.

property lazy

Map evaluation mode – lazy or greedy?

load(*, map=None, maps=None, spectrum=None, spectra=None, cube=None, smoothing=None, fac=1.0, eps=1e-12)

Load a spatial and/or spectral representation of a stellar surface.

Users may load spatial maps and/or spectra: one of each is required per map component (nc). Alternatively, users may provide a single data cube containing the full spectral-spatial representation of the surface. This cube should have shape (nlat, nlon, nw0), where nlat is the number of pixels in latitude, nlon is the number of pixels in longitude, and nw0 is the number of wavelength bins in the rest frame spectrum.

This routine performs a simple spherical harmonic transform (SHT) to compute the spherical harmonic expansion given maps or a data cube. If a data cube is provided, this routine performs singular value decomposition (SVD) to compute the nc component surface maps and spectra (the “eigen” components) that best approximate the input.

Parameters
  • map (str or ndarray, optional) – A list or ndarray of surface maps on a rectangular latitude-longitude grid. This may also be a list of strings corresponding to the names of PNG image files representing the surface maps.

  • maps (str or ndarray, optional) – Alias for map.

  • spectrum (ndarray, optional) – A list or ndarray of vectors containing the spectra corresponding to each map component.

  • spectra (ndarray, optional) – Aliass for spectrum.

  • cube (ndarray, optional) – A 3-dimensional ndarray of shape (nlat, nlon, nw0) containing the spectra at each position on a latitude-longitude grid spanning the entire surface.

  • smoothing (float, optional) – Gaussian smoothing strength. Increase this value to suppress ringing or explicitly set to zero to disable smoothing. Default is 2/self.ydeg.

  • fac (float, optional) – Factor by which to oversample the image when applying the SHT. Default is 1.0. Increase this number for higher fidelity (at the expense of increased computational time).

  • eps (float, optional) – Regularization strength for the spherical harmonic transform. Default is 1e-12.

property nc

Number of spectro-spatial map components. Read-only

property nt

Number of spectral epochs. Read-only

property nw

Length of the user-facing flux wavelength grid wav. Read-only

property nw0

Length of the user-facing rest frame spectrum wavelength grid wav0. Read-only

property nw0_

Length of the internal rest frame spectrum wavelength grid wav0_. Read-only

property nw_

Length of the internal flux wavelength grid wav_. Read-only

property obl

The obliquity of the rotation axis in units of angle_unit.

optimize(model=None, start=None, niter=10000, lr=0.0001, quiet=False, **kwargs)

Optimizes the log likelihood function within a pymc3 model context using the Adam solver.

Parameters
  • model (optional) – The pymc3 model to optimize. Default is None, in which case the current model context is used.

  • start (optional) – The starting point for the optimization. Default is model.test_point.

  • niter (int, optional) – The number of iterations. Default is 10000.

  • lr (float, optional) – The learning rate. Default is 1e-4.

  • quiet (bool, optional) – If True, disables the progress bar. Default is False.

  • kwargs – Additional kwargs passed directly to pymc3_ext.optim.Adam.

Below is an example of how to use this method.

import pymc3 as pm
import starry

# Define a `pymc3` model
with pm.Model() as model:

    # Instantiate a Doppler map
    # Assuming ``nt=16`` spectral epochs
    map = starry.DopplerMap(ydeg=15, nt=16)

    # SHT matrix: converts from pixels to Ylms
    A = map.sht_matrix(smoothing=0.075)
    npix = A.shape[1]

    # Prior on the map: intensity uniform in [0, 1]
    p = pm.Uniform("pixels", lower=0.0, upper=1.0, shape=(npix,))
    amp = pm.Uniform("amp", lower=0.0, upper=1.0)
    map[:, :] = amp * tt.dot(A, p)

    # Prior on the spectrum: Gaussian about unity
    map.spectrum = pm.Normal(
        "spectrum", mu=1.0, sigma=1e-1, shape=(map.nw0,)
    )

    # Compute the model
    # Assuming `theta` is a vector of 16 rotational phases at
    # which each of the spectra were observed
    flux_model = map.flux(theta=theta)

    # Likelihood term
    # Assuming `flux` is the data and `flux_err` is the
    # (scalar) uncertainty
    pm.Normal(
        "obs",
        mu=tt.reshape(flux_model, (-1,)),
        sd=flux_err,
        observed=flux.reshape(-1,)
    )

    # Optimize (this function)
    map_soln, loss = map.optimize()
property oversample

Spectrum oversampling factor. Read-only

render(n=0, **kwargs)

Render a component surface map. See the documentation for starry.Map().render() for details.

Parameters

n (int, optional) – The index of the surface map component. Default is 0.

reset(**kwargs)

Reset all map coefficients and attributes.

show(n=0, **kwargs)

Show a component surface map. See the documentation for starry.Map().show() for details.

Parameters

n (int, optional) – The index of the surface map component. Default is 0.

sht_matrix(inverse=False, return_grid=False, smoothing=None, oversample=2, lam=1e-06)

Return the Spherical Harmonic Transform (SHT) matrix.

This matrix dots into a vector of pixel intensities defined on a set of latitude-longitude points, resulting in a vector of spherical harmonic coefficients that best approximates the image.

This can be useful for transforming between the spherical harmonic and pixel representations of the surface map, such as when specifying priors during inference or optimization. A common prior in Doppler imaging problems is the assumption of sparsity, often in the form of maximum entropy regularization, where the likelihood penalty has the form (e.g., Vogt et al. 1987)

In this case, the prior must be applied to the vector p of pixel intensities, not the vector of spherical harmonic coefficients. The former is obtained from the latter as follows:

(Note, importantly, that p is not guaranteed to be positive, so one should be careful when taking the log in the expression above!)

An alternative way to approach the problem is to sample in pixels p directly; in this case, one must compute y from p so that starry can compute the model for the flux:

Parameters
  • inverse (bool, optional) – which transforms from spherical harmonic coefficients to pixels. Default is False.

  • return_grid (bool, optional) – shape (npix, 2) corresponding to the latitude-longitude points (in units of angle_unit) on which the SHT is evaluated. Default is False.

  • smoothing (float, optional) – Gaussian smoothing strength. Increase this value to suppress ringing (forward SHT only) or explicitly set to zero to disable smoothing. Default is 2/self.ydeg.

  • oversample (int, optional) – Factor by which to oversample the pixelization grid. Default 2.

  • lam (float, optional) – Regularization parameter for the inverse pixel transform. Default 1e-6.

Returns

A matrix of shape (Ny, npix) or (npix, Ny) (if inverse is True), where npix is the number of pixels on the grid. This number is determined from the degree of the map and the oversample keyword. If return_grid is True, also returns the latitude-longitude points (in units of angle_unit) on which the SHT is evaluated, a matrix of shape (npix, 2).

solve(flux, solver='bilinear', **kwargs)

Iteratively solves the bilinear or nonlinear problem for the spatial and/or spectral map given a spectral timeseries.

Parameters
  • flux (matrix) – The observed spectral timeseries. Must be an ndarray of shape (py:attr:nt, py:attr:nw).

  • solver (str, optional) – Which solver to use. Options are “bilinear” or “nonlinear”. Default is “bilinear”.

  • flux_err (float, vector, or matrix, optional) – The data uncertainty. If a scalar, the data is assumed to be homoscedastic. If a vector, the data for each epoch is assumed to be homoscedastic; in this case, must have length py:attr:nt. If a matrix, must have shape (py:attr:nt, py:attr:nw). Note that correlated errors (i.e., non-diagonal data covariances) are not currently supported. Please reach out if this is something you’d like to see implemented. Default is 1e-4.

  • theta (vector, optional) – The angular phase(s) at which the spectra were observed, in units of angle_unit. This must be a vector of size nt. Default is uniformly spaced values in the range [0, 2 * pi).

  • spatial_mean (float, vector, or list, optional) – The prior mean on the spherical harmonic coefficients of the map. If a scalar, the same mean is assumed for all coefficients. If a vector of length Ny, the same mean vector is assumed for all map components. Users can also provide a list of length nc containing scalars or vectors corresponding to the prior mean for each component. Default is the vector [1.0, 0.0, …, 0.0].

  • spatial_cov (float, vector, matrix, or list, optional) – The prior (co)variance on the spherical harmonic coefficients of the map. If a scalar, assumes the same variance for all map coefficients (with no covariance across them). If a vector of length Ny, sets the prior variance for each spherical harmonic coefficient individually and assumes the same variance vector for all map components. If a matrix of shape (Ny, Ny), sets the prior covariance of each map component equal to this matrix. Users can also provide a list of length nc containing scalars, vectors, or matrices corresponding to the prior covariance for each component. Default is 1e-4.

  • spectral_mean (float, vector, or list, optional) – The prior mean on the spectral components of the map. If a scalar, the same mean is assumed for all spectral elements. If a vector of length nw0, the same mean vector is assumed for all spectral components. Users can also provide a list of length nc containing scalars or vectors corresponding to the prior mean for each component. Default is 1.0.

  • spectral_cov (float, vector, matrix, or list, optional) – The prior (co)variance on the spectra. If a scalar, assumes the same variance for all spectral elements (with no covariance across them). If a vector of length nw0, sets the prior variance for each spectral element individually and assumes the same variance vector for all spectral components. If a matrix of shape (nw0, nw0), sets the prior covariance of each spectral component equal to this matrix. Users can also provide a list of length nc containing scalars, vectors, or matrices corresponding to the prior covariance for each component. Default is 1e-3.

  • spectral_guess (float, vector, or list, optional) – The guess for the spectral components of the map. If a scalar, the same value is assumed for all spectral elements. If a vector of length nw0, the same guess vector is assumed for all spectral components. Users can also provide a list of length nc containing scalars or vectors corresponding to the guess for each component. Default is to compute the guess based on a deconvolution of the mean observed spectrum.

  • spectral_lambda (float, optional) – The regularization parameter for the L1 solver. Increasing this value increases the sparsity of the solution. Default is 1e5.

  • spectral_maxiter (int, optional) – Maximum number of iterations in the L1 solver. Default is 100.

  • spectral_eps (float, optional) – Small parameter added to the diagonal of the spectral covariance matrix for stability. Default is 1e-12.

  • spectral_tol (float, optional) – Tolerance for termination of the L1 iterative solver. Default is 1e-8.

  • spectral_method (str, optional) – Regularization method when solving for the spectrum. Options are “L1” or “L2” (default). When solving for both the spectrum and the map, the L1 solver is used to obtain an initial guess for the spectrum, regardless of this setting.

  • normalized (bool, optional) – Whether the flux dataset is continuum-normalized. Default is True. If it is normalized, the solution for the map is non-linear, but typically converges with an iterative tempered solver. See the nlogT tempering parameter below.

  • baseline (vector, optional) – If normalized is True, users may provide a vector of length nt corresponding to the photometric baseline at each epoch. This is used to un-normalize the data, restoring the linearity of the problem. If this quantity is not known exactly (as is almost always the case), do not provide a value for it; instead, set the baseline_var parameter, which controls the variance of the prior on the baseline; the code will then marginalize over it. Default is None.

  • baseline_var (float, optional) – Variance of the prior on the baseline, when baseline is not provided and the data is normalized. Default is 1e-2.

  • fix_spectrum (bool, optional) – If True, fixes the spectrum at the current value and solves only for the map. Default is False.

  • fix_map (bool, optional) – If True, fixes the map at the current value and solves only for the spectrum. Default is False.

  • logT0 (float, optional) – Initial log temperature for the tempering scheme. Default is 12.0. See the nlogT tempering parameter below for more details.

  • logTf (float, optional) – Final log temperature for the tempering scheme. Default is 0.0 (untempered). See the nlogT tempering parameter below for more details.

  • nlogT (int, optional) – Number of steps in the tempering scheme. Default is 50. The tempering scheme is used to solve for the map when the baseline is not known, or to solve for both the map and the specturm when neither are known. At each step, the square of the flux uncertainty (i.e., the variance) is multiplied by the current value of the “temperature”, which is slowly decreased from an initial large number controlled by logT0 to a small number, controlled by logTf. In practice, the correct tempering schedule can greatly improve the odds that the solver converges to the global minimum. Different problems in general require different settings for these parameters, so fine tuning is usually required.

  • quiet (bool, optional) – Suppress messages and progress bars? Default is False.

property spectral_map

The spectral-spatial map vector. Read only

This is equal to the (unrolled) outer product of the spherical harmonic decompositions and their corresponding spectral components. Dot the design matrix into this quantity to obtain the observed spectral timeseries (the flux()).

property spectrum

The rest frame spectrum for each component.

This quantity is defined on the wavelength grid wav0. Shape must be (nc, nw0). If nc is unity, a one-dimensional array of length nw0 is also accepted.

property spectrum_

The internal rest frame spectrum for each component. Read only

This is defined on the wavelength grid wav0_.

property u

The vector of limb darkening coefficients. Read-only

To set this vector, index the map directly using one index: map[n] = ... where n is the degree of the limb darkening coefficient. This may be an integer or an array of integers. Slice notation may also be used.

property udeg

Degree of the limb darkening applied to the star. Read-only

property velocity_unit

An astropy.units unit defining the velocity metric for this map.

property veq

The equatorial velocity of the body in units of velocity_unit.

visualize(backend='bokeh', **kwargs)

Display or save a visualization of the star with optional interactivity.

If backend is set to bokeh, this method uses the bokeh package to render an interactive visualization of the spectro-spatial stellar surface and the model for the spectral timeseries. The output is an HTML page that is either saved to disk (if file is provided) or displayed in a browser window or inline (if calling this method from within a Jupyter notebook).

Users can interact with the visualization by moving the mouse over the map to show the emergent, rest frame spectrum at different points on the surface. Users can also scroll (with the mouse wheel or track pad) to change the wavelength at which the map is visualized (in the left panel) or to rotate the orthographic projection of the map (in the right panel).

If instead backend is set to matplotlib, this method shows static plots of the component surface maps and/or spectra.

Parameters
  • backend (str, optional) – The visualization backend. Options are bokeh or matplotlib. Default is bokeh.

  • theta (vector, optional) – The angular phase(s) at which to compute the design matrix, in units of angle_unit. If backend is bokeh, this must be a vector of size nt. If backend` is ``matplotlib, this may be a scalar. Default is uniformly spaced values in the range [0, 2 * pi) (bokeh) or 0.0 (matplotlib).

  • res (int, optional) – Resolution of the map image in pixels on a side. Default is 150.

  • file (str, optional) – Path to an HTML file (bokeh backend) or PDF, JPG, PNG, etc. figure file (matplotlib backend) to which the visualization will be saved. Default is None, in which case the visualization is displayed.

If the backend is set to matplotlib, additional keyword arguments are allowed:

Parameters
  • show_maps (bool, optional) – Show the component surface maps? Default is True.

  • show_spectra (bool, optional) – Show the component spectra? Default is True.

  • projection (str, optional) – Cartographic projection to plot. Options are rect (equirectangular), moll (Mollweide), and ortho (orthographic). Default is moll.

  • vmin (float, optional) – Minimum value in the color scale. Default is None.

  • vmax (float, optional) – Maximum value in the color scale. Default is None.

  • nmax (int, optional) – Maximum number of components to show. Default is 5.

Any other keywords accepted by the show method are also allowed.

Note

The bokeh visualization can be somewhat memory-intensive! Try decreasing the map resolution or switch to the matplotlib backend if you experience issues.

property vsini

The projected equatorial radial velocity in units of velocity_unit. Read-only

property wav

The wavelength grid for the spectral timeseries model. Read-only

This is the wavelength grid on which quantities like the flux() and design_matrix() are defined.

property wav0

The rest-frame wavelength grid. Read-only

This is the wavelength grid on which the spectrum is defined.

property wav0_

The internal rest frame spectrum wavelength grid. Read-only

This is the wavelength grid on which the spectrum_ is defined.

property wav_

The internal model wavelength grid. Read-only

property wavc

A wavelength corresponding to the continuum level. Read only

This is specified when instantiating a DopplerMap and should be a wavelength at which there are no spectral features in the observed spectrum. This is used for normalizing the model at each epoch.

property y

The spherical harmonic coefficient matrix. Read-only

property ydeg

Spherical harmonic degree of the map. Read-only