iltpy.fit
iltpy.fit.inversion
- class iltpy.fit.inversion.IltData(data_dict)[source]
Bases:
object
Main data handling class for ILTpy.
Attributes
- datandarray
The input data to be inverted.
- tlist of ndarray
t (input sampling) vectors for each dimension
- taulist of ndarray
Tau vectors for each dimension.
- dataset_typestr
Refers to how the data was loaded, ndarray or txt
- data_dictdict
Dictionary containing ‘dataset_type’, ‘data’, ‘dims’ and ‘filepath’ as keys
- kernellist of functions
Kernel function object used for inversion.
- fit_dictdict
Dictionary of fit parameters before initialization as provided by user.
- fit_dict_initdict
Dictionary of fit parameters after initialization.
- gndarray
The computed distribution after inversion.
- g_guessndarray, optional
User-provided initial guess for the distribution.
- fitndarray
The fitted data after inversion.
- residualsndarray
Difference between fitted data and input data.
- statisticsdict
Statistical analysis results after running iltstats().
- invertedbool
Indicates whether inversion has been performed.
- initializedbool
Indicates whether the object has been initialized.
- convergencelist
List of convergence values for each iteration.
- Gamma2sparse matrix
Regularization matrix used in inversion.
- Vbcsparse matrix
Boundary condition regularization matrix.
- nnbool
If True, non-negativity is enforced.
- maxloopint
Maximum number of inversion iterations.
- conv_limitfloat
Convergence threshold for stopping iterations.
- alt_gint
Alternative strategy for calculating g (0, 2, 4, or 6).
- force_sparsebool
If True, forces sparse matrix operations, mainly applicable for alt_g = 0.
- sparse_thresholdfloat
Threshold for treating kernel matrices as sparse.
- iltstats(n=10, solver=None, n_jobs=-1, verbose=1, random_seed=None)[source]
Calculates statistical parameters for a distribution by adding random noise to the initial fit of the data and inverting n times.
Parameters
- nint, optional
The number of times the data is refit with different noise vectors, must be at least 2. If None, then n is set equal to 10.
- solverfunction, optional
The solver used for inversion.
- n_jobsint, default -1
The number of parallel jobs to create, by default -1. This will create as many jobs as cpu cores present on the system.
- verboseint, default 0
Controls joblib’s verbosity.
- random_seedint, optional
Random seed for reproducible noise vector initialization.
- init(tau=None, kernel=None, parameters=None, dim_ndim=None, sigma=None, g_guess=None, verbose=0, reinitialize=False)[source]
Initializes IltData. The initialization steps are as follows:
init_ndim : Sets dimension-related properties
init_fit_dict : Sets the default parameters
set_parameters : Sets values of fit_dict as class attributes
init_parameters : Initializes default and any user-defined parameters, consistency checks, etc.
init_kernelmatrix : Computes the kernel matrix, LHS, RHS of the Ax=B equation and initializes the weighting matrix if provided.
init_regmatrix : Computes the initial regularization matrices.
Parameters
- taundarray, optional
The tau vector used by the inversion,if None, a logarithmically spaced tau is generated (only valid for relaxation data)
- kernelobject, parent class IltKernel
A class defining the kernel function.
- dim_ndimndarray, optional
The array defining the n-th dimension of n-D data in case of an uninverted but regularized dimension.
- parametersdict, optional
User-defined parameters, provided as a dict. Note that data, kernel, tau, g_guess, dim_ndim or sigma can’t be provided as parameters and must be provided using keyword arguments.
- sigmandarray
Array used for weighted inversions. It must have the same shape as the data.
- g_guessndarray
User-defined guess for the initial g. It must have the same shape as the expected distribution.
- verboseint
Controls verbosity.
- reinitializebool
If IltData is already initialized and if this argument is set to True, default parameters will be reloaded again.
- invert(solver=None, verbose=1)[source]
Perform Inversion with Regularization
The method iteratively solves the system (W + gamma_2)g = Y, where: - W is a sparse matrix representing the system, - g is the spectrum to be solved for, - Y is the observed data, and - gamma_2 is the regularization matrix, which includes terms for uniform penalty, zero-crossing penalty, and non-negativity.
The regularization penalties are updated at each iteration, and the system is solved using either a provided solver or a default solver.
Parameters
- solverfunction, optional
A function that solves sparse linear systems. The solver should take two arguments: - Sparse matrix A (representing W+gamma2) - Dense vector b (Y) and return the spectrum g as a flattened 1D numpy array of shape (Nxtot,). If not provided, the default solver is used.
- verboseint, default 1
Controls verbosity. Setting this to zero will not display the progress bar.
- store_gbool, default False
If True, g (distribution or spectrum) after each iteration is stored in the attrbiute g_list
Raises
- RuntimeError
If the IltData object has not been initialized with IltData.init(tau), an error is raised.
- ValueError
If the provided initial guess for g does not match the expected size.
- report(filepath=None, save_arrays=False, level=0, parameters_to_file=False)[source]
Generates .txt files from results of inversion. Requires the dataset to be inverted.
Parameters
- filepathstr, optional
Path where the results should be saved. If None, results are saved in the current working directory, by default None
- save_arraysbool, by default False
Whether to also output a file with numpy arrays in npz format.
- levelint, optional
Controls the depth of saved data. Default is 0. - 0: Save attributes of type np.ndarray. - 1: Additionally, save lists containing only np.ndarray elements. - 2: Additionally, save sparse matrices (converted to dense arrays).
- parameters_to_filebool, optional
If True, saves the fit parameters to a separate file along with the results.
- class iltpy.fit.inversion.IltSolver(iltdata)[source]
Bases:
object
Class for handling alternative strategies for calculating g.
- alt_g_solve(iltdata)[source]
Switches the strategy for calculating g based on the value of alt_g.
Parameters
- iltdataIltData
An instance of IltData class which is already initialized.
Returns
- ndarray
g as a flattened array
Raises
- ValueError
If alt_g parameter is not one of the accepted values : [0,2,4,6]
- get_default_solver()[source]
Selects the default solver in case of alt_g=0 based on kernel sparsity. This function either returns scipy’s spsolve or numpy’s solve.
The permc_spec=’NATURAL’ parameter specifies the natural ordering of the matrix for sparse LU decomposition. This ordering is chosen as the default because it avoids additional reordering overhead and is generally efficient for well-structured matrices, which are common in this application.
Returns
- function
A solver function as a lambda function based on kernel sparsity. For dense kernels, scipy.linalg.solve and for sparse kernels, scipy.sparse.linalg.spsolve is returned.
iltpy.fit.kernels
- class iltpy.fit.kernels.IltKernel[source]
Bases:
object
Parent class for kernels. For defining kernels:
- class MyKernel(IltKernel):
- def __init__(self):
self.name = ‘my_kernel’ self.kernel = lambda t,tau : kernel_func(t,tau) self.kernel_str = “lambda t,tau : kernel_func(t,tau)”
- iltpy.fit.kernels.init_kernel_dense(iltdata, verbose=0)[source]
This function initializes kernel matrices if sparsity of kernel matrices is less than iltdata.sparse_threshold. It uses the dense form of matrices and uses corresponding functions. The sparse version of this function can be forced by setting iltdata.force_sparse to True during initialization.
Parameters
- iltdataIltData
An instance of IltData class
- verboseint
Controls verbosity
- iltpy.fit.kernels.init_kernel_sparse(iltdata, verbose=0)[source]
This function initializes kernel matrices if sparsity of kernel matrices is larger than iltdata.sparse_threshold. It uses the sparse form of matrices and uses corresponding functions.
Parameters
- iltdataIltData
An instance of IltData class
- verboseint
Controls verbosity
iltpy.fit.regularization
- iltpy.fit.regularization.amp_reg(iltdata)[source]
Amplitude regularization
Parameters
- iltdataIltData
An instance of IltData class.
- iltpy.fit.regularization.init_g_red(iltdata)[source]
Initialize g_red
Parameters
- iltdataIltData
An instance of IltData class.
- iltpy.fit.regularization.init_regmatrix(iltdata, verbose=0)[source]
This function generates and initializes the regularization matrices.
Parameters
- iltdataIltData
An instance of IltData class.
- verboseint, optional
Controls verbosity, by default 0
- iltpy.fit.regularization.nn_reg(iltdata)[source]
Non negativity constraint.
Parameters
- iltdataIltData
An instance of IltData class.
- iltpy.fit.regularization.up_reg(iltdata)[source]
Uniform Penalty regularization.
Parameters
- iltdataIltData
An instance of IltData class.
- iltpy.fit.regularization.zc_reg(iltdata, it)[source]
Zero crossing regularization.
This function computes the zero-crossing (ZC) penalty matrix. The ZC penalty discourages oscillatory solutions in the spectrum g by penalizing zero-crossings. After zc_down number of iterations, zc_ratio is updated.
Parameters
- iltdataIltData
An instance of IltData class.
- itint
Number of iterations