Pre-computed Values
Motivation -- Nodes First, Values Later
In production environments -- HPC clusters, distributed pricing engines, cloud compute -- the function to be interpolated often cannot be called from within C#. The evaluation may run on a separate process, a remote node, or inside a proprietary pricing library.
ChebyshevSharp's standard workflow (constructor + Build()) requires a C#
delegate. The pre-computed values workflow decouples node generation from
function evaluation:
- Generate nodes -- call
Nodes()to get the Chebyshev grid points. - Evaluate externally -- feed the points to your own pricing engine.
- Load values -- call
FromValues()to construct the interpolant.
The resulting object is usable for stored-tensor operations: evaluation,
derivatives, integration, rootfinding, optimization, algebra,
extrusion/slicing, serialization, and error estimation. It cannot call
Build() again because it does not retain a callable function.
Workflow
ChebyshevApproximation
using ChebyshevSharp;
// Step 1: Get nodes (no function needed)
NodeInfo info = ChebyshevApproximation.Nodes(
numDimensions: 2,
domain: new[] { new[] { -1.0, 1.0 }, new[] { 0.0, 2.0 } },
nNodes: new[] { 15, 11 }
);
// info.NodesPerDim -> double[2][] (15 nodes, 11 nodes)
// info.FullGrid -> double[165][] (each row is one [x0, x1] point)
// info.Shape -> int[] { 15, 11 }
// Step 2: Evaluate externally
// Each row of FullGrid is one (x0, x1) point.
// Feed these to your pricing engine, collect results.
double[] values = new double[info.FullGrid.Length];
for (int i = 0; i < info.FullGrid.Length; i++)
values[i] = MyPricingEngine(info.FullGrid[i]);
// Step 3: Build from values
var cheb = ChebyshevApproximation.FromValues(
values,
numDimensions: 2,
domain: new[] { new[] { -1.0, 1.0 }, new[] { 0.0, 2.0 } },
nNodes: new[] { 15, 11 }
);
// Now use it like any other interpolant
double price = cheb.VectorizedEval(new[] { 0.5, 1.0 }, new[] { 0, 0 });
double delta = cheb.VectorizedEval(new[] { 0.5, 1.0 }, new[] { 1, 0 });
double integral = (double)cheb.Integrate();
ChebyshevSpline
For piecewise interpolation with knots, the same pattern applies per piece:
using ChebyshevSharp;
// Step 1: Get per-piece nodes
SplineNodeInfo info = ChebyshevSpline.Nodes(
numDimensions: 1,
domain: new[] { new[] { -1.0, 1.0 } },
nNodes: new[] { 15 },
knots: new[] { new[] { 0.0 } }
);
// info.Pieces -> SplinePieceNodeInfo[] (2 pieces)
// info.NumPieces -> 2
// info.PieceShape -> int[] { 2 }
// Step 2: Evaluate each piece externally
double[][] pieceValues = new double[info.NumPieces][];
for (int p = 0; p < info.NumPieces; p++)
{
SplinePieceNodeInfo piece = info.Pieces[p];
pieceValues[p] = new double[piece.FullGrid.Length];
for (int i = 0; i < piece.FullGrid.Length; i++)
pieceValues[p][i] = MyPricingEngine(piece.FullGrid[i]);
}
// Step 3: Build
var spline = ChebyshevSpline.FromValues(
pieceValues,
numDimensions: 1,
domain: new[] { new[] { -1.0, 1.0 } },
nNodes: new[] { 15 },
knots: new[] { new[] { 0.0 } }
);
ChebyshevTT
ChebyshevTT.FromValues() compresses a dense value tensor with TT-SVD. It is
useful when another system has already produced a full Chebyshev grid and you
want TT storage afterward. It is not a sparse alternative to TT-Cross; the input
array must still contain every row-major grid value.
using ChebyshevSharp;
var (nodes, shape) = ChebyshevTT.Nodes(2,
new[] { new[] { -1.0, 1.0 }, new[] { -1.0, 1.0 } },
new[] { 8, 8 });
double[] values = new double[shape[0] * shape[1]];
for (int i = 0; i < shape[0]; i++)
for (int j = 0; j < shape[1]; j++)
values[i * shape[1] + j] = Math.Sin(nodes[0][i]) * Math.Cos(nodes[1][j]);
var tt = ChebyshevTT.FromValues(
values,
numDimensions: 2,
domain: new[] { new[] { -1.0, 1.0 }, new[] { -1.0, 1.0 } },
nNodes: shape,
maxRank: 6,
tolerance: 1e-10);
For high-dimensional models, prefer new ChebyshevTT(...).Build(method: "cross")
unless the dense tensor is intentionally small enough to materialize.
Indexing Convention
The tensor entries must satisfy:
For ChebyshevApproximation.Nodes(), the rows of FullGrid follow row-major
(C-order) enumeration, so iterating over FullGrid in order and storing the
results in a flat double[] produces the correct tensor layout automatically.
For ChebyshevTT.Nodes(), construct the same row-major order explicitly by
making the last index vary fastest, as shown in the TT example above.
Warning: Use row-major (C-order) layout. The flat
double[]passed toFromValues()must be in row-major order (last index varies fastest). Do not use column-major (Fortran-order) layout, which would silently produce incorrect tensor entries and lead to wrong interpolation results.
For spline pieces, the list order follows row-major multi-index enumeration
over PieceShape. In 2D with knots = { { 0.0 }, { 1.0 } } on domain
[[-1,1], [0,2]], the four pieces are:
| Index | PieceIndex | SubDomain |
|---|---|---|
| 0 | (0, 0) | [(-1, 0), (0, 1)] |
| 1 | (0, 1) | [(-1, 0), (1, 2)] |
| 2 | (1, 0) | [(0, 1), (0, 1)] |
| 3 | (1, 1) | [(0, 1), (1, 2)] |
Mathematical Justification
Everything the interpolant needs -- except the function values themselves -- depends only on the node positions:
| Pre-computed data | Depends on |
|---|---|
| Chebyshev nodes | domain, nNodes |
| Barycentric weights | nodes |
| Differentiation matrices | nodes, weights |
| Fejer quadrature weights | nodes |
The function values appear only in TensorValues. The Chebyshev nodes are
the zeros of \(T_n(x)\), mapped to each dimension's domain (Trefethen 2013,
Ch. 3). Since FromValues() computes all numerical pre-computation from the
same node formula as Build(), evaluation and derivative results match a
Build()-based interpolant from the same values. Runtime metadata differs:
Function is null, and build timing or evaluation-count metadata may differ.
Examples
1-D: sin(x) on [0, pi]
using ChebyshevSharp;
NodeInfo info = ChebyshevApproximation.Nodes(1,
new[] { new[] { 0.0, Math.PI } }, new[] { 20 });
double[] values = new double[info.FullGrid.Length];
for (int i = 0; i < info.FullGrid.Length; i++)
values[i] = Math.Sin(info.FullGrid[i][0]);
var cheb = ChebyshevApproximation.FromValues(values, 1,
new[] { new[] { 0.0, Math.PI } }, new[] { 20 });
// Integral of sin(x) on [0, pi] ~ 2.0
double integral = (double)cheb.Integrate();
Console.WriteLine(integral); // 1.9999999...
Combine with Algebra
// Two proxies built externally
var chebA = ChebyshevApproximation.FromValues(
valsA, 2, domain, nNodes);
var chebB = ChebyshevApproximation.FromValues(
valsB, 2, domain, nNodes);
// Portfolio-level proxy
var portfolio = 0.6 * chebA + 0.4 * chebB;
Save and Load
cheb.Save("my_proxy.json");
var loaded = ChebyshevApproximation.Load("my_proxy.json");
Error Handling
| Error | Raised when |
|---|---|
ArgumentException (shape) |
tensorValues.Length does not equal the product of nNodes |
ArgumentException (NaN/Inf) |
tensorValues contains non-finite values |
ArgumentException (dimensions) |
domain.Length != numDimensions or nNodes.Length != numDimensions |
ArgumentException (domain) |
lo >= hi for any dimension |
InvalidOperationException (build) |
Calling Build() on a FromValues object (Function is null) |
Calling
Build()on aFromValuesresult: Objects created viaFromValues()haveFunction = null. To re-build with a different set of values, create a new object viaFromValues(). To re-build from a callable, create a newChebyshevApproximationorChebyshevSpline, orChebyshevTTwith that function and callBuild()on the new object.
Combining with Other Features
Objects created via FromValues() support the stored-interpolant operations for
their class:
- Evaluation & derivatives -- for example
VectorizedEval()/VectorizedEvalMulti()on dense interpolants, orEval()/EvalMulti()onChebyshevTT - Calculus --
Integrate(),Roots(),Minimize(),Maximize() - Algebra --
+,-,*,/(including withBuild()-based objects) - Extrusion & slicing --
Extrude(),Slice() - Serialization --
Save(),Load() - Error estimation --
ErrorEstimate()
Methods that require the original callable function, such as Build() or TT
RunCompletion(), are not available on FromValues() results.
Note:
ChebyshevSliderdoes not supportNodes()orFromValues(). The sliding technique builds a separate low-dimensional interpolant for each partition group, each requiring a different subset of dimensions evaluated at pivot values for the remaining dimensions. There is no single "nodes first, values later" grid that covers this decomposition. To use pre-computed values with a slider-like workflow, build aChebyshevApproximationfrom values for each group independently.
API Reference
ChebyshevApproximation.Nodes()
| Parameter | Type | Description |
|---|---|---|
numDimensions |
int |
Number of dimensions |
domain |
double[][] |
Bounds per dimension (new[] { lo, hi } each) |
nNodes |
int[] |
Nodes per dimension |
Returns NodeInfo with properties:
| Property | Type | Description |
|---|---|---|
NodesPerDim |
double[][] |
Chebyshev nodes for each dimension, sorted ascending |
FullGrid |
double[][] |
Full Cartesian product grid, shape (totalPoints, numDimensions) |
Shape |
int[] |
Expected tensor shape (equal to nNodes) |
ChebyshevApproximation.FromValues()
| Parameter | Type | Description |
|---|---|---|
tensorValues |
double[] |
Function values, flat array of length Product(nNodes) |
numDimensions |
int |
Number of dimensions |
domain |
double[][] |
Bounds per dimension |
nNodes |
int[] |
Nodes per dimension |
maxDerivativeOrder |
int (default 2) |
Maximum derivative order |
Returns ChebyshevApproximation with Function = null.
ChebyshevSpline.Nodes()
| Parameter | Type | Description |
|---|---|---|
numDimensions |
int |
Number of dimensions |
domain |
double[][] |
Bounds per dimension |
nNodes |
int[] |
Nodes per dimension per piece |
knots |
double[][] |
Knot positions per dimension |
Returns SplineNodeInfo with properties:
| Property | Type | Description |
|---|---|---|
Pieces |
SplinePieceNodeInfo[] |
Per-piece node info |
NumPieces |
int |
Total number of pieces |
PieceShape |
int[] |
Per-dimension piece counts |
Each SplinePieceNodeInfo contains:
| Property | Type | Description |
|---|---|---|
PieceIndex |
int[] |
Multi-index of this piece |
SubDomain |
double[][] |
Sub-domain bounds for this piece |
NodesPerDim |
double[][] |
Per-dimension node arrays |
FullGrid |
double[][] |
Full Cartesian product grid |
Shape |
int[] |
Tensor shape |
ChebyshevSpline.FromValues()
| Parameter | Type | Description |
|---|---|---|
pieceValues |
double[][] |
Per-piece values in row-major order |
numDimensions |
int |
Number of dimensions |
domain |
double[][] |
Bounds per dimension |
nNodes |
int[] |
Nodes per dimension per piece |
knots |
double[][] |
Knot positions per dimension |
maxDerivativeOrder |
int (default 2) |
Maximum derivative order |
Returns ChebyshevSpline with Function = null.
ChebyshevTT.Nodes()
| Parameter | Type | Description |
|---|---|---|
numDimensions |
int |
Number of dimensions |
domain |
double[][] |
Bounds per dimension |
nNodes |
int[] |
Nodes per dimension |
Returns (double[][] nodesPerDim, int[] shape).
ChebyshevTT.FromValues()
| Parameter | Type | Description |
|---|---|---|
tensorValues |
double[] |
Dense row-major values, length Product(nNodes) |
numDimensions |
int |
Number of dimensions |
domain |
double[][] |
Bounds per dimension |
nNodes |
int[] |
Nodes per dimension |
maxRank |
int (default 10) |
Maximum TT rank |
tolerance |
double (default 1e-6) |
Non-negative TT-SVD truncation tolerance |
Returns ChebyshevTT with no callable function and Method = "svd".
See Also
- Advanced Usage -- arithmetic operators and algebra
- Calculus -- integrate, find roots, optimize
- Extrusion & Slicing -- add/fix dimensions
- Serialization & Construction -- persist FromValues objects
References
- Trefethen, L. N. (2013). Approximation Theory and Approximation Practice. SIAM. Chapter 3.