Construct 1D interpolator with extrapolation config.
#pragma once
#include <algorithm>
#include <casadi/casadi.hpp>
#include <cmath>
#include <optional>
#include <vector>
};
};
struct ExtrapolationConfig {
static ExtrapolationConfig
clamp() {
return {}; }
static ExtrapolationConfig
linear(std::optional<double> lower = std::nullopt,
std::optional<double> upper = std::nullopt) {
}
};
namespace detail {
switch (method) {
return "linear";
return "bspline";
return "linear";
default:
return "linear";
}
}
return "Interpolator: Hermite/Catmull-Rom is numeric-only because interval selection "
"requires runtime comparisons. Use BSpline for symbolic or optimization queries.";
}
template <typename Derived>
using Scalar = typename Derived::Scalar;
Eigen::Index idx = 0;
for (Eigen::Index j = 0; j < values.cols(); ++j) {
for (Eigen::Index i = 0; i < values.rows(); ++i) {
flattened(idx++) = values(i, j);
}
}
return flattened;
}
switch (method) {
return "interpn: symbolic table values are not supported for Hermite/Catmull-Rom "
"because interval selection remains numeric-only. Use BSpline for symbolic or "
"optimization queries.";
return "interpn: symbolic table values are not supported for Nearest because the "
"method is non-differentiable. Use Linear or BSpline.";
default:
return "interpn: symbolic table values are only supported for Linear and BSpline.";
}
}
struct InterpnGridData {
std::vector<std::vector<double>>
grid;
};
const char *context) {
if (points.empty()) {
throw InterpolationError(std::string(context) + ": points cannot be empty");
}
const int n_dims = static_cast<int>(points.size());
InterpnGridData data;
data.grid.resize(n_dims);
for (int d = 0; d < n_dims; ++d) {
if (points[d].size() < 2) {
throw InterpolationError(std::string(context) + ": points[" + std::to_string(d) +
"] must have at least 2 points");
}
data.grid[d].resize(points[d].size());
Eigen::VectorXd::Map(data.grid[d].data(), points[d].size()) = points[d];
for (size_t i = 0; i + 1 < data.grid[d].size(); ++i) {
if (data.grid[d][i + 1] <= data.grid[d][i]) {
throw InterpolationError(std::string(context) + ": points[" + std::to_string(d) +
"] must be strictly increasing");
}
}
data.expected_size *= static_cast<int>(points[d].size());
}
return data;
}
const char *context) {
for (size_t d = 0; d < grid.size(); ++d) {
if (grid[d].size() < 4) {
throw InterpolationError(std::string(context) +
": BSpline requires at least 4 grid points in dimension " +
std::to_string(d));
}
}
}
template <JanusScalar Scalar>
const char *context) {
if (xi.cols() != n_dims && xi.rows() != n_dims) {
throw InterpolationError(std::string(context) +
": xi must have shape (n_points, n_dims) or (n_dims, n_points)");
}
const bool need_transpose = (xi.rows() == n_dims && xi.cols() != n_dims);
if (need_transpose) {
return xi.transpose().eval();
}
return xi.eval();
}
template <JanusScalar Scalar>
const std::vector<std::vector<double>> &grid) {
for (int d = 0; d < point.size(); ++d) {
std::is_same_v<Scalar, SymbolicScalar> ? point(d) :
SymbolicScalar(point(d));
value = SymbolicScalar::fmax(value, grid[d].front());
value = SymbolicScalar::fmin(value, grid[d].back());
clamped(d) = value;
}
return clamped;
}
template <JanusScalar Scalar>
const std::vector<std::vector<double>> &grid) {
static_assert(std::is_floating_point_v<Scalar>);
for (int d = 0; d < point.size(); ++d) {
if (point(d) < grid[d].front() || point(d) > grid[d].back()) {
return true;
}
}
return false;
}
const std::vector<std::vector<double>> &grid) {
for (int d = 0; d < point.size(); ++d) {
out_of_bounds = SymbolicScalar::logic_or(out_of_bounds, point(d) < grid[d].front());
out_of_bounds = SymbolicScalar::logic_or(out_of_bounds, point(d) > grid[d].back());
}
return out_of_bounds;
}
}
}
casadi::Dict opts;
opts["inline"] = true;
}
inline double catmull_rom_slope(
const std::vector<double> &x,
const std::vector<double> &y,
size_t i) {
size_t n = x.size();
if (n < 2)
return 0.0;
if (i == 0) {
return (y[1] - y[0]) / (x[1] - x[0]);
} else if (i == n - 1) {
return (y[n - 1] - y[n - 2]) / (x[n - 1] - x[n - 2]);
} else {
return (y[i + 1] - y[i - 1]) / (x[i + 1] - x[i - 1]);
}
}
double query) {
size_t n = x.size();
if (n < 2) {
return y[0];
}
query = std::max(query, x.front());
query = std::min(query, x.back());
auto it = std::upper_bound(x.begin(), x.end(), query);
size_t idx = 0;
if (it == x.begin()) {
idx = 0;
} else if (it == x.end()) {
idx = n - 2;
} else {
idx = static_cast<size_t>(std::distance(x.begin(), it)) - 1;
}
if (idx >= n - 1)
idx = n - 2;
double x0 = x[idx];
double x1 = x[idx + 1];
double y0 = y[idx];
double y1 = y[idx + 1];
double h = x1 - x0;
double t = (query - x0) / h;
t = std::max(0.0, std::min(1.0, t));
double t2 = t * t;
double t3 = t2 * t;
double h00 = 2.0 * t3 - 3.0 * t2 + 1.0;
double h10 = t3 - 2.0 * t2 + t;
double h01 = -2.0 * t3 + 3.0 * t2;
double h11 = t3 - t2;
return h00 * y0 + h10 * (m0 * h) + h01 * y1 + h11 * (m1 * h);
}
std::vector<int> multi(dims.size());
int remaining = flat_idx;
for (size_t d = 0; d < dims.size(); ++d) {
multi[d] = remaining % dims[d];
remaining /= dims[d];
}
return multi;
}
int flat = 0;
int stride = 1;
for (size_t d = 0; d < dims.size(); ++d) {
flat += multi[d] * stride;
stride *= dims[d];
}
return flat;
}
const std::vector<double> &values,
const std::vector<double> &query) {
int n_dims = static_cast<int>(grid.size());
std::vector<int> dims(n_dims);
for (int d = 0; d < n_dims; ++d) {
dims[d] = static_cast<int>(grid[d].size());
}
std::vector<double> current_values = values;
std::vector<int> current_dims = dims;
for (int d = 0; d < n_dims; ++d) {
const std::vector<double> &x = grid[d];
int n_x = current_dims[0];
int remaining_size = 1;
for (size_t dd = 1; dd < current_dims.size(); ++dd) {
remaining_size *= current_dims[dd];
}
std::vector<double> new_values(remaining_size);
std::vector<double> y_slice(n_x);
for (int i = 0; i < n_x; ++i) {
int flat_idx = i +
slice * n_x;
y_slice[i] = current_values[flat_idx];
}
}
current_values = std::move(new_values);
current_dims.erase(current_dims.begin());
}
return current_values[0];
}
}
class Interpolator {
private:
std::vector<std::vector<double>> m_grid;
std::vector<double> m_values;
casadi::Function m_casadi_fn;
int m_dims = 0;
bool m_valid = false;
bool m_use_custom_hermite = false;
ExtrapolationConfig m_extrap_config;
std::vector<double> m_slopes_left;
std::vector<double> m_slopes_right;
public:
if (points.empty()) {
throw InterpolationError("Interpolator: points cannot be empty");
}
m_dims = static_cast<int>(points.size());
m_grid.resize(m_dims);
int expected_size = 1;
for (int d = 0; d < m_dims; ++d) {
if (points[d].size() < 2) {
throw InterpolationError("Interpolator: need at least 2 grid points in dimension " +
std::to_string(d));
}
m_grid[d].resize(points[d].size());
Eigen::VectorXd::Map(m_grid[d].data(), points[d].size()) = points[d];
if (!std::is_sorted(m_grid[d].begin(), m_grid[d].end())) {
throw InterpolationError("Interpolator: points[" + std::to_string(d) +
"] must be sorted");
}
expected_size *= static_cast<int>(points[d].size());
}
if (values.size() != expected_size) {
throw InterpolationError("Interpolator: values size mismatch. Expected " +
std::to_string(expected_size) + ", got " +
std::to_string(values.size()));
}
m_values.resize(values.size());
Eigen::VectorXd::Map(m_values.data(), values.size()) = values;
}
if (x.size() != y.size()) {
throw InterpolationError("Interpolator: x and y must have same size");
}
if (x.size() < 2) {
throw InterpolationError("Interpolator: need at least 2 grid points");
}
m_dims = 1;
m_grid.resize(1);
m_grid[0].resize(x.size());
Eigen::VectorXd::Map(m_grid[0].data(), x.size()) = x;
if (!std::is_sorted(m_grid[0].begin(), m_grid[0].end())) {
throw InterpolationError("Interpolator: x grid must be sorted");
}
m_values.resize(y.size());
Eigen::VectorXd::Map(m_values.data(), y.size()) = y;
}
m_extrap_config = extrap;
compute_boundary_slopes();
}
ExtrapolationConfig extrap)
m_extrap_config = extrap;
compute_boundary_slopes();
}
int dims()
const {
return m_dims; }
bool valid()
const {
return m_valid; }
template <JanusScalar Scalar> Scalar
operator()(
const Scalar &query)
const {
if (!m_valid)
throw InterpolationError("Interpolator: not initialized");
if (m_dims != 1)
throw InterpolationError("Interpolator: scalar query only valid for 1D");
if constexpr (std::is_floating_point_v<Scalar>) {
return eval_numeric_scalar(query);
} else {
if (m_use_custom_hermite) {
}
return eval_symbolic_scalar(query);
}
}
if (!m_valid)
throw InterpolationError("Interpolator: not initialized");
if (query.size() != m_dims) {
throw InterpolationError("Interpolator: query size must match dims. Expected " +
std::to_string(m_dims) + ", got " +
std::to_string(query.size()));
}
if constexpr (std::is_floating_point_v<Scalar>) {
return eval_numeric_point(query);
} else {
if (m_use_custom_hermite) {
}
return eval_symbolic_point(query);
}
}
template <typename Derived>
auto operator()(
const Eigen::MatrixBase<Derived> &queries)
const
using Scalar = typename Derived::Scalar;
if (!m_valid)
throw InterpolationError("Interpolator: not initialized");
int n_points;
if (m_dims == 1) {
n_points = static_cast<int>(queries.size());
} else {
bool is_transposed = (queries.rows() == m_dims && queries.cols() != m_dims);
n_points =
is_transposed ? static_cast<int>(queries.cols()) : static_cast<int>(queries.rows());
}
if constexpr (std::is_floating_point_v<Scalar>) {
if (m_dims == 1) {
for (int i = 0; i < n_points; ++i) {
Scalar val;
if (queries.cols() == 1) {
val = queries(i, 0);
} else {
val = queries(0, i);
}
result(i) = eval_numeric_scalar(val);
}
} else {
bool is_transposed = (queries.rows() == m_dims && queries.cols() != m_dims);
for (int i = 0; i < n_points; ++i) {
for (int d = 0; d < m_dims; ++d) {
point(d) = is_transposed ? queries(d, i) : queries(i, d);
}
result(i) = eval_numeric_point(point);
}
}
} else {
if (m_use_custom_hermite) {
}
if (m_dims == 1) {
for (int i = 0; i < n_points; ++i) {
Scalar val;
if (queries.cols() == 1) {
val = queries(i, 0);
} else {
val = queries(0, i);
}
result(i) = eval_symbolic_scalar(val);
}
} else {
bool is_transposed = (queries.rows() == m_dims && queries.cols() != m_dims);
for (int i = 0; i < n_points; ++i) {
for (int d = 0; d < m_dims; ++d) {
point(d) = is_transposed ? queries(d, i) : queries(i, d);
}
result(i) = eval_symbolic_point(point);
}
}
}
return result;
}
private:
m_use_custom_hermite = true;
m_valid = true;
for (int d = 0; d < m_dims; ++d) {
if (m_grid[d].size() < 4) {
throw InterpolationError(
"Interpolator: BSpline requires at least 4 grid points in dimension " +
std::to_string(d));
}
}
m_casadi_fn = casadi::interpolant("interp", "bspline", m_grid, m_values);
m_valid = true;
m_casadi_fn = casadi::interpolant("interp", "linear", m_grid, m_values);
m_valid = true;
} else {
m_casadi_fn = casadi::interpolant("interp", "linear", m_grid, m_values);
m_valid = true;
}
}
void compute_boundary_slopes() {
return;
}
m_slopes_left.resize(m_dims, 0.0);
m_slopes_right.resize(m_dims, 0.0);
for (int d = 0; d < m_dims; ++d) {
const auto &grid_d = m_grid[d];
size_t n_d = grid_d.size();
if (n_d < 2)
continue;
double x_min = grid_d.front();
double x_max = grid_d.back();
double h_left = grid_d[1] - grid_d[0];
double h_right = grid_d[n_d - 1] - grid_d[n_d - 2];
std::vector<double> query_base(m_dims);
for (int dd = 0; dd < m_dims; ++dd) {
query_base[dd] = 0.5 * (m_grid[dd].front() + m_grid[dd].back());
}
query_base[d] = x_min;
std::vector<casadi::DM> args0 = {casadi::DM(query_base)};
double f0 = static_cast<double>(m_casadi_fn(args0)[0]);
query_base[d] = x_min + h_left;
std::vector<casadi::DM> args1 = {casadi::DM(query_base)};
double f1 = static_cast<double>(m_casadi_fn(args1)[0]);
m_slopes_left[d] = (f1 - f0) / h_left;
query_base[d] = x_max - h_right;
std::vector<casadi::DM> args2 = {casadi::DM(query_base)};
double f2 = static_cast<double>(m_casadi_fn(args2)[0]);
query_base[d] = x_max;
std::vector<casadi::DM> args3 = {casadi::DM(query_base)};
double f3 = static_cast<double>(m_casadi_fn(args3)[0]);
m_slopes_right[d] = (f3 - f2) / h_right;
}
}
double apply_output_bounds(double value) const {
if (m_extrap_config.lower_bound.has_value()) {
value = std::max(value, m_extrap_config.lower_bound.value());
}
if (m_extrap_config.upper_bound.has_value()) {
value = std::min(value, m_extrap_config.upper_bound.value());
}
return value;
}
if (m_extrap_config.lower_bound.has_value()) {
result = SymbolicScalar::fmax(result, m_extrap_config.lower_bound.value());
}
if (m_extrap_config.upper_bound.has_value()) {
result = SymbolicScalar::fmin(result, m_extrap_config.upper_bound.value());
}
return result;
}
double eval_numeric_scalar(double query) const {
const double x_min = m_grid[0].front();
const double x_max = m_grid[0].back();
if (query < x_min) {
double result = m_values.front() + m_slopes_left[0] * (query - x_min);
return apply_output_bounds(result);
} else if (query > x_max) {
double result = m_values.back() + m_slopes_right[0] * (query - x_max);
return apply_output_bounds(result);
}
}
double clamped = std::max(query, x_min);
clamped = std::min(clamped, x_max);
double result;
if (m_use_custom_hermite) {
auto it = std::lower_bound(m_grid[0].begin(), m_grid[0].end(), clamped);
size_t idx = std::distance(m_grid[0].begin(), it);
if (idx > 0 && (idx == m_grid[0].size() || std::abs(clamped - m_grid[0][idx - 1]) <
std::abs(clamped - m_grid[0][idx]))) {
idx--;
}
result = m_values[idx];
} else {
std::vector<casadi::DM> args = {casadi::DM(clamped)};
std::vector<casadi::DM> res = m_casadi_fn(args);
result = static_cast<double>(res[0]);
}
return apply_output_bounds(result);
}
double extrap_correction = 0.0;
bool needs_extrapolation = false;
for (int d = 0; d < m_dims; ++d) {
double val = query(d);
double x_min = m_grid[d].front();
double x_max = m_grid[d].back();
if (val < x_min) {
extrap_correction += m_slopes_left[d] * (val - x_min);
needs_extrapolation = true;
} else if (val > x_max) {
extrap_correction += m_slopes_right[d] * (val - x_max);
needs_extrapolation = true;
}
}
}
std::vector<double> clamped(m_dims);
for (int d = 0; d < m_dims; ++d) {
double val = query(d);
val = std::max(val, m_grid[d].front());
val = std::min(val, m_grid[d].back());
clamped[d] = val;
}
double result;
if (m_use_custom_hermite) {
std::vector<casadi::DM> args = {casadi::DM(clamped)};
std::vector<casadi::DM> res = m_casadi_fn(args);
result = static_cast<double>(res[0]);
} else {
std::vector<casadi::DM> args = {casadi::DM(clamped)};
std::vector<casadi::DM> res = m_casadi_fn(args);
result = static_cast<double>(res[0]);
}
if (needs_extrapolation) {
result += extrap_correction;
}
return apply_output_bounds(result);
}
const double x_min = m_grid[0].front();
const double x_max = m_grid[0].back();
clamped = SymbolicScalar::fmin(clamped, x_max);
std::vector<casadi::MX> args = {clamped};
std::vector<casadi::MX> res = m_casadi_fn(args);
SymbolicScalar left_extrap = m_values.front() + m_slopes_left[0] * (query - x_min);
SymbolicScalar right_extrap = m_values.back() + m_slopes_right[0] * (query - x_max);
query < x_min, left_extrap,
casadi::MX::if_else(query > x_max, right_extrap, interp_result));
return apply_output_bounds_sym(result);
}
return apply_output_bounds_sym(interp_result);
}
casadi::MX clamped(m_dims, 1);
for (int d = 0; d < m_dims; ++d) {
casadi::MX val = query(d);
val = casadi::MX::fmax(val, m_grid[d].front());
val = casadi::MX::fmin(val, m_grid[d].back());
clamped(d) = val;
}
std::vector<casadi::MX> args = {clamped};
std::vector<casadi::MX> res = m_casadi_fn(args);
for (int d = 0; d < m_dims; ++d) {
double x_min = m_grid[d].front();
double x_max = m_grid[d].back();
casadi::MX val = query(d);
SymbolicScalar left_contrib = casadi::MX::if_else(val < x_min, left_corr, 0.0);
SymbolicScalar right_contrib = casadi::MX::if_else(val > x_max, right_corr, 0.0);
extrap_correction = extrap_correction + left_contrib + right_contrib;
}
return apply_output_bounds_sym(interp_result + extrap_correction);
}
return apply_output_bounds_sym(interp_result);
}
};
template <typename Scalar>
std::optional<Scalar> fill_value = std::nullopt) {
const int n_dims = static_cast<int>(points.size());
const int n_points = static_cast<int>(xi_work.rows());
if (values_flat.size() != grid_data.expected_size) {
throw InterpolationError("interpn: values_flat size mismatch. Expected " +
std::to_string(grid_data.expected_size) + ", got " +
std::to_string(values_flat.size()));
}
Interpolator interp(points, values_flat, method);
if (fill_value.has_value()) {
for (int i = 0; i < n_points; ++i) {
bool out_of_bounds = false;
if constexpr (std::is_floating_point_v<Scalar>) {
}
if (out_of_bounds) {
result(i) = fill_value.value();
} else {
result(i) = interp(point);
}
}
} else {
for (int i = 0; i < n_points; ++i) {
result(i) = interp(point);
}
}
return result;
}
template <typename Scalar>
std::optional<SymbolicScalar> fill_value = std::nullopt) {
const int n_dims = static_cast<int>(points.size());
const int n_points = static_cast<int>(xi_work.rows());
if (values_flat.size() != grid_data.expected_size) {
throw InterpolationError("interpn: values_flat size mismatch. Expected " +
std::to_string(grid_data.expected_size) + ", got " +
std::to_string(values_flat.size()));
}
for (int i = 0; i < n_points; ++i) {
SymbolicScalar interp_value = interp(std::vector<SymbolicScalar>{clamped_point, coeffs})[0];
if (fill_value.has_value()) {
if constexpr (std::is_floating_point_v<Scalar>) {
? fill_value.value()
: interp_value;
} else {
result(i) =
SymbolicScalar::if_else(out_of_bounds, fill_value.value(), interp_value);
}
} else {
result(i) = interp_value;
}
}
return result;
}
}
C++20 concepts constraining valid Janus scalar types.
Custom exception hierarchy for Janus framework.
Core type aliases for numeric and symbolic Eigen/CasADi interop.
Linear algebra operations (solve, inverse, determinant, eigendecomposition, norms).
Scalar operator()(const Scalar &query) const
Evaluate interpolant at a scalar point (1D only).
Definition Interpolate.hpp:669
InterpolationMethod method() const
Get the interpolation method.
Definition Interpolate.hpp:651
bool valid() const
Check if interpolator is valid (initialized).
Definition Interpolate.hpp:656
Interpolator()=default
Default constructor (invalid state).
int dims() const
Get number of dimensions.
Definition Interpolate.hpp:646
auto slice(const Eigen::MatrixBase< Derived > &m, Eigen::Index start, Eigen::Index end)
Definition IntegrateDiscrete.hpp:32
JanusVector< typename Derived::Scalar > flatten_fortran_order(const Eigen::MatrixBase< Derived > &values)
Flatten N-D values array in Fortran order for CasADi.
Definition Interpolate.hpp:115
std::string method_to_casadi_string(InterpolationMethod method)
Convert InterpolationMethod enum to CasADi string.
Definition Interpolate.hpp:86
double hermite_interp_1d_numeric(const std::vector< double > &x, const std::vector< double > &y, double query)
Evaluate 1D Catmull-Rom cubic Hermite spline (numeric only).
Definition Interpolate.hpp:306
double hermite_interpn_numeric(const std::vector< std::vector< double > > &grid, const std::vector< double > &values, const std::vector< double > &query)
N-dimensional Catmull-Rom interpolation via tensor product (numeric only).
Definition Interpolate.hpp:390
int multi_to_flat_index(const std::vector< int > &multi, const std::vector< int > &dims)
Compute flat index from multi-dimensional index (Fortran order).
Definition Interpolate.hpp:374
casadi::Function make_parametric_interpolant(const std::vector< std::vector< double > > &grid, InterpolationMethod method)
Definition Interpolate.hpp:242
bool point_out_of_bounds_numeric(const JanusVector< Scalar > &point, const std::vector< std::vector< double > > &grid)
Definition Interpolate.hpp:221
double catmull_rom_slope(const std::vector< double > &x, const std::vector< double > &y, size_t i)
Compute Catmull-Rom slope at a point using neighboring values.
Definition Interpolate.hpp:271
void validate_bspline_grid(const std::vector< std::vector< double > > &grid, const char *context)
Definition Interpolate.hpp:180
InterpnGridData build_interpn_grid(const std::vector< NumericVector > &points, const char *context)
Definition Interpolate.hpp:148
std::vector< int > flat_to_multi_index(int flat_idx, const std::vector< int > &dims)
Compute multi-dimensional index from flat index (Fortran order).
Definition Interpolate.hpp:361
const char * hermite_symbolic_error_message()
Definition Interpolate.hpp:101
JanusMatrix< Scalar > normalize_query_matrix(const JanusMatrix< Scalar > &xi, int n_dims, const char *context)
Definition Interpolate.hpp:192
SymbolicScalar point_out_of_bounds_symbolic(const SymbolicVector &point, const std::vector< std::vector< double > > &grid)
Definition Interpolate.hpp:232
SymbolicScalar clamp_query_point(const JanusVector< Scalar > &point, const std::vector< std::vector< double > > &grid)
Definition Interpolate.hpp:207
std::string symbolic_values_method_error_message(InterpolationMethod method)
Definition Interpolate.hpp:129
Definition Diagnostics.hpp:19
InterpolationMethod
Supported interpolation methods.
Definition Interpolate.hpp:29
@ Linear
Piecewise linear (C0 continuous) - fast, simple.
Definition Interpolate.hpp:30
@ Nearest
Nearest neighbor - fast, non-differentiable.
Definition Interpolate.hpp:33
@ Hermite
Cubic Hermite/Catmull-Rom (C1 continuous) - smooth gradients.
Definition Interpolate.hpp:31
@ BSpline
Cubic B-spline (C2 continuous) - smoothest, good for optimization.
Definition Interpolate.hpp:32
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > JanusMatrix
Dynamic-size matrix for both numeric and symbolic backends.
Definition JanusTypes.hpp:43
ExtrapolationMode
Controls behavior when query points fall outside grid bounds.
Definition Interpolate.hpp:43
@ Clamp
Clamp queries to grid bounds (default, safe).
Definition Interpolate.hpp:44
@ Linear
Linear extrapolation from boundary slope (1D only).
Definition Interpolate.hpp:45
JanusVector< Scalar > interpn(const std::vector< NumericVector > &points, const NumericVector &values_flat, const JanusMatrix< Scalar > &xi, InterpolationMethod method=InterpolationMethod::Linear, std::optional< Scalar > fill_value=std::nullopt)
N-dimensional interpolation on regular grids (backwards compatibility).
Definition Interpolate.hpp:1127
@ Hermite
Definition AutoDiff.hpp:31
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > JanusVector
Dynamic-size column vector for both numeric and symbolic backends.
Definition JanusTypes.hpp:49
casadi::MX SymbolicScalar
CasADi MX symbolic scalar.
Definition JanusTypes.hpp:70
SymbolicScalar as_mx(const SymbolicVector &v)
Get the underlying MX representation of a SymbolicVector.
Definition JanusTypes.hpp:173
std::vector< std::vector< double > > grid
Definition Interpolate.hpp:144
int expected_size
Definition Interpolate.hpp:145
JanusVector< SymbolicScalar > SymbolicVector
Eigen vector of MX elements.
Definition JanusTypes.hpp:72
JanusVector< NumericScalar > NumericVector
Eigen::VectorXd equivalent.
Definition JanusTypes.hpp:67
casadi::MX SymbolicScalar
CasADi MX symbolic scalar.
Definition JanusTypes.hpp:70