solve_fits¶
- sofia_redux.toolkit.resampling.solve_fits(sample_indices, sample_coordinates, sample_phi_terms, sample_data, sample_error, sample_mask, fit_coordinates, fit_phi_terms, order, alpha, adaptive_alpha, is_covar=False, mean_fit=False, cval=nan, fit_threshold=0.0, error_weighting=True, estimate_covariance=False, order_algorithm_idx=1, order_term_indices=None, derivative_term_map=None, derivative_term_indices=None, edge_algorithm_idx=1, edge_threshold=None, minimum_points=None, get_error=True, get_counts=True, get_weights=True, get_distance_weights=True, get_rchi2=True, get_cross_derivatives=True, get_offset_variance=True)[source]¶
Solve all fits within one intersection block.
This function is a wrapper for
solve_fit()
over all data sets and fit points. The main computations here involve:Creating and populating the output arrays.
Selecting the correct samples within the region of each fitting window.
Calculating the full weighting factors for the fits.
For further details on the actual fitting, please see
solve_fit()
.- Parameters:
- sample_indicesnumba.typed.List
A list of 1-dimensional numpy.ndarray (dtype=int) of length n_fits. Each list element
sample_indices[i]
, contains the indices of samples within the “window” region offit_indices[i]
.- sample_coordinatesnumpy.ndarray (n_dimensions, n_samples)
The independent coordinates for each sample in n_dimensions
- sample_phi_termsnumpy.ndarray (n_terms, n_samples)
The polynomial terms of
sample_coordinates
. Please seepolynomial_terms()
for further details.- sample_datanumpy.ndarray (n_sets, n_samples)
The dependent values of the samples for n_sets, each containing n_samples.
- sample_errornumpy.ndarray (n_sets, n_samples)
The associated 1-sigma error values for each sample in each set. The user may also supply an array of shape (n_sets, 1) in which case all samples in a set will share the same associated error value. If the shape is set to (n_sets, 0), this indicates that no error values are available for the samples.
- sample_masknumpy.ndarray (n_sets, n_samples)
A mask where
False
indicates that the associated sample should be excluded from all fits.- fit_coordinatesnumpy.ndarray (n_dimensions, n_fits)
The independent variables at each fit coordinate in d_dimensions.
- fit_phi_termsnumpy.ndarray (n_terms, n_fits)
The polynomial terms of
fit_coordinates
. Please seepolynomial_terms()
for further details.- ordernumpy.ndarray
The desired order of the fit as a (1,) or (n_dimensions,) array. If only a single value is supplied, it will be applied over all dimensions.
- alphanumpy.ndarray (n_dimensions,)
A distance weighting scaling factor per dimension. The weighting kernel is applied equally to all sets and samples. For further details, please see
calculate_distance_weights()
. Will be overridden byadaptive_alpha
ifadaptive_alpha.size > 0
.- adaptive_alphanumpy.ndarray
Shape = (n_samples, n_sets, [1 or n_dimensions], n_dimensions). Defines a weighting kernel for each sample in each set. The function
calculate_adaptive_distance_weights_scaled()
will be used for kernels of shape (1, n_dimensions), andcalculate_adaptive_distance_weights_shaped()
will be used for kernels of shape (n_dimensions, n_dimensions).adaptive_alpha
is a required parameter due to Numba constraints, and will override thealpha
parameter unless it has a size of 0. Therefore, to disable, please set the size of any dimension to zero.- is_covarbool, optional
If
True
, indicates thatsample_data
contains covariance values that should be propagated through algorithm. If this is the case, polynomial fitting is disabled, and a weighted variance is calculated instead.- mean_fitbool, optional
If
True
, a weighted mean is performed instead of calculating a polynomial fit.- cvalfloat, optional
In a case that a fit is unable to be calculated at certain location,
cval
determines the fill value for the outputfit
array at those locations.- fit_thresholdfloat, optional
If fit_threshold is non-zero, perform a check on the goodness of the fit. When the reduced-chi statistic is greater than abs(fit_threshold), the fit is determined to be a failure, and a replacement value is used. If
fit_threshold
< 0, failed fit values will be set tocval
. Iffit_threshold
> 0, failed fit values will be replaced by the weighted mean.- error_weightingbool, optional
If
True
, weight the samples in the fit by the inverse variance (1 / sample_error^2) in addition to distance weighting.- estimate_covariancebool, optional
If True, calculate the covariance of the fit coefficients using
estimated_covariance_matrix_inverse()
. Otherwise, usecovariance_matrix_inverse()
.- order_algorithm_idxint, optional
An integer specifying which polynomial order validation algorithm to use. The default (1), will always be the more robust of all available options. For further information, please see
check_edges()
.- order_term_indicesnumpy.ndarray (> max(order) + 1,), optional
A 1-dimensional lookup array for use in determining the correct phi terms to use for a given polynomial order. The order validation algorithm ensures a fit of the requested order is possible. If not, and the orders are equal in all dimensions, it may also optionally return a suggested order. In this case,
order_term_indices
is used to select the correctsample_phi_terms
andfit_phi_terms
for a given order (k), where terms are extracted viaphi[order_term_indices[k]:order_term_indices[k + 2]]
.- derivative_term_mapnumpy.ndarray, optional
A mapping array for the determination of derivatives from the coefficients of the fit, and available terms in “phi”. The shape of the array is (n_dimensions, 3, n_derivative_terms). This is only required if the gradient is required as an output. For a full description of the derivative map, please see
polynomial_derivative_map()
.- derivative_term_indicesnumpy.ndarray (max(order) + 1,), optional
If the fit order is allowed to vary, gives the indices in
derivative_term_map
for a given symmetrical order. The correctderivative_term_map
mapping for order k is given asderivative_term_map[:, :, indices[k]:indices[k + 2]]
.- edge_algorithm_idxint, optional
Integer specifying the algorithm used to determine whether a fit should be attempted with respect to the sample distribution. Please see
check_edges()
for further information. The default (1), is always the most robust of the available algorithms.- edge_thresholdnumpy.ndarray (n_dimensions,)
A threshold parameter determining how close an edge should be to the center of the distribution during
check_edges()
. Higher values result in an edge closer to the sample mean. A value should be provided for each dimension. A zero value in any dimension will result in an infinite edge for that dimension.- minimum_pointsint, optional
Certain order validation algorithms check the number of available samples as a means to determine what order of fit is appropriate. If pre-calculated for the base
order
, it may be passed in here for a slight speed advantage.- get_errorbool, optional
If
True
, return the error on the fit.- get_countsbool, optional
If
True
, return the number of samples used when determining the fit at each fitting point.- get_weightsbool, optional
If
True
, return the sum of all sample weights used in determining the fit at each point.- get_distance_weightsbool, optional
If
True
, return the sum of only the distance weights used in determining the fit at each point.- get_rchi2bool, optional
If
True
, return the reduced chi-squared statistic for each of the fitted points.- get_cross_derivativesbool, optional
If
True
, return the derivative mean-squared-cross-products of the samples for each of the fitted points. Seederivative_mscp()
for further information.- get_offset_variancebool optional
If
True
, return the offset of the fitting point from the sample distribution. Seeoffset_variance()
for further information.
- Returns:
- fit_results8-tuple of numpy.ndarray
fit_results[0]: Fitted values. fit_results[1]: Error on the fit. fit_results[2]: Number of samples in each fit. fit_results[3]: Weight sums. fit_results[4]: Distance weight sums. fit_results[5]: Reduced chi-squared statistic. fit_results[6]: Derivative mean-squared-cross-products. fit_results[7]: Offset variances from the sample distribution center.
All arrays except for fit_results[6] have the shape (n_sets, n_fits) or (0, 0) depending on whether
get_
isTrue
orFalse
respectively. The derivative MSCP is of shape (n_sets, n_fits, n_dimensions, n_dimensions) if requested, and (1, 0, 0, 0) otherwise.