astrobase.lcfit.utils module¶
This contains utilities for fitting routines in the rest of this subpackage.
-
astrobase.lcfit.utils.
get_phased_quantities
(stimes, smags, serrs, period)[source]¶ Does phase-folding for the mag/flux time-series given a period.
Given finite and sigma-clipped times, magnitudes, and errors, along with the period at which to phase-fold the data, perform the phase-folding and return the phase-folded values.
Parameters: - stimes,smags,serrs (np.array) – The sigma-clipped and finite input mag/flux time-series arrays to operate on.
- period (float) – The period to phase the mag/flux time-series at. stimes.min() is used as the epoch value to fold the times-series around.
Returns: (phase, pmags, perrs, ptimes, mintime) – The tuple returned contains the following items:
- phase: phase-sorted values of phase at each of stimes
- pmags: phase-sorted magnitudes at each phase
- perrs: phase-sorted errors
- ptimes: phase-sorted times
- mintime: earliest time in stimes.
Return type: tuple
-
astrobase.lcfit.utils.
make_fit_plot
(phase, pmags, perrs, fitmags, period, mintime, magseriesepoch, plotfit, magsarefluxes=False, wrap=False, model_over_lc=True, fitphase=None)[source]¶ This makes a plot of the LC model fit.
Parameters: - phase,pmags,perrs (np.array) – The actual mag/flux time-series.
- fitmags (np.array) – The model fit time-series.
- period (float) – The period at which the phased LC was generated.
- mintime (float) – The minimum time value.
- magseriesepoch (float) – The value of time around which the phased LC was folded.
- plotfit (str) – The name of a file to write the plot to.
- magsarefluxes (bool) – Set this to True if the values in pmags and fitmags are actually fluxes.
- wrap (bool) – If True, will wrap the phased LC around 0.0 to make some phased LCs easier to look at.
- model_over_lc (bool) – Usually, this function will plot the actual LC over the model LC. Set this to True to plot the model over the actual LC; this is most useful when you have a very dense light curve and want to be able to see how it follows the model.
- fitphase (optional np.array) – If passed, use this as x values for fitmags
Returns: Return type: Nothing.
-
astrobase.lcfit.utils.
iterative_fit
(data_x, data_y, init_coeffs, objective_func, objective_args=None, objective_kwargs=None, optimizer_func=<function least_squares>, optimizer_kwargs=None, optimizer_needs_scalar=False, objective_residualarr_func=None, fit_iterations=5, fit_reject_sigma=3.0, verbose=True, full_output=False)[source]¶ This is a function to run iterative fitting based on repeated sigma-clipping of fit outliers.
Parameters: - data_x (np.array) – Array of the independent variable.
- data_y (np.array) – Array of the dependent variable.
- init_coeffs – The initial values of the fit function coefficients.
- objective_func (Python function) –
A function that is used to calculate residuals between the model and the data_y array. This should have a signature similar to:
def objective_func(fit_coeffs, data_x, data_y, *objective_args, **objective_kwargs)
and return an array of residuals or a scalar value indicating some sort of sum of residuals (depending on what the optimizer function requires).
If this function returns a scalar value, you must set optimizer_needs_scalar to True, and provide a Python function in objective_residualarr_func that returns an array of residuals for each value of data_x and data_y given an array of fit coefficients.
- objective_args (tuple or None) – A tuple of arguments to pass into the objective_func.
- objective_kwargs (dict or None) – A dict of keyword arguments to pass into the objective_func.
- optimizer_func (Python function) –
The function that minimizes the residual between the model and the data_y array using the objective_func. This should have a signature similar to one of the optimizer functions in scipy.optimize, i.e.:
def optimizer_func(objective_func, initial_coeffs, args=(), kwargs={}, ...)
and return a scipy.optimize.OptimizeResult. We’ll rely on the
.success
attribute to determine if the EPD fit was successful, and the.x
attribute to get the values of the fit coefficients. - optimizer_kwargs (dict or None) – A dict of kwargs to pass into the optimizer_func function.
- optimizer_needs_scalar (bool) – If True, this indicates that the optimizer requires a scalar value to be returned from the objective_func. This is the case for scipy.optimize.minimize. If this is True, you must also provide a function in objective_residual_func.
- objective_residualarr_func (Python function) –
This is used in conjunction with optimizer_needs_scalar. The function provided here must return an array of residuals for each value of data_x and data_y given an array of fit coefficients. This is then used to calculate which points are outliers after a fit iteration. The function here must have the following signature:
def objective_residualarr_func(coeffs, data_x, data_y, *objective_args, **objective_kwargs)
- fit_iterations (int) – The number of iterations of the fit to perform while throwing out outliers to the fit.
- fit_reject_sigma (float) – The maximum deviation allowed to consider a data_y item as an outlier to the fit and to remove it from consideration in a successive iteration of the fit.
- verbose (bool) – If True, reports per iteration on the cost function value and the number of items remaining in data_x and data_y after sigma-clipping outliers.
- full_output (bool) – If True, returns the full output from the optimizer_func along with the resulting fit function coefficients.
Returns: result – If full_output was True, will return the fit coefficients np.array as the first element and the optimizer function fit output from the last iteration as the second element of a tuple. If full_output was False, will only return the final fit coefficients as an np.array.
Return type: np.array or tuple