Low-level Solver Interface¶
This is the low-level interface to the C++ implementation of the celerite algorithm. These methods do most of the heavy lifting but most users shouldn’t need to call these directly. This interface was built using pybind11.
-
celerite.solver.
get_kernel_value
(arg0: numpy.ndarray[float64[m, 1]], arg1: numpy.ndarray[float64[m, 1]], arg2: numpy.ndarray[float64[m, 1]], arg3: numpy.ndarray[float64[m, 1]], arg4: numpy.ndarray[float64[m, 1]], arg5: numpy.ndarray[float64[m, 1]], arg6: numpy.ndarray[float]) → object¶ Get the value of the kernel for given parameters and lags
Parameters: - alpha_real (array[j_real]) – The coefficients of the real terms.
- beta_real (array[j_real]) – The exponents of the real terms.
- alpha_complex_real (array[j_complex]) – The real part of the coefficients of the complex terms.
- alpha_complex_imag (array[j_complex]) – The imaginary part of the coefficients of the complex terms.
- beta_complex_real (array[j_complex]) – The real part of the exponents of the complex terms.
- beta_complex_imag (array[j_complex]) – The imaginary part of the exponents of the complex terms.
- tau (array[n]) – The time lags where the kernel should be evaluated.
Returns: The kernel evaluated at
tau
.Return type: array[n]
-
celerite.solver.
get_psd_value
(arg0: numpy.ndarray[float64[m, 1]], arg1: numpy.ndarray[float64[m, 1]], arg2: numpy.ndarray[float64[m, 1]], arg3: numpy.ndarray[float64[m, 1]], arg4: numpy.ndarray[float64[m, 1]], arg5: numpy.ndarray[float64[m, 1]], arg6: numpy.ndarray[float]) → object¶ Get the PSD of the kernel for given parameters and angular frequencies
Parameters: - alpha_real (array[j_real]) – The coefficients of the real terms.
- beta_real (array[j_real]) – The exponents of the real terms.
- alpha_complex_real (array[j_complex]) – The real part of the coefficients of the complex terms.
- alpha_complex_imag (array[j_complex]) – The imaginary part of the coefficients of the complex terms.
- beta_complex_real (array[j_complex]) – The real part of the exponents of the complex terms.
- beta_complex_imag (array[j_complex]) – The imaginary part of the exponents of the complex terms.
- omega (array[n]) – The frequencies where the PSD should be evaluated.
Returns: The PSD evaluated at
omega
.Return type: array[n]
-
celerite.solver.
check_coefficients
(arg0: numpy.ndarray[float64[m, 1]], arg1: numpy.ndarray[float64[m, 1]], arg2: numpy.ndarray[float64[m, 1]], arg3: numpy.ndarray[float64[m, 1]], arg4: numpy.ndarray[float64[m, 1]], arg5: numpy.ndarray[float64[m, 1]]) → bool¶ Apply Sturm’s theorem to check if parameters yield a positive PSD
Parameters: - alpha_real (array[j_real]) – The coefficients of the real terms.
- beta_real (array[j_real]) – The exponents of the real terms.
- alpha_complex_real (array[j_complex]) – The real part of the coefficients of the complex terms.
- alpha_complex_imag (array[j_complex]) – The imaginary part of the coefficients of the complex terms.
- beta_complex_real (array[j_complex]) – The real part of the exponents of the complex terms.
- beta_complex_imag (array[j_complex]) – The imaginary part of the exponents of the complex terms.
Returns: True
if the PSD is everywhere positive.Return type:
-
class
celerite.solver.
Solver
¶ A thin wrapper around the C++ BandSolver class
The class provides all of the computation power for the
celerite
module. The key methods are listed below but thesolver.Solver.compute()
method must always be called first.-
compute
(self: celerite.solver.Solver, arg0: numpy.ndarray[float64[m, 1]], arg1: numpy.ndarray[float64[m, 1]], arg2: numpy.ndarray[float64[m, 1]], arg3: numpy.ndarray[float64[m, 1]], arg4: numpy.ndarray[float64[m, 1]], arg5: numpy.ndarray[float64[m, 1]], arg6: numpy.ndarray[float64[m, 1]], arg7: numpy.ndarray[float64[m, 1]]) → int¶ Assemble the extended matrix and perform the banded LU decomposition
Parameters: - alpha_real (array[j_real]) – The coefficients of the real terms.
- beta_real (array[j_real]) – The exponents of the real terms.
- alpha_complex_real (array[j_complex]) – The real part of the coefficients of the complex terms.
- alpha_complex_imag (array[j_complex]) – The imaginary part of the coefficients of the complex terms.
- beta_complex_real (array[j_complex]) – The real part of the exponents of the complex terms.
- beta_complex_imag (array[j_complex]) – The imaginary part of the exponents of the complex terms.
- x (array[n]) – The _sorted_ array of input coordinates.
- diag (array[n]) – An array that should be added to the diagonal of the
matrix. This often corresponds to measurement uncertainties and in
that case,
diag
should be the measurement _variance_ (i.e. sigma^2).
Returns: 1
if the dimensions are inconsistent and0
otherwise. No attempt is made to confirm that the matrix is positive definite. If it is not positive definite, thesolve
andlog_determinant
methods will return incorrect results.Return type:
-
computed
(self: celerite.solver.Solver) → bool¶ A flag that indicates if
compute
has been executedReturns: True
ifsolver.Solver.compute()
was previously executed successfully.Return type: bool
-
dot
(self: celerite.solver.Solver, arg0: numpy.ndarray[float64[m, 1]], arg1: numpy.ndarray[float64[m, 1]], arg2: numpy.ndarray[float64[m, 1]], arg3: numpy.ndarray[float64[m, 1]], arg4: numpy.ndarray[float64[m, 1]], arg5: numpy.ndarray[float64[m, 1]], arg6: numpy.ndarray[float64[m, 1]], arg7: numpy.ndarray[float64[m, n]]) → numpy.ndarray[float64[m, n]]¶ Compute the dot product of a
celerite
matrix and another arbitrary matrixThis method computes
A.b
whereA
is defined by the parameters andb
is an arbitrary matrix of the correct shape.Parameters: - alpha_real (array[j_real]) – The coefficients of the real terms.
- beta_real (array[j_real]) – The exponents of the real terms.
- alpha_complex_real (array[j_complex]) – The real part of the coefficients of the complex terms.
- alpha_complex_imag (array[j_complex]) – The imaginary part of the coefficients of the complex terms.
- beta_complex_real (array[j_complex]) – The real part of the exponents of the complex terms.
- beta_complex_imag (array[j_complex]) – The imaginary part of the exponents of the complex terms.
- x (array[n]) – The _sorted_ array of input coordinates.
- b (array[n] or array[n, neq]) – The matrix
b
described above.
Returns: The dot product
A.b
as described above.Return type: Raises: ValueError
– For mismatched dimensions.
-
dot_solve
(self: celerite.solver.Solver, arg0: numpy.ndarray[float64[m, n]]) → float¶ Solve the system
b^T . A^-1 . b
A previous call to
solver.Solver.compute()
defines a matrixA
and this method solvesb^T . A^-1 . b
for a vectorb
.Parameters: b (array[n]) – The right hand side of the linear system. Returns: The solution of b^T . A^-1 . b
.Return type: float Raises: ValueError
– For mismatched dimensions.
-
log_determinant
(self: celerite.solver.Solver) → float¶ Get the log-determinant of the matrix defined by
compute
Returns: The log-determinant of the matrix defined by solver.Solver.compute()
.Return type: float
-
solve
(self: celerite.solver.Solver, arg0: numpy.ndarray[float64[m, n]]) → numpy.ndarray[float64[m, n]]¶ Solve a linear system for the matrix defined in
compute
A previous call to
solver.Solver.compute()
defines a matrixA
and this method solves forx
in the matrix equationA.x = b
.Parameters: b (array[n] or array[n, nrhs]) – The right hand side of the linear system. Returns: The solution of the linear system. Return type: array[n] or array[n, nrhs] Raises: ValueError
– For mismatched dimensions.
-
use_lapack
(self: celerite.solver.Solver) → bool¶
-