Class ChebyshevApproximation
- Namespace
- ChebyshevSharp
- Assembly
- ChebyshevSharp.dll
Multi-dimensional Chebyshev tensor interpolation with analytical derivatives. Uses barycentric interpolation with pre-computed weights.
public class ChebyshevApproximation
- Inheritance
-
ChebyshevApproximation
- Inherited Members
Constructors
ChebyshevApproximation(Func<double[], object?, double>, int, double[][], int[], int)
Create a new ChebyshevApproximation.
public ChebyshevApproximation(Func<double[], object?, double> function, int numDimensions, double[][] domain, int[] nNodes, int maxDerivativeOrder = 2)
Parameters
functionFunc<double[], object, double>Function to approximate: f(point, data) -> double.
numDimensionsintNumber of input dimensions.
domaindouble[][]Bounds for each dimension as double[ndim][2].
nNodesint[]Number of Chebyshev nodes per dimension.
maxDerivativeOrderintMaximum derivative order to support (default 2).
Properties
BuildTime
Time taken by Build() in seconds.
public double BuildTime { get; }
Property Value
DiffMatrices
Spectral differentiation matrices per dimension.
public double[][,]? DiffMatrices { get; }
Property Value
- double[][,]
Domain
Domain bounds for each dimension, as list of [lo, hi].
public double[][] Domain { get; }
Property Value
- double[][]
Function
The function to approximate. Null after load or from_values.
public Func<double[], object?, double>? Function { get; }
Property Value
MaxDerivativeOrder
Maximum supported derivative order.
public int MaxDerivativeOrder { get; }
Property Value
NEvaluations
Number of function evaluations during Build().
public int NEvaluations { get; }
Property Value
NNodes
Number of Chebyshev nodes per dimension.
public int[] NNodes { get; }
Property Value
- int[]
NodeArrays
Chebyshev nodes per dimension, each sorted ascending.
public double[][] NodeArrays { get; }
Property Value
- double[][]
NumDimensions
Number of input dimensions.
public int NumDimensions { get; }
Property Value
TensorValues
Flat tensor of function values at all node combinations (C-order).
public double[]? TensorValues { get; }
Property Value
- double[]
Weights
Barycentric weights per dimension.
public double[][]? Weights { get; }
Property Value
- double[][]
Methods
Build(bool)
Evaluate the function at all node combinations and pre-compute weights.
public void Build(bool verbose = true)
Parameters
verboseboolIf true, print build progress.
ChebyshevCoefficients1D(double[])
Get Chebyshev coefficients for a 1D array of values at Type I nodes. Public for testing.
public static double[] ChebyshevCoefficients1D(double[] values)
Parameters
valuesdouble[]
Returns
- double[]
ErrorEstimate()
Estimate the supremum-norm interpolation error using Chebyshev coefficient decay.
public double ErrorEstimate()
Returns
- double
Estimated maximum interpolation error.
Eval(double[], int[])
Evaluate using dimensional decomposition with barycentric interpolation. Loop-based implementation matching Python eval().
public double Eval(double[] point, int[] derivativeOrder)
Parameters
pointdouble[]Query point, one coordinate per dimension.
derivativeOrderint[]Derivative order per dimension.
Returns
- double
Interpolated value or derivative at the query point.
Extrude(params (int dimIndex, double[] bounds, int nNodes)[])
Add new dimensions where the function is constant.
public ChebyshevApproximation Extrude(params (int dimIndex, double[] bounds, int nNodes)[] extrudeParams)
Parameters
Returns
FromValues(double[], int, double[][], int[], int)
Create an interpolant from pre-computed function values.
public static ChebyshevApproximation FromValues(double[] tensorValues, int numDimensions, double[][] domain, int[] nNodes, int maxDerivativeOrder = 2)
Parameters
Returns
GetDerivativeId(int[])
Return derivative order as-is (for API compatibility).
public int[] GetDerivativeId(int[] derivativeOrder)
Parameters
derivativeOrderint[]
Returns
- int[]
Integrate(int[]?, (double lo, double hi)[]?)
Integrate the interpolant over one or more dimensions.
public object Integrate(int[]? dims = null, (double lo, double hi)[]? bounds = null)
Parameters
dimsint[]Dimensions to integrate out. Null = all.
bounds(double lo, double hi)[]Sub-interval bounds per dim. Null = full domain.
Returns
- object
Scalar if all dims integrated, otherwise a lower-dimensional interpolant.
Load(string)
Load a previously saved interpolant from a file.
public static ChebyshevApproximation Load(string path)
Parameters
pathstringPath to the saved file.
Returns
- ChebyshevApproximation
The restored interpolant.
Maximize(int?, Dictionary<int, double>?)
Find the maximum value of the interpolant along a dimension.
public (double value, double location) Maximize(int? dim = null, Dictionary<int, double>? fixedDims = null)
Parameters
dimint?fixedDimsDictionary<int, double>
Returns
Minimize(int?, Dictionary<int, double>?)
Find the minimum value of the interpolant along a dimension.
public (double value, double location) Minimize(int? dim = null, Dictionary<int, double>? fixedDims = null)
Parameters
dimint?fixedDimsDictionary<int, double>
Returns
Nodes(int, double[][], int[])
Generate Chebyshev nodes without evaluating any function.
public static NodeInfo Nodes(int numDimensions, double[][] domain, int[] nNodes)
Parameters
numDimensionsintNumber of dimensions.
domaindouble[][]Lower and upper bounds for each dimension.
nNodesint[]Number of Chebyshev nodes per dimension.
Returns
- NodeInfo
Dictionary with "NodesPerDim", "FullGrid", and "Shape".
Roots(int?, Dictionary<int, double>?)
Find all roots of the interpolant along a specified dimension.
public double[] Roots(int? dim = null, Dictionary<int, double>? fixedDims = null)
Parameters
dimint?fixedDimsDictionary<int, double>
Returns
- double[]
Save(string)
Save the built interpolant to a file using JSON serialization.
public void Save(string path)
Parameters
pathstringDestination file path.
Slice(params (int dimIndex, double value)[])
Fix one or more dimensions at given values, reducing dimensionality.
public ChebyshevApproximation Slice(params (int dimIndex, double value)[] sliceParams)
Parameters
Returns
ToReprString()
Compact representation of the interpolant.
public string ToReprString()
Returns
ToString()
Returns a string that represents the current object.
public override string ToString()
Returns
- string
A string that represents the current object.
VectorizedEval(double[], int[])
Fully vectorized evaluation using matrix operations. Replaces the Python loop with BLAS-style matrix-vector products.
public double VectorizedEval(double[] point, int[] derivativeOrder)
Parameters
pointdouble[]Query point, one coordinate per dimension.
derivativeOrderint[]Derivative order per dimension.
Returns
- double
Interpolated value or derivative.
VectorizedEvalBatch(double[][], int[])
Evaluate at multiple points.
public double[] VectorizedEvalBatch(double[][] points, int[] derivativeOrder)
Parameters
pointsdouble[][]Points as double[N][numDimensions].
derivativeOrderint[]Derivative order per dimension.
Returns
- double[]
Results array of length N.
VectorizedEvalMulti(double[], int[][])
Evaluate multiple derivative orders at the same point, sharing barycentric weights.
public double[] VectorizedEvalMulti(double[] point, int[][] derivativeOrders)
Parameters
pointdouble[]Query point.
derivativeOrdersint[][]Each inner array specifies derivative order per dimension.
Returns
- double[]
One result per derivative order.
Operators
operator +(ChebyshevApproximation, ChebyshevApproximation)
Add two interpolants with the same grid.
public static ChebyshevApproximation operator +(ChebyshevApproximation a, ChebyshevApproximation b)
Parameters
Returns
operator /(ChebyshevApproximation, double)
Divide interpolant by a scalar.
public static ChebyshevApproximation operator /(ChebyshevApproximation a, double scalar)
Parameters
aChebyshevApproximationscalardouble
Returns
operator *(ChebyshevApproximation, double)
Multiply interpolant by a scalar.
public static ChebyshevApproximation operator *(ChebyshevApproximation a, double scalar)
Parameters
aChebyshevApproximationscalardouble
Returns
operator *(double, ChebyshevApproximation)
Multiply scalar by interpolant.
public static ChebyshevApproximation operator *(double scalar, ChebyshevApproximation a)
Parameters
scalardoubleaChebyshevApproximation
Returns
operator -(ChebyshevApproximation, ChebyshevApproximation)
Subtract two interpolants with the same grid.
public static ChebyshevApproximation operator -(ChebyshevApproximation a, ChebyshevApproximation b)
Parameters
Returns
operator -(ChebyshevApproximation)
Negate interpolant.
public static ChebyshevApproximation operator -(ChebyshevApproximation a)