Code Documentation
Here we give detailed code documentation of the classes, constructors, and functions within APPRENTICE. Here is an index of all the code documentation presented in this page.
Surrogate model
Predition function
Function Optimization
Surrogate models
Base class
SurrogateModel (apprentice/surrogatemodel.py)
- class apprentice.surrogatemodel.SurrogateModel(dim, fnspace=None, **kwargs)
Bases:
objectSurrogate model base class
- __init__(dim, fnspace=None, **kwargs)
Initialize surrogate models
- Parameters
dim (int) – parameter dimension
fnspace (apprentice.space.Space) – function space
- classmethod from_interpolation_points(X=None, Y=None, **kwargs)
A class method to construct surrogate model from interpolation points
- Parameters
X (list) – a 2-D array of size \(dim imes N_p\) and it is the x data values to fit where \(dim\) is the parameter dimension and \(N_p\) is the the number of data points
Y (list) – an array of size \(N_p\) and it is the y data values to fit where \(N_p\) is the the number of data points
- classmethod from_data_structure(data_structure, **kwargs)
A class method to construct surrogate model from data structure. This function needs to be implemented in a class that inherits this class
- Parameters
data_structure (dict) – previously fit surrogate model saved as a data structure
- classmethod from_file(filepath, **kwargs)
A class method to construct surrogate model from data structure saved in a file. This function needs to be implemented in a class that inherits this class
- Parameters
filepath (str) – file path of a previously fit surrogate model saved as a data structure
- fit(X, Y)
Surrogate model fitting method. This function needs to be implemented in a class that inherits this class
- Parameters
X (list) – a 2-D array of size \(dim imes N_p\) and it is the x data values to fit where \(dim\) is the parameter dimension and \(N_p\) is the the number of data points
Y (list) – an array of size \(N_p\) and it is the y data values to fit where \(N_p\) is the the number of data points
- f_x(x)
Calculate the surrogate model at a new point. This function needs to be implemented in a class that inherits this class
- Parameters
x (list) – a new x point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
surrogate model value at the new point
- Return type
float
- f_X(x)
Calculate the surrogate model at a multiple date points. This function needs to be implemented in a class that inherits this class
- Parameters
X (list) – multiple x point, an araay of size \(dim imes n\) where \(dim\) is the parameter dimension and \(n\) is the number of new data points
- Returns
surrogate model value at the new points an araay of size \(n\) where \(n\) is the number of new data points
- Return type
list
- property training_size
Get the training size
- Returns
the training size \(N_p\)
- Return type
int
- property function_space
Get the function space object
- Returns
function space object
- Return type
- property has_gradient
Check if gradient is implemented
- Returns
true if an implementation of gradient is found.
- Return type
bool
- property has_hessian
Check if hessian is implemented
- Returns
true if an implementation of hessian is found.
- Return type
bool
- property fnspace
Get the function space object
- Returns
function space object
- Return type
- property bounds
Get function space bounds
- Returns
2-D array of size \(dim ime 2\) containing bounds of parameter space
- Return type
list
- property dim
Get the dimension of the parameter space
- Returns
dimension value
- Return type
int
Available implementation
PolynomialApproximation (apprentice/polynomialapproximation.py)
- class apprentice.polynomialapproximation.PolynomialApproximation(dim, fnspace=None, **kwargs: dict)
Bases:
SurrogateModelPolynomial approximation surrogate model
- __init__(dim, fnspace=None, **kwargs: dict)
Polynomial approximation surrogate model construction function
- Parameters
dim (int) – parameter dimension
fnspace (apprentice.space.Space) – function space object
- property dim
Get the parameter dimension value
- Returns
parameter dimension
- Return type
int
- property pnames
Get names of the parameter dimension names
- Returns
array of parameter dimension string names
- Return type
list
- property order_numerator
Get the numerator order
- Returns
numerator order
- Return type
int
- property fit_strategy
Get the fit strategy
- Returns
fit strategy
- Return type
int
- property to_compute_covariance
Get the choice of whether covariance should be computed
- Returns
true if covariance should be computed
- Return type
bool
- property coeff_numerator
Get the numerator coefficients. The order of coefficients is as in https://people.math.sc.edu/Burkardt/py_src/polynomial/polynomial.html
- Returns
list of numerator coefficients.
- Return type
list
- set_structures()
Set monomial structures into self. The order is as in https://people.math.sc.edu/Burkardt/py_src/polynomial/polynomial.html
- coeff_solve(VM, Y)
Solve coefficients using Singular Value Decomposition (SVD)
- Parameters
VM (np.array) – Vandermonde matrix
Y (np.array) – an array of size \(N_p\) and it is the y data values to fit where \(N_p\) is the the number of data points
- coeff_solve2(VM, Y)
Solve coefficients using least squares regression
- Parameters
VM (np.array) – Vandermonde matrix
Y (np.array) – an array of size \(N_p\) and it is the y data values to fit where \(N_p\) is the the number of data points
- fit(X, Y)
Surrogate model fitting method.
- Parameters
X (list) – a 2-D array of size \(dim imes N_p\) and it is the x data values to fit where \(dim\) is the parameter dimension and \(N_p\) is the the number of data points
Y (np.array) – an array of size \(N_p\) and it is the y data values to fit where \(N_p\) is the the number of data points
- f_x_slow(x)
Calculate the surrogate model at a new point. This is a slower version.
- Parameters
x (list) – a new x point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
surrogate model value at the new point
- Return type
float
- f_x(x)
Calculate the surrogate model at a new point.
- Parameters
x (list) – a new x point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
surrogate model value at the new point
- Return type
float
- f_X(X)
Calculate the surrogate model at a multiple date points.
- Parameters
X (list) – multiple x point, an araay of size \(dim imes n\) where \(dim\) is the parameter dimension and \(n\) is the number of new data points
- Returns
surrogate model value at the new points an araay of size \(n\) where \(n\) is the number of new data points
- Return type
list
- property as_dict
Get the polynomial approximation fit as a dictionary
- Returns
polynomial approximation fit dictionary
- Return type
dict
- save(fname)
Save the polynomial approximation fit into a file
- Parameters
fname (str) – file path to save polynomial approximation fit
- classmethod from_data_structure(data_structure, **kwargs)
A class method to construct surrogate model from data structure.
- Parameters
data_structure (dict) – previously fit surrogate model saved as a data structure
- classmethod from_file(filepath, **kwargs)
A class method to construct surrogate model from data structure saved in a file.
- Parameters
filepath (str) – file path of a previously fit surrogate model saved as a data structure
- property coeff_norm
Get 1-norm of the polynomial approximation fit coefficients
- Returns
1-norm of the polynomial approximation fit coefficients
- Return type
float
- property coeff2_norm
Get 2-norm of the polynomial approximation fit coefficients
- Returns
2-norm of the polynomial approximation fit coefficients
- Return type
float
- gradient(X)
Get gradient of the polynomial approximation
- Parameters
X (list) – a new point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
gradient of the polynomial approximation at a new point
- Return type
list
- hessian(X)
Get hessian of the polynomial approximation
- Parameters
X (list) – a new point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
hessian of the polynomial approximation at a new point
- Return type
list
RationalApproximation (apprentice/rationalapproximation.py)
- class apprentice.rationalapproximation.RationalApproximation(dim, fnspace=None, **kwargs: dict)
Bases:
SurrogateModelRational approximation surrogate model
- __init__(dim, fnspace=None, **kwargs: dict)
Rational approximation surrogate model construction function
- Parameters
dim (int) – parameter dimension
fnspace (apprentice.space.Space) – function space object
- property max_restarts
Get max restarts
- Returns
max restarts
- Return type
int
- property threshold
Get threshold value
- Returns
threshold value
- Return type
float
- property ftol
Get function tolerance
- Returns
function tolerance
- Return type
float
- property iterations
Get number of iterations
- Returns
number of iterarations
- Return type
int
- property tmpdir
Get PYOMO temp directory
- Returns
temp directory path
- Return type
str
- property debug
Get debug value
- Returns
true if debug mode is on, false otherwise
- Return type
bool
- property use_abstract_model
Whether to use PYOMO abstract model
- Returns
true if to use PYOMO abstract model
- Return type
bool
- property dim
Get the parameter dimension value
- Returns
parameter dimension
- Return type
int
- property pnames
Get names of the parameter dimension names
- Returns
array of parameter dimension string names
- Return type
list
- property order_numerator
Get the numerator order
- Returns
numerator order
- Return type
int
- property order_denominator
Get the denominator order
- Returns
denominator order
- Return type
int
- property fit_solver
Get the fit solver
- Returns
fit solver
- Return type
str
- property local_solver
Get the local order
- Returns
local solver
- Return type
str
- property fit_strategy
Get the fit strategy
- Returns
fit strategy
- Return type
str
- property coeff_numerator
Get the numerator coefficients. The order of coefficients is as in https://people.math.sc.edu/Burkardt/py_src/polynomial/polynomial.html
- Returns
list of numerator coefficients.
- Return type
list
- property coeff_denominator
Get the denominator coefficients. The order of coefficients is as in https://people.math.sc.edu/Burkardt/py_src/polynomial/polynomial.html
- Returns
list of denominator coefficients.
- Return type
list
- set_structures()
Set monomial structures into self. The order is as in https://people.math.sc.edu/Burkardt/py_src/polynomial/polynomial.html
- property M
Get number of numerator coefficients
- Returns
number of numerator coefficients
- Return type
int
- property N
Get number of denominator coefficients
- Returns
number of denominator coefficients
- Return type
int
- coeff_solve(VM, VN, Y)
Find coefficients of the rational approximation. The resulting rational approximation may contain poles using the linear algebra technique proposed in our paper https://www.sciencedirect.com/science/article/abs/pii/S0010465520303222 [https://arxiv.org/abs/1912.02272]
- Parameters
VM (np.array) – Vandermonde matrix of numerator
VN (np.array) – Vandermonde matrix of denominator
Y (np.array) – an array of size \(N_p\) and it is the y data values to fit where \(N_p\) is the the number of data points
- coeff_solve2(VM, VN, Y)
Find coefficients of the rational approximation. F = p/q is reformulated as 0 = p - qF using the Vandermonde matrices. That defines the problem Ax = b and we solve for x using using Singular Value Decomposition (SVD), exploiting A = U x S x V.T. There is an additional manipulation exploiting on setting the constant coefficient in q to 1. The resulting rational approximation may contain poles.
- Parameters
VM (np.array) – Vandermonde matrix of numerator
VN (np.array) – Vandermonde matrix of denominator
Y (np.array) – an array of size \(N_p\) and it is the y data values to fit where \(N_p\) is the the number of data points
- coeff_solve3(VM, VN, Y)
Find coefficients of the rational approximation F = p/q is reformulated as 0 = p - qF using the Vandermonde matrices. That defines the problem Ax = b and we solve for x using using Singular Value Decomposition (SVD), exploiting A = U x S x V.T. We get the solution as the last column in V (corresponds to the smallest singular value). The resulting rational approximation may contain poles.
- Parameters
VM (np.array) – Vandermonde matrix of numerator
VN (np.array) – Vandermonde matrix of denominator
Y (np.array) – an array of size \(N_p\) and it is the y data values to fit where \(N_p\) is the the number of data points
- fit_order_denominator_one(X, Y)
Find pole-free coefficients of the rational approximation with numerator order \(m\) and denominator order \(1\)
- Parameters
X (list) – a 2-D array of size \(dim imes N_p\) and it is the x data values to fit where \(dim\) is the parameter dimension and \(N_p\) is the the number of data points
Y (np.array) – an array of size \(N_p\) and it is the y data values to fit where \(N_p\) is the the number of data points
- fit_sip(X, Y)
Find coefficients of the rational approximation with numerator order \(m\) and denominator order \(n\) with minimal number of poles. The number of poles is reduced using the semi-infinite programming technique proposed in our paper https://www.sciencedirect.com/science/article/abs/pii/S0010465520303222 [https://arxiv.org/abs/1912.02272]
- Parameters
X (list) – a 2-D array of size \(dim imes N_p\) and it is the x data values to fit where \(dim\) is the parameter dimension and \(N_p\) is the the number of data points
Y (np.array) – an array of size \(N_p\) and it is the y data values to fit where \(N_p\) is the the number of data points
- fit(X, Y)
Surrogate model fitting method.
- Parameters
X (np.array) – a 2-D array of size \(dim imes N_p\) and it is the x data values to fit where \(dim\) is the parameter dimension and \(N_p\) is the the number of data points
Y (np.array) – an array of size \(N_p\) and it is the y data values to fit where \(N_p\) is the the number of data points
- Q(X)
Evaluation of the denominator polynomial at new point.
- Parameters
X (list) – a new point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
denominator polynomial at new point
- Return type
np.array
- denom(X)
Evaluation of the denominator polynomial at new point.
- Parameters
X (list) – a new point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
denominator polynomial at new point
- Return type
np.array
- P(X)
Evaluation of the numerator polynomial at new point.
- Parameters
X (list) – a new point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
numerator polynomial at new point
- Return type
np.array
- predict(X)
Evaluate the rational approximation at a new point
- Parameters
X (list) – a new point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
surrogate model value at the new point
- Return type
float
- gradient(X)
Get gradient of the rational approximation
- Parameters
X (list) – a new point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
gradient of the rational approximation at a new point
- Return type
list
- f_x(x)
Calculate the surrogate model at a new point.
- Parameters
x (list) – a new x point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
surrogate model value at the new point
- Return type
float
- f_X(X)
Calculate the surrogate model at a multiple date points.
- Parameters
X (list) – multiple x point, an araay of size \(dim imes n\) where \(dim\) is the parameter dimension and \(n\) is the number of new data points
- Returns
surrogate model value at the new points an araay of size \(n\) where \(n\) is the number of new data points
- Return type
list
- property as_dict
Get the rational approximation fit as a dictionary
- Returns
rational approximation fit dictionary
- Return type
dict
- save(fname)
Save the rational approximation fit into a file
- Parameters
fname (str) – file path to save rational approximation fit
- classmethod from_data_structure(data_structure, **kwargs)
A class method to construct surrogate model from data structure.
- Parameters
data_structure (dict) – previously fit surrogate model saved as a data structure
- classmethod from_file(filepath, **kwargs)
A class method to construct surrogate model from data structure saved in a file.
- Parameters
filepath (str) – file path of a previously fit surrogate model saved as a data structure
- property coeff_norm
Get 1-norm of the rational approximation fit coefficients
- Returns
1-norm of the rational approximation fit coefficients
- Return type
float
- property coeff2_norm
Get 2-norm of the polynomial approximation fit coefficients
- Returns
2-norm of the polynomial approximation fit coefficients
- Return type
float
GaussianProcess (apprentice/gaussianprocess.py)
- class apprentice.gaussianprocess.GaussianProcess(dim, fnspace=None, **kwargs: dict)
Bases:
SurrogateModelGaussian process surrogate model
- __init__(dim, fnspace=None, **kwargs: dict)
Gaussian process surrogate model construction function
- Parameters
dim (int) – parameter dimension
fnspace (apprentice.space.Space) – function space object
- property use_mpi_tune
Check if to use MPI tune
- Returns
true if to use MPI tune
- Return type
bool
- property stopping_bound
Get the stopping bound
- Returns
stopping bound
- Return type
float
- property polynomial_order
Find the order of the polynomial
- Returns
polynomial order
- Return type
int
- property sample_size
Get the sample size
- Returns
sample size
- Return type
int
- property fit_strategy
Get the fit strategy
- Parameters
print_descr (true if to print description, false otherwise) – whether to print description
- Returns
fit strategy
- Return type
str
- property debug
Get debug value
- Returns
true if debug mode is on, false otherwise
- Return type
bool
- property max_restarts
Get the maximum restarts
- Returns
max restarts
- Return type
int
- property mean_surrogate_model
Get the mean surrogate model
- Returns
mean surrogate model
- Return type
- property error_surrogate_model
Get the error surrogate model
- Returns
error surrogate model
- Return type
- property dim
Get the parameter dimension value
- Returns
parameter dimension
- Return type
int
- property pnames
Get names of the parameter dimension names
- Returns
array of parameter dimension string names
- Return type
list
- property kernel
Get the kernel
- Returns
kernel
- Return type
str
- property seed
Get the seed
- Returns
seed
- Return type
int
- property keepout_percentage
Get the keepout percentage
- Returns
keepout percentage
- Return type
float
- fit(X, Y)
Surrogate model fitting method.
- Parameters
X (np.array) – a 2-D array of size \(dim imes N_p\) and it is the x data values to fit where \(dim\) is the parameter dimension and \(N_p\) is the the number of data points
Y (np.array) – an array of size \(N_p\) and it is the y data values to fit where \(N_p\) is the the number of data points
- save(fname)
Save the gaussian process fit into a file
- Parameters
fname (str) – file path to save gaussian process fit
- property as_dict
Get the gaussian process fit as a dictionary
- Returns
gaussian process fit dictionary
- Return type
dict
- static get_kernel(kernel_str, dim, poly_order)
Get kernel
- Parameters
kernel_str (str) – kernel
dim (int) – parameter dimension
poly_order (int) – polynomial order
- Returns
GPy kernel
- Return type
GPy.kern
- static mpi_tune(model, num_restarts, use_mpi=False, robust=True, debug=False)
Do MPI tuning
- Parameters
model (GPy.core.GP or GPy.models) – GPy model
num_restarts (int) – number of restarts
use_mpi (bool) – whether to use MPI
robust (bool) – whether the tuning should be robust
debug (bool) – debug value
- Returns
return value
- Return type
int
- build_MLHGP_model_from_data(X, Y)
Build GP model using Most Likely Heteroscedastic Gaussian Process (HeGP-ML)
- Parameters
X (np.array) – a 2-D array of size \(dim imes N_p\) and it is the x data values to fit where \(dim\) is the parameter dimension and \(N_p\) is the the number of data points
Y (np.array) – an array of size \(N_p\) and it is the y data values to fit where \(N_p\) is the the number of data points
- Returns
Model for values (model Y), model for heteroscedastic varaince (model Z), and database of model Y & model Z information
- Return type
GPy.core.GP/GPy.models, GPy.core.GP/GPy.models, dict
- build_SK_model_from_data(X, Y)
Build GP model using Heteroscedastic Gaussian Process using Stochastic Kriging (HeGP-SK)
- Parameters
X (np.array) – a 2-D array of size \(dim imes N_p\) and it is the x data values to fit where \(dim\) is the parameter dimension and \(N_p\) is the the number of data points
Y (np.array) – an array of size \(N_p\) and it is the y data values to fit where \(N_p\) is the the number of data points
- Returns
Model for values (model Y), model for heteroscedastic varaince (model Z), and database of model Y & model Z information
- Return type
GPy.core.GP/GPy.models, GPy.core.GP/GPy.models, dict
- build_GP_model_from_data(X, Y)
Build GP model using Homoscedastic Gaussian Process (HoGP)
- Parameters
X (np.array) – a 2-D array of size \(dim imes N_p\) and it is the x data values to fit where \(dim\) is the parameter dimension and \(N_p\) is the the number of data points
Y (np.array) – an array of size \(N_p\) and it is the y data values to fit where \(N_p\) is the the number of data points
- Returns
Model for values (model Y) and database of model Y & model Z information
- Return type
GPy.core.GP/GPy.models, dict
- classmethod from_file(filepath, **kwargs)
A class method to construct surrogate model from data structure saved in a file.
- Parameters
filepath (str) – file path of a previously fit surrogate model saved as a data structure
- classmethod from_data_structure(data_structure, **kwargs)
A class method to construct surrogate model from data structure.
- Parameters
data_structure (dict) – previously fit surrogate model saved as a data structure
- static predict_static_homoscedastic(Xte, Mte, modely)
Predict static homoscedastic at test data points
- Parameters
Xte (np.array) – test data points
Mte (list) – test mean values
modely (GPy.core.GP) – Y model
- Returns
mean value and standard deviation
- Return type
np.array, np.array
- static predict_static_heteroscedastic(Xte, Mte, modely, modelz)
Predict static heteroschedastic at test data points
- Parameters
Xte (np.array) – test data points
Mte (list) – test mean values
modely (GPy.core.GP/GPy.models) – Y model
modelz (GPy.core.GP/GPy.models) – Z model
- Returns
mean value and standard deviation
- Return type
np.array, np.array
- static get_metrics(Xte, MCte, DeltaMCte, Mte, modely, modelz)
Get metrics
- Parameters
Xte (np.array) – test data points
MCte (list) – test mean values
DeltaMCte (list) – test standard deviation values
Mte (list) – test mean values
modely (GPy.core.GP/GPy.models) – Y model
modelz (GPy.core.GP/GPy.models) – Z model
- Returns
mean MSE metric, standard deviation MSE metric, \(\chi^2\)
- Return type
np.array, np.array, np.array
- predict_heteroscedastic(X)
Predict using heteroschedastic GP at new point
- Parameters
X (list) – multiple x point, an araay of size \(dim imes n\) where \(dim\) is the parameter dimension and \(n\) is the number of new data points
- Returns
surrogate model mean and standard deviation value at the new points
- Return type
list,list
- predict_homoscedastic(X)
Predict using homoscedastic GP at new point
- Parameters
X (list) – multiple x point, an araay of size \(dim imes n\) where \(dim\) is the parameter dimension and \(n\) is the number of new data points
- Returns
surrogate model mean and standard deviation value at the new points
- Return type
list,list
- f_x(x)
Calculate the surrogate model at a new point.
- Parameters
x (list) – a new x point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
surrogate model value at the new point
- Return type
float, float
- f_X(X)
Calculate the surrogate model at a multiple date points.
- Parameters
X (list) – multiple x point, an araay of size \(dim imes n\) where \(dim\) is the parameter dimension and \(n\) is the number of new data points
- Returns
surrogate model value at the new points an araay of size \(n\) where \(n\) is the number of new data points
- Return type
list
Prediction function
Base class
Function (apprentice/function.py)
- class apprentice.function.Function(dim, fnspace=None, **kwargs)
Bases:
objectFunction base class
- __init__(dim, fnspace=None, **kwargs)
Initialize the Function class
- Parameters
dim (int) – parameter dimension
fnspace (apprentice.space.Space) – function space object
- set_bounds(bmin, bmax)
Set minimum and maximum bounds
- Parameters
bmin (list) – bound minimum
bmax (list) – bound maximum
- set_fixed_parameters(fixed)
Set fixed dimension onf the function
- Parameters
fixed (list) – fixed dimensions of the function
- classmethod mk_empty(dim)
A class function to make an empty function
- Parameters
dim (int) – parameter dimension
- classmethod from_space(spc, **kwargs)
A class function to make a function from parameter space
- Parameters
spc (apprentice.space.Space) – parameter space
- classmethod from_surrogates(surrogates)
A class function to make a function from surrogates
- Parameters
surrogates (list) – surrogate models for all terms for the function
- property dim
Get the parameter dimension value
- Returns
parameter dimension
- Return type
int
- property bounds
Get bounds of the parameter space
- Returns
parameter bounds
- Return type
np.array
- property has_gradient
Check if gradient is implemented
- Returns
true if an implementation of gradient is found.
- Return type
bool
- property has_hessian
Check if hessian is implemented
- Returns
true if an implementation of hessian is found.
- Return type
bool
- objective(x)
Compute the function objective at a new data point. This function needs to be implemented in a class that inherits this class
- Parameters
x (list) – a new x point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
objective value of the function at the new point
- Return type
float
Available implementation
LeastSquares (apprentice/leastsquares.py)
- class apprentice.leastsquares.LeastSquares(dim, fnspace, data, s_val, errors, prefactors, e_val=None, **kwargs)
Bases:
FunctionLease squares objective function
- __init__(dim, fnspace, data, s_val, errors, prefactors, e_val=None, **kwargs)
Initialize the least squares function
- Parameters
dim (int) – parameter dimension
fnspace (apprentice.space.Space) – function space object
data (np.array) – data values
s_val (apprentice.polyset.PolySet or list) – surrogate polyset or list of SurrogateModel
errors (np.array) – data errors
prefactors (np.array) – prefactors (weights)
e_val (apprentice.polyset.PolySet or list) – surrogate polyset or list of SurrogateModel
- objective(x)
Compute least squares objective function value at new data point x
- Parameters
x (list) – a new x point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
least squares objective function value at new data point x
- Return type
float
- gradient(x)
Get gradient of the least squares function
- Parameters
x (list) – a new point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
gradient of the least squares function
- Return type
np.array
- hessian(x)
Get hessian of the least squares function
- Parameters
x (list) – a new point, an araay of size \(dim\) where \(dim\) is the parameter dimension
- Returns
hessian of the least squares function
- Return type
np.array
GeneratorTuning (apprentice/generatortuning.py)
- class apprentice.generatortuning.GeneratorTuning(dim, fnspace, data, errors, s_val, e_val, weights, binids, bindn, binup, **kwargs)
Bases:
LeastSquaresGenerator Tuning function
- __init__(dim, fnspace, data, errors, s_val, e_val, weights, binids, bindn, binup, **kwargs)
- Parameters
dim (int) – parameter dimension
fnspace (apprentice.space.Space) – function space object
data (np.array) – data values
s_val (apprentice.polyset.PolySet or list) – surrogate polyset or list of SurrogateModel
errors (np.array) – data errors
weights (np.array) – weights
e_val (apprentice.polyset.PolySet or list) – surrogate polyset or list of SurrogateModel
binids (list) – list of binids
bindn (list) – bin down
binup (list) – bin up
- initWeights(fname, hnames, bnums, blows, bups)
Initialize weights
- Parameters
fname (str) – read weight file
hnames (list) – observable names
bnums (list) – bin ids
blows (np.array) – lower end of bins
bups (lower end of bins) – upper end of bins
- Returns
list of weights
- Return type
np.array
- setWeights(wdict)
Convenience function to update the bins weights. NOTE that hnames is in fact an array of strings repeating the histo name for each corresp bin
- Parameters
wdict (dict) – weight dictionary
- mkReduced(keep, **kwargs)
Make reduced function
- Parameters
keep (list) – terms to keep, true if to keep and false if not to keep
- Returns
reduced tuning objective object
- Return type
apprentice.appset.TuningObjective2
Function Optimization
Base class
Minimizer (apprentice/minimizer.py)
- class apprentice.minimizer.Minimizer(function, **kwargs)
Bases:
objectMinimizer (Optimization) base class
- __init__(function, **kwargs)
Initialize minimizer object
- Parameters
function (apprentice.function.Function) – function
- set_bounds(bounds)
Set bounds of the optimization
- Parameters
bounds (np.array) – bounds of optimization
- minimize(x0)
Minimize. This function needs to be implemented in a class that inherits this class
- Parameters
x0 (np.array) – starting point
- sample(npoints, **kwargs)
Sample the domain
- Parameters
npoints (int) – number of points in the sample
- Returns
sample of the domain
- Return type
np.array
Available implementation
ScipyMinimizer (apprentice/minimizer.py)
- class apprentice.scipyminimizer.ScipyMinimizer(function, **kwargs)
Bases:
MinimizerSciPy minimizer
- __init__(function, **kwargs)
Initialize minimizer object
- Parameters
function (apprentice.function.Function) – function
- minimize(x0=None, nrestart=1, method='tnc', tol=1e-06)
Minimize
- Parameters
x0 (np.array) – starting point
nrestart (int) – number of restarts
method (str) – solver method name
tol (float) – tolerance
- Returns
minimization result
- Return type
dict
- minimiseTNC(x0, tol=1e-06)
Minimize using Truncated Newton (TNC) algorithm
- Parameters
x0 (np.array) – starting point
tol (float) – tolerance
- Returns
minimization result
- Return type
dict
- minimiseLBFGSB(x0, tol=1e-06)
Minimize using L-BFGS-B algorithm
- Parameters
x0 (np.array) – starting point
tol (float) – tolerance
- Returns
minimization result
- Return type
dict
- minimiseSLSQP(x0, tol=1e-06)
Minimize using Sequential Least Squares Programming (SLSQP)
- Parameters
x0 (np.array) – starting point
tol (float) – tolerance
- Returns
minimization result
- Return type
dict
- minimizeNCG(x0, sel=slice(None, None, None), tol=1e-06)
Minimize using Newton conjugate gradient trust-region algorithm
- Parameters
x0 (np.array) – starting point
sel (slice) – selected bins
tol (float) – tolerance
- Returns
minimization result
- Return type
dict
Monomial
Utility functions for building and working with monomials
- apprentice.monomial.mono_next_grlex(x)
Return next monomial. This is a deterministic procedure.
- Parameters
x (list) – current monomial
- Returns
next monomial
- Return type
list
- Example
mono_next_grlex([0,0,0])returns[0,0,1][0,0,0]corresponds to \(x^0y^0z^0\)[0,0,1]corresponds to \(x^0y^0z^1\)
- apprentice.monomial.genStruct(mnm)
Generator for mono_next_grlex. Allows for a deterministic generation of monomial terms.
- Parameters
mnm (list) – zero array of length dimension
- Returns
generated structure
- Return type
list
- Example
g = genStruct(np.zeros(10)) for i in range(10): print(next(g)) [0. 0. 0. 0. 0. 0. 0. 0. 0. 0.] [0. 0. 0. 0. 0. 0. 0. 0. 0. 1.] [0. 0. 0. 0. 0. 0. 0. 0. 1. 0.] [0. 0. 0. 0. 0. 0. 0. 1. 0. 0.] [0. 0. 0. 0. 0. 0. 1. 0. 0. 0.] [0. 0. 0. 0. 0. 1. 0. 0. 0. 0.] [0. 0. 0. 0. 1. 0. 0. 0. 0. 0.] [0. 0. 0. 1. 0. 0. 0. 0. 0. 0.] [0. 0. 1. 0. 0. 0. 0. 0. 0. 0.] [0. 1. 0. 0. 0. 0. 0. 0. 0. 0.]
- apprentice.monomial.monomialStructure(dim, order)
Generate the monomial structure of a dim-dimensional polynomial of order order. The rows are the terms. The columns are the dimensions.
- Parameters
dim (int) – dimension monomial
order (int) – order of monomial
- Returns
array of terms (rows)
- Return type
numpy.ndarray
- Example
[2,3,4]means \(x^2y^3z^4\)[[0 0],[0 1],[1 0],...]means \(x^0y^0, x^0y^1, x^1y^0,...\)
- apprentice.monomial.recurrence1D(x, structure)
Evaluate recurrence relation for a one-dimensional polynomial structure.
- Parameters
x (float) – parameter value
structure (numpy.ndarray) – monomial structure (1D)
- Returns
individual terms of f(x) without multiplication by the coefficients with \(f(x) = ax^0 + bx^1 + cx^2 +dx^3\)
- Return type
numpy.ndarray
- Exaxmple
let
x = 0.1and the polynomial order be 3, thenstructure = [0,1,2,3]and therefore return[0.1^0, 0.1^1, 0.1^2, 0.1^3]i.e. the individual terms of f(x) without multiplication by the coefficients with \(f(x) = ax^0 + bx^1 + cx^2 +dx^3\)
- apprentice.monomial.recurrence(x, structure)
Generalised version of recurrence1D for multidimensional problems.
- Parameters
x (list) – array of parameter values
structure (numpy.ndarray) – monomial structure (>=2D)
- Returns
individual terms of f(x) WITHOUT multiplication by the coefficients
- Return type
numpy.ndarray
- apprentice.monomial.recurrence2(x, structure, nnz=slice(None, None, None))
Efficient version of recurrence with beforehand knowledge of nonzero (nnz) structure elements.
- Parameters
x (list) – array of parameter values
structure (numpy.ndarray) – monomial structure (>=2D)
nnz (list) – nonzero structure elements
- Returns
individual terms of f(x) WITHOUT multiplication by the coefficients
- Return type
numpy.ndarray
- Example
struct = monomialStructure(5,3) nnz = struct>0 recurrence2([0.1,0.2,0.3,0.4,0.5], struct, nnz)
- apprentice.monomial.vandermonde(params, order)
Construct the Vandermonde matrix for polynomial of order order
- Parameters
params (list) – array of parameter values
order (int) – order of monomial
- Returns
array of shape (len(params), numcoeffs(dim, order). The row i contains the recurrence relations obtained for param[i].
- Return type
numpy.ndarray
Function parameter space
- class apprentice.space.Space(dim, a, b, sa=None, sb=None, pnames=None)
Bases:
objectParameter space
- __init__(dim, a, b, sa=None, sb=None, pnames=None)
Initialize parameter space
- Parameters
dim (int) – dimension
a (list) – lower bound
b (list) – upper bound
sa (list) – scaled lower bound
sb (list) – scaled upper bound
pnames (list) – parameter names
- classmethod fromList(lst, pnames=None)
Create parameter space from list
- Parameters
lst (list) – list of parameter bounds
pnames (list) – parameter names
- property box
The parameter box
- Returns
lower and upper bounds of parameter box in the unscaled world
- Return type
np.array
- property box_scaled
The parameter box in scaled world
- Returns
lower and upper bounds of parameter box in the scaled world
- Return type
np.array
- property dim
Get the parameter dimension
- Returns
parameter dimension
- Return type
int
- property center
Get the center of the parameter space in unscaled world
- Returns
center of the parameter space
- Return type
list
- property center_scaled
Get the center of the parameter space in scaled world
- Returns
center of the parameter space
- Return type
list
- property jacfac
Get the jacobian factor
- Returns
jacobian factor
- Return type
list
- property as_dict
Get the parameter space as dictionary
- Returns
parameter space as dictionary
- Return type
dict
- scale(x)
Scale the point x from the observed range _Xmin, _Xmax to the interval _interval (newmax-newmin)/(oldmax-oldmin)*(x-oldmin)+newmin
- Parameters
x (np.array) – point to scale
- Returns
scaled point
- Return type
np.array
- unscale(x)
Convert a point from the scaled world back to the real world.
- Parameters
x (np.array) – point to unscale
- Returns
unscaled point
- Return type
np.array
- property pnames
Get the parameter names
- Returns
parameter names
- Return type
list
- mkSubSpace(dims)
Return a Space using only the dimensions specified in dims. Useful when fixing parameters.
- Parameters
dims (list) – dimensions
- Returns
new space (subspace)
- Return type
- static sample_main(b_min, b_max, npoints, method='uniform', seed=None)
Sample npoints dim-dimensional pointsrandomly from within this space’s bounds. Provided methods: uniform,lhs,sobol With seed=None, it is guaranteed that successive calls yield different points.
- Parameters
b_min (np.array) – bound minimum
b_max (np.array) – bound maximum
npoints (int) – number of points
method (str) – method to use
seed (int) – seed to use
- Returns
sample
- Return type
np.array
- sample(npoints, method='uniform', seed=None)
Sample in unsclaed box
- Parameters
npoints (int) – number of points
method (str) – method to use
seed (int) – seed to use
- Returns
sample in unsclaed world
- Return type
np.array
- sample_scaled(npoints, method='uniform', seed=None)
Sample in sclaed box
- Parameters
npoints (int) – number of points
method (str) – method to use
seed (int) – seed to use
- Returns
sample in sclaed world
- Return type
np.array
Utility
- class apprentice.util.Util
Bases:
objectUtility class
- static inherits_from(child, parent_name)
Check if child inherits from parent
- Parameters
child (any) – child object
parent_name (any) – parent object
- Returns
true if child inherits from parent, false otherwise
- Return type
bool
- static num_coeffs_poly(dim, order)
Find number of coefficients a dim-dimensional polynomial of order order has.
- Parameters
dim (int) – dimension value
order (int) – order of polynomial
- Returns
number of coeffecients in polynomial
- Return type
int
- static num_coeffs_to_order(dim, ncoeffs)
Infer the polynomial order from the number of coefficients and the dimension
- Parameters
dim (int) – dimension
ncoeffs (int) – number of coeffecients
- Returns
polynomial order
- Return type
int
- static num_coeffs_rapp(dim, order)
Number of coefficients a dim-dimensional rational approximation of order (m,n) has.
- Parameters
dim (int) – dimension value
order (int) – order of polynomial
- Returns
number of coeffecients in rational approximation
- Return type
int
- static gradient_recurrence(X, struct, jacfac, NNZ, sred)
Get array suitable for multiplication with coefficient vector
- Parameters
X – scaled point
struct – polynomial structure
jacfac – jacobian factor
NNZ – list of np.where results
sred – reduced structure
- Returns
array suitable for multiplication with coefficient vector
- Return type
np.array
- static prime(GREC, COEFF, dim, NNZ)
- static doubleprime(dim, xs, NSEL, HH, HNONZ, EE, COEFF)
- static hreduction(xs, ee)
- static calcSpans(spans1, DIM, G1, G2, H2, H3, grads, egrads)