Janus 2.0.0
High-performance C++20 dual-mode numerical framework
Loading...
Searching...
No Matches
Pseudospectral

janus::Pseudospectral implements global polynomial optimal-control transcription using spectral differentiation matrices. It enforces dynamics globally via D * X = (dt / 2) * F(X, U, t), where D is a differentiation matrix on Lobatto nodes. For smooth problems this gives spectral convergence, so high accuracy is often possible with fewer nodes than local collocation. This works in symbolic mode via the janus::Opti interface. The class lives in <janus/optimization/Pseudospectral.hpp>.

Quick Start

#include <janus/janus.hpp>
opts.n_nodes = 31;
auto [x, u, tau] = ps.setup(3, 1, 0.0, 2.0, opts);
ps.set_dynamics([](const janus::SymbolicVector& x,
dxdt(0) = x(2) * janus::sin(u(0));
dxdt(1) = -x(2) * janus::cos(u(0));
dxdt(2) = 9.81 * janus::cos(u(0));
return dxdt;
});
ps.add_dynamics_constraints();
ps.set_initial_state(janus::NumericVector{{0.0, 10.0, 0.001}});
ps.set_final_state(0, 10.0);
ps.set_final_state(1, 5.0);
opti.minimize(/* objective */);
auto sol = opti.solve();
Main optimization environment class.
Definition Opti.hpp:167
void minimize(const SymbolicScalar &objective)
Set objective to minimize.
Definition Opti.hpp:555
OptiSol solve(const OptiOptions &options={})
Solve the optimization problem.
Definition Opti.hpp:603
Pseudospectral (Gauss-Lobatto) transcription.
Definition Pseudospectral.hpp:35
Umbrella header that includes the entire Janus public API.
JanusVector< SymbolicScalar > SymbolicVector
Eigen vector of MX elements.
Definition JanusTypes.hpp:72
JanusVector< NumericScalar > NumericVector
Eigen::VectorXd equivalent.
Definition JanusTypes.hpp:67
@ LGL
Legendre-Gauss-Lobatto.
Definition Pseudospectral.hpp:18
T cos(const T &x)
Computes cosine of x.
Definition Trig.hpp:46
T sin(const T &x)
Computes sine of x.
Definition Trig.hpp:21
casadi::MX SymbolicScalar
CasADi MX symbolic scalar.
Definition JanusTypes.hpp:70
Options for Pseudospectral setup.
Definition Pseudospectral.hpp:23
int n_nodes
Number of collocation nodes (including endpoints).
Definition Pseudospectral.hpp:25
PseudospectralScheme scheme
Definition Pseudospectral.hpp:24

Core API

Method Description
Pseudospectral(opti) Construct with a janus::Opti instance
setup(n_states, n_controls, t0, tf, opts) Create decision variables and time grid
set_dynamics(ode) Set the ODE function: (x, u, t) -> dxdt
add_dynamics_constraints() Apply spectral differentiation matrix constraints
add_defect_constraints() Alias for add_dynamics_constraints()
set_initial_state(x0) Set initial boundary condition
set_final_state(xf) Set final boundary condition
n_nodes() Number of collocation nodes
time_grid() Normalized time grid [0, 1]
diff_matrix() Access the differentiation matrix D
quadrature_weights() Access the quadrature weight vector
quadrature(integrand) Compute weighted integral of an integrand vector

Supported schemes:

enum class PseudospectralScheme {
LGL, // Legendre-Gauss-Lobatto
CGL // Chebyshev-Gauss-Lobatto
};

Both include endpoints, so boundary-condition setup matches DirectCollocation.

Usage Patterns

Free Final Time

auto T = opti.variable(2.0, std::nullopt, 0.1, 10.0);
auto [x, u, tau] = ps.setup(3, 1, 0.0, T, opts);
ps.set_dynamics(my_ode);
ps.add_dynamics_constraints();
ps.set_initial_state(x0);
ps.set_final_state(0, xf_x);
ps.set_final_state(1, xf_y);
opti.minimize(T);
auto sol = opti.solve();
SymbolicScalar variable(double init_guess=0.0, std::optional< double > scale=std::nullopt, std::optional< double > lower_bound=std::nullopt, std::optional< double > upper_bound=std::nullopt)
Create a scalar decision variable.
Definition Opti.hpp:200

Quadrature for Running Costs

Use quadrature() for Lagrange objectives:

janus::SymbolicVector integrand(ps.n_nodes());
for (int k = 0; k < ps.n_nodes(); ++k) {
integrand(k) = u(k, 0) * u(k, 0);
}
opti.minimize(ps.quadrature(integrand));

This computes J = (dt / 2) * sum_i w_i * L_i where w_i are the scheme quadrature weights.

Unified Comparison Example

The file examples/optimization/transcription_comparison_demo.cpp runs Direct Collocation, Multiple Shooting, Pseudospectral, and Birkhoff Pseudospectral on the same brachistochrone problem and prints a side-by-side comparison with performance metrics and a convergence study.

Advanced Usage

Accessing Internal Matrices

const auto &D = ps.diff_matrix();
const auto &w = ps.quadrature_weights();

Notes

  • time_grid() returns normalized [0, 1] values for API consistency with DirectCollocation.
  • Internally, nodes and D are built on [-1, 1] and scaled correctly in constraints/objectives.
  • Pseudospectral methods are strongest for smooth trajectories; sharp bang-bang controls can need mesh refinement or alternative transcription.

See Also