Tinyopt’s documentation

Overview

Tinyopt Builds Tinyopt-example Builds

Tinyopt

Is your optimization’s convergence rate rivaling the speed of continental drift?

Tinyopt, the header-only C++ hero, swoops in to save the day! It’s like a tiny, caffeinated mathematician living in your project, ready to efficiently tackle those small-to-large optimization beasties, including unconstrained and non-linear least squares puzzles. Perfect for when your science or engineering project is about to implode from too much math.

Tinyopt provides high-accuracy and computationally efficient optimization capabilities, supporting both dense and sparse problem structures. The library integrates a collection of iterative solvers including Gradient Descent, Gauss-Newton and Levenberg-Marquardt algorithms (more are coming).

Furthermore, to facilitate the computation of derivatives, Tinyopt seamlessly integrates the automatic differentiation capabilities which empowers users to effortlessly compute accurate gradients.

Tinyopt is open-source, licensed under the Apache 2.0 License. 🧾

Table of Contents

Installation

Usage

Benchmarks

Roadmap

Get Involved

Installation πŸ“₯

Simply clone the repo, configure and install.

git clone https://github.com/julien-michot/tinyopt
cd tinyopt && mkdir build && cd build
cmake .. -DTINYOPT_BUILD_TESTS=OFF
sudo make install

Header files will be copied to /usr/local/include.

Usage πŸ‘¨πŸ»β€πŸ’»οƒ

Tinyopt: The Easy Way πŸ˜Žοƒ

Tinyopt is inspired by the simple syntax of python so it is very developer friendly*, just call Optimize and give it something to optimize, say x and something to minimize.

Optimize performs automatic differentiation so you just have to specify the residual(s), no jacodians/derivatives to calculate because you know the pain, right? No pain, vanished, thank you Julien.

* but not compiler friendly, sorry gcc/clang but you’ll have to work double because it’s all templated.

Example: What’s the square root of 2? πŸ€“οƒ

Beause using std::sqrt is over hyped, let’s try to recover it using Tinyopt, here is how to do:

// Import the optimizer you want, say the default NLLS one
using namespace tinyopt::nlls;
// Define 'x', the parameter to optimize, initialized to '1' (yeah, who doesn't like 1?)
double x = 1;
Optimize(x, [](auto &x) { return x * x - 2.0; }); // Let's minimize Ξ΅ = x*x - 2
// 'x' is now √2, amazing.

That’s it. Is it too verbose? Well remove the comments then. Come on, it’s just two lines, I can’t do better.

Running this will give you x = √2 at the end as well as some nerdy info:

tinyopt# make run_tinyopt_test_sqrt2
πŸ’‘ #0: Ο„:0.00ms x:{1} |Ξ΄x|:5.00e-01 Ξ»:1.00e-04 Ξ΅:1.00e+00 n:1 dΞ΅:-3.403e+38 |βˆ‡|:4.000e+00
βœ… #1: Ο„:0.06ms x:{1.49995} |Ξ΄x|:8.33e-02 Ξ»:3.33e-05 Ξ΅:2.50e-01 n:1 dΞ΅:-7.502e-01 |βˆ‡|:5.618e-01
βœ… #2: Ο„:0.07ms x:{1.41667} |Ξ΄x|:2.45e-03 Ξ»:1.11e-05 Ξ΅:6.94e-03 n:1 dΞ΅:-2.429e-01 |βˆ‡|:3.871e-04
βœ… #3: Ο„:0.08ms x:{1.41422} |Ξ΄x|:2.11e-06 Ξ»:3.70e-06 Ξ΅:5.96e-06 n:1 dΞ΅:-6.938e-03 |βˆ‡|:2.842e-10
βœ… #4: Ο„:0.08ms x:{1.41421} |Ξ΄x|:4.21e-08 Ξ»:1.23e-06 Ξ΅:1.19e-07 n:1 dΞ΅:-5.841e-06 |βˆ‡|:1.137e-13
🌞 Reached minimal gradient (success)

Tinyopt: The β€œJust Works” Example (Minimal Edition)

Feeling lost? Fear not! We’ve crafted a delightful, teeny-tiny CMake project in tinyopt-example that’ll have you parsing options faster than you can say β€œcommand-line arguments.” It’s so simple, even your pet rock could probably figure it out. (Though, we haven’t tested that rigorously.)

API Documentation πŸ“šοƒ

Have a look at our API doc or delve into the full doc at ReadTheDocs.

Benchmarks: How fast is Tinyopt? πŸš€οƒ

Tinyopt is fast, really fast, one of the fastest optimization library out there!

Setup

We’re currently evaluating small dense problems (<50 dimensions) with one cost function on a Ubuntu GNU/Linux 2024.04 64b. We’re showing without Automatic Differentiation as there’s some time increase with it, but not that much. The script benchmarks/scripts/run.sh was called after making sure the CPU powermodes were all in β€˜performance’. Plotting is done using the notebook benchmarks/scripts/results.ipynb.

Benchmarks Results πŸŽ―οƒ

Check this out, we’re comparing Tinyopt against the well known Ceres solver.

Benchmarks Results

It’s Super Green πŸ’š, Korben!

Benchmarks Plot πŸ“ˆοƒ

Benchmarks Plot

Dang, this thing’s got some serious pep in its step, making other optimization libraries look like they’re enjoying a leisurely Sunday drive. We put it in the ring with the well-respected Ceres solver, and the results were… eye-opening! My beloved old Cadillac Eldorado Biarritz computer (a true classic, if a bit slow by modern standards) practically sputtered trying to process Tinyopt’s sheer velocity.

I’m still double-checking the numbers to make sure my vintage machine wasn’t just having a particularly enthusiastic day. Is it actually that quick? Well, Tinyopt’s certainly making me wonder if my computer needs a pit stop for some performance upgrades!

Note that the current benchmarks are only for somewhat small problems, I expect the timings difference will reduce as the problem size (not residuals!) increases as Ceres-Solver has nice tricks that Tinyopt hasn’t…just yet!

Roadmap πŸ—ΊοΈοƒ

Here is what is coming up. Don’t trust too much the versions as I go with the flow.

v1 (stable API + Armadillo)

  • [ ] Add l-BFGS for large sparse problems

  • [ ] Native support of Armadillo (as alternative to Eigen)

  • [ ] Refactor Numerical/Automatic differentiation to support NLLS and general optimizations

  • [ ] Update all docs

v1.x (Bindings)

  • [ ] Add C API

  • [ ] Add python binding

  • [ ] Add Rust binding Also,

  • [ ] A wider array of benchmarks

v2 (Refactoring, Speed-ups & Many Solvers)

  • [ ] Speed-up compilation (e.g. c++20 Concepts)

  • [ ] Refactor Solvers

  • [ ] Add various more solvers (CG, Adam, …) and backend (e.g. Cuda)

  • [ ] Speed-up large problems (e.g. AMD)

Ah ah, you thought I would use Jira for this list? No way.

Get Involved & Get in Touch! πŸ€οƒ

Citation πŸ“‘οƒ

If you find yourself wanting to give us a scholarly nod, feel free to use this BibTeX snippet:

@misc{michot2025,
    author = {Julien Michot},
    title = {tinyopt: A tiny optimization library},
    howpublished = "\url{https://github.com/julien-michot/tinyopt}",
    year = {2025}
}

Fancy Lending a Hand? (We’d Love That!) πŸ€©οƒ

Feel free to contribute to the project, there’s plenty of things to add, from bindings to various languages to adding more solvers, examples and code optimizations in order to make Tinyopt, truly the fastest optimization library!

Otherwise, have fun using Tinyopt ;)

Got Big Ideas (or Just Want to Chat Business)?

If your business needs a super fast πŸ”₯ Bundle Adjustment (BA), a multi-sensor SLAM or if Tinyopt is still taking its sweet time with your application and you’re finding yourself drumming your fingers impatiently, don’t despair!

Feel free to give me a shout, I can probably help!

Usage

API Documentation

Full Documentation

Dive into the glorious depths of our documentation on ReadTheDocs. It’s packed with all the juicy details, and maybe a few hidden jokes if you look hard enough.

Simple API

tinyopt is inspired by the simple syntax of python so it is very developer friendly*, just call Optimize and give it something to optimize, say x and something to minimize.

Optimize performs automatic differentiation so you just have to specify the residual(s), no jacodians/derivatives to calculate because you know the pain, right? No pain, vanished, thank you Julien.

* but not compiler friendly, sorry gcc/clang but you’ll have to work double because it’s all templated.

Example: What’s the square root of 2?

Beause using std::sqrt is over hyped, let’s try to recover it using tinyopt, here is how to do:

// Define 'x', the parameter to optimize, initialized to '1' (yeah, who doesn't like 1?)
double x = 1;
Optimize(x, [](const auto &x) { return x * x - 2.0; }); // Let's minimize Ξ΅ = x*x - 2
// 'x' is now √2, amazing.

That’s it. Is it too verbose? Well remove the comments then. Come on, it’s just two lines, I can’t do better.

Example: Fitting a circle to a set of points

In this use case, you’re given n 2D points your job is to fit a circle to them. Today is your lucky day, tinyopt is here to help! Let’s see how.

Mat2Xf obs(2, n); // fill the observations (n 2D points)

// loss is the sum of || ||p - center||Β² - radiusΒ² ||
auto loss = [&obs]<typename Derived>(const MatrixBase<Derived> &x) {
  using T = typename Derived::Scalar;
  const auto &center = x.template head<2>(); // the first two elements are the cicle position
  const auto radius2 = x.z() * x.z(); // the last one is its radius, taking the square to avoid a sqrt later on
  // Here we compute the squared distances of each point to the center
  auto residuals = (obs.cast<T>().colwise() - center).colwise().squaredNorm();
  // Here we compute the difference of of squared distances are the circle's squared radius
  return (residuals.array() - radius2).eval();
  // Make sure the returned type is a scalar or Eigen::Matrix (thus the .eval())
};


// Define 'x', the parameter to optimize, using the following parametrization: x = {center (x, y), radius}
Vec3 x(0, 0, 1);
Optimize(x, loss); // Optimize!
// 'x' should have converged to a reasonable cicle
std::cout << "Residuals: " << loss(x) << "\n"; // Let's print the final residuals

Ok so this is quite more verbose than in the first example but I’m trying to help you understand the syntax, it’s easy no?

Supported Types

With tinyopt, you can directly optimize several types of parameters x, namely

  • float or double scalar types

  • std:array of scalars or another type

  • std::vector of a scalars or another type

  • Eigen::Vector of fixed or dynamic size

  • Eigen::Matrix of fixed or dynamic size

  • Your custorm class or struct, see below

You can also use a one level nesting of types as long as the dimensions of the nested type are known at compile time, e.g. std::array<Vec2f, 2> or std::vector<Vec3>.

tinyopt will detect whether the size is known at compile time and use optimized data structs to make the optimization faster.

Residuals of the following types can also be returned

  • float or double (or the typename T, typically a Jet<S,N>)

  • Eigen::Vector or Eigen::Matrix

Advanced API

Direct Accumulation, a faster - but manual - way.

When working with more whan one residuals, tinyopt allows you to avoid storing a full vector of residuals. You can directly accumulate the residuals and manually update the full (real or approx.) Hessian and gradient matrices this way:

// Define 'x', the parameter to optimize, initialized to '1'
double x = 1;

// Define the residuals / loss function, here Ρ² = ||x*x - 2||²
auto loss = [](const auto &x, auto &grad, auto &H) {
  float res = x * x - 2; // since we want x to be sqrt(2), x*x should be 2
  float J   = 2 * x; // residual's jacobian/derivative w.r.t x
  // Manually update the Hessian H ~= Jt * J and Gradient = Jt * residuals
  if constexpr (!traits::is_nullptr_v<decltype(grad)>) { // Gradient requested?
    H(0, 0) = J * J;   // normal matrix (initialized to 0s before so only update what is needed)
    grad(0) = J * res; // gradient (half of it actually)
  }
  // Return both the norm and the number of residuals (here, we have only one)
  return std::make_pair(std::sqrt(res*res), 1);
};

// Setup optimizer options (optional)
Options options;
// Optimize!
const auto &out = Optimize(x, loss, options);
// 'x' is now std::sqrt(2.0), you can check the convergence with out.Converged()

For second order solvers, H and grad are the only things you need to update for LM to solve the normal equations and optimize x. It looks a bit rustic I know but we can’t all live in a fancy city with sleek buidlings, sometimes it’s good to go back to your (square) roots, so take your boots and start coding.

NOTE that you only need to fill the upper triangle part only since H is assumed to be symmetric.

Sparse Systems

Ok so you have quite large and sparse systems? Just tell Tinyopt about it by simply defining you loss to have a SparseMatrix<double> &H type instead of a Dense Matrix or auto.

NOTE that automatic differentation is not supported for sparse Hessian matrices but is for first order solvers. In that case, simply use this loss signature auto loss = []<typename T>(const auto &x, SparseMatrix<T> &gradient).

auto loss = [](const auto &x, auto &grad, SparseMatrix<double> &H) {
  // Define your residuals
  const VecX res = 10 * x.array() - 2; // the residuals
  // Update the full gradient matrix, using the Jacobian J of the residuals w.r.t 'x'
  MatX J = Matx::Zero(res.rows(), x.size());
  for (int i = 0; i < x.size(); ++i) J(i, i) = 10;
  grad = J.transpose() * res;
  // Update the (approx. or real) Hessian
  std::vector<Eigen::Triplet<double>> triplets;
  triplets.reserve(x.size());
  for (int i = 0; i < x.size(); ++i) triplets.emplace_back(i, i, 10 * 10);
  H.setFromTriplets(triplets.begin(), triplets.end());
  // Returns the norm + number of residuals
  return std::make_pair(res.norm(), res.size());
};

As an alternative, you can use the Optimizer<SparseMatrix<MatX>> class instead of the Optimize function.

There are many ways to fill H in Eigen, have a look at tests/sparse.cpp for some examples.

User defined parameters

Let’s say you have a fancy class/struct, say a Rectangle and you want to optimize its dimensions. For instance,

template <typename T> // Template is only needed if you need automatic differentiation
struct Rectangle {
  T area() const { return (p2 - p1).norm(); }
  T width() const { return p2.x() - p1.x(); }
  T height() const { return p2.y() - p1.y(); }
  Eigen::Vector<T, 2> center() const { return T(0.5) * (p1 + p2); }
  Eigen::Vector<T, 2> p1, p2; // top left and bottom right positions
};

Now if you wan to simply call Optimize(rectangle, loss) on your rectangle struct, you need to either add specific members and methods or use a trait specialization of params_trait for you object type:

namespace tinyopt::traits { // must be defined in tinyopt::traits

template <typename T>
struct params_trait<Rectangle<T>> {
  using Scalar = T;              // The scalar type
  static constexpr Index Dims = 4; // Compile-time parameters dimensions (use Eigen::Dynamic if unknown)
  // Execution-time parameters dimensions [OPTIONAL, if Dims is known)
  static int dims(const Rectangle<T> &) { return Dims; }

  // Convert a Rectangle to another type 'T2', e.g. T2 = Jet<T> [OPTIONAL, if no Jet]
  // Not needed if you use manual Jacobians instead of automatic differentiation
  template <typename T2> static Rectangle<T2> cast(const Rectangle<T> &rect) {
    return Rectangle<T2>(rect.p1.template cast<T2>(),
                         rect.p2.template cast<T2>());
  }

  // Define how to update the object members (parametrization and manifold)
  static void PlusEq(Rectangle<T> &rect,
                     const Eigen::Vector<Scalar, Dims> &delta) {
    // In this case delta is defined as 2 deltas for p1 and p2: [dx1, dy1, dx2, dy2]
    rect.p1 += delta.template head<2>();
    rect.p2 += delta.template tail<2>();
  }
};

} // namespace tinyopt::traits

You can also skip the trait all together if you add these members and methods directly to the class. Now you can start optimizing your custom object, e.g.

// Let's say I want the rectangle area to be 10*20, the width = 2 * height and
// the center at (1, 2).
auto loss = [&]<typename T>(const Rectangle<T> &rect) {
  using std::max;
  Eigen::Vector<T, 4> residuals;
  residuals[0] = rect.area() - 10.0f * 20.0f;
  residuals[1] = 100.0f * (rect.width() / max(rect.height(), T(1e-8f)) -
                            2.0f); // the 1e-8 is to prevent division by 0
  residuals.template tail<2>() = rect.center() - Eigen::Vector<T, 2>(1, 2);
  return residuals;
};

Rectangle rectangle;
Optimize(rectangle, loss);
// That's it, rectangle is now fitted to your loss

Numerical Differentiation

Not all cost functions are the same. By default, tinyopt will try to use automatic differentiation when the function has only one parameter x but if your function does not allow it, you can use a numerical differentiation one. Here is an exmaple

auto original_loss = [&](const auto &x) -> Vec3 { return 2 * (x - y_prior); };
auto new_loss = CreateNumDiffFunc1(x, original_loss);
// you can now pass this 'new_loss' to an optimizer, e.g. Optimize(x, new_loss);

NOTE CreateNumDiffFunc1 is when using first order optimizers which use the gradient only and CreateNumDiffFunc2 for second or pseudo-second order methods, which use both gradient and Hessian.

Losses and Norms

You can play with different losses, robust norms and M-estimators, have a look at losses.h.

Here is an example of a loss that uses a Mahalanobis distance with a covariance C.

auto loss = [&]<typename T>(const Eigen::Vector<T, 2> &x) {
  const Matrix<T, 2, 2> C_ = C.template cast<T>();
  return MahaSquaredNorm(x - y, C_);  // return res.T * C.inv() * res
};

API

struct DefaultOptions : public tinyopt::Options2

Default Options struct in case _Options is a nullptr_t.

Public Functions

inline DefaultOptions(const Options2 options = {})

Public Members

SolverType::Options solver
struct ErrorScaling

Public Members

bool use_squared_norm = true

Use squared norm instead of norm (faster)

bool downscale_by_2 = false

Rescale the cost by 0.5

bool normalize = false

Normalize the final error by the number of residuals (after use_squared_norm)

struct Export

Subclassed by tinyopt::Options2::Export

Public Members

bool acc_dx = false

Saves the accumulated delta dx as part of the output results.

struct Export : public tinyopt::Options1::Export

Public Members

bool H = true

Saves the last Hessian H as part of the output results.

template<typename T>
struct is_bool : public std::is_same<std::decay_t<T>, bool>
template<typename T>
struct is_jet_type

Public Static Attributes

static constexpr bool value = false
template<typename T, int N>
struct is_jet_type<diff::Jet<T, N>>

Public Static Attributes

static constexpr bool value = true
template<typename T>
struct is_matrix_or_array : public std::disjunction<std::is_base_of<MatrixBase<T>, T>, std::is_base_of<ArrayBase<T>, T>>
template<typename T>
struct is_nullptr_t : public std::is_same<std::decay_t<T>, std::nullptr_t>
template<typename T>
struct is_pair : public std::false_type
template<typename T, typename U>
struct is_pair<std::pair<T, U>> : public std::true_type
template<typename T>
struct is_scalar : public std::is_scalar<std::decay_t<T>>
template<typename T, int N>
struct is_scalar<diff::Jet<T, N>>

Public Static Attributes

static constexpr bool value = true
template<typename T>
struct is_sparse_matrix : public std::false_type
template<typename T>
struct is_sparse_matrix<SparseMatrix<T>> : public std::true_type
template<typename T, typename = void>
struct is_streamable : public std::false_type
template<typename T>
struct is_streamable<T, typename std::enable_if<std::is_convertible<decltype(std::declval<std::ostream&>() << std::declval<T>()), std::ostream&>::value>::type> : public std::true_type
template<typename T>
struct jet_details

Public Types

using Scalar = double

Public Static Attributes

static constexpr Index Dims = 0
template<typename T, int N>
struct jet_details<diff::Jet<T, N>>

Public Types

using Scalar = T

Public Static Attributes

static constexpr Index Dims = N
template<typename SolverType, typename _Options = std::nullptr_t>
class Optimizer

Public Types

using Scalar = typename SolverType::Scalar
using OutputType = Output<typename SolverType::H_t>
using Options = std::conditional_t<std::is_null_pointer_v<_Options>, DefaultOptions, _Options>

Public Functions

inline Optimizer(const Options &_options = {})
template<int FO = SolverType::FirstOrder, std::enable_if_t<!FO, int> = 0>
inline void InitWith(const auto &g, const auto &h)

Initialize solver with specific gradient and hessian.

template<int FO = SolverType::FirstOrder, std::enable_if_t<FO, int> = 0>
inline void InitWith(const auto &g)

Initialize solver with specific gradient.

inline void reset()

Reset the optimization and solver.

template<typename X_t>
inline std::variant<StopReason, bool> ResizeIfNeeded(X_t &x, OutputType &out)
template<typename X_t, typename AccFunc>
inline OutputType operator()(X_t &x, const AccFunc &acc, int max_iters = -1)

Main optimization function.

template<typename X_t, typename AccFunc>
inline OutputType Optimize(X_t &x, const AccFunc &acc, int max_iters = -1)

Main optimization loop.

Public Static Attributes

static constexpr Index Dims = SolverType::Dims

Protected Functions

template<typename X_t, typename AccFunc>
inline std::pair<bool, std::optional<Vector<Scalar, Dims>>> Step(X_t &x, const AccFunc &acc, OutputType &out)

Run one optimization iteration, return the estimated next step (solve + decreased error)

Protected Attributes

const Options options_

Optimization options.

SolverType solver_

Linear solver.

struct Options : public tinyopt::Options1
#include <gd.h>

Gradient Descent Optimization Options.

Public Functions

inline Options(const Options1 options = {})

Public Members

gd::SolverOptions solver
struct Options : public tinyopt::Options2
#include <gn.h>

Gauss-Newton Optimization Options.

Public Functions

inline Options(const Options2 options = {})

Public Members

gn::SolverOptions solver
struct Options : public tinyopt::Options2
#include <lm.h>

Levenberg-Marquardt Optimization Options.

Public Functions

inline Options(const Options2 options = {})

Public Members

lm::SolverOptions solver
struct Options1

Subclassed by tinyopt::Options2, tinyopt::gd::Options

Optimization options

bool check_final_err = false

Recompute the current error with latest state to eventually roll back. Only performed at the very last iteration as a safety measure (to prevent unlucky divergence at the very end…).

bool use_step_quality_approx = false

Use relative error decrease as step quality, other 0.0 will be used.

Stop criteria

uint16_t max_iters = 50

Maximum number of outter iterations.

float min_error = 1e-6

Minimum error.

float min_rerr_dec = 1e-6

Minimum relative error (Ξ΅_rel = (Ξ΅_prev-Ξ΅_new)/Ξ΅_prev)

float min_step_norm2 = 1e-16

Minimum step (dx) squared norm.

float min_grad_norm2 = 1e-18f

Minimum gradient squared norm.

uint8_t max_total_failures = 0

Overall max failures to decrease error.

uint8_t max_consec_failures = 5

Maximum consecutive failures to decrease error.

double max_duration_ms = 0

Maximum optimization duration in milliseconds (ms)

std::function<bool(double, double, double)> stop_callback

User defined callback. It will be called with the current error, step size and the gradient norm, i.e. stop = stop_callback(Ξ΅, |Ξ΄x|Β², βˆ‡). The user returns true to stop the optimization iterations early.

std::function<bool(float, const VecXf&, const VecXf&)> stop_callback2

User defined callback. It will be called with the current error, step vector and the gradient, i.e. stop = stop_callback(Ξ΅, Ξ΄x, βˆ‡). The user returns true to stop the optimization iterations early.

Logging Options

struct tinyopt::Options1 log

Export Options

struct tinyopt::Options1::Export save

Public Members

bool enable = true

Whether to enable the logging.

std::string e = "Ρ²"

Symbol used when logging the error, e.g Ρ, Ρ² or √Ρ etc.

bool print_emoji = true

Whether to show the emoji or not.

bool print_x = true

Log the value of β€˜x’.

bool print_t = true

Log the duration (in ms)

bool print_J_jet = false

Log the value of β€˜J’ from the Jet.

bool print_failure = true

Log the value of β€˜H’ and β€˜grad’ from the Jet.

bool print_max_stdev = false

Log the maximum of all standard deviations (sqrt((co-)variance)) (need to invert H)

struct Options1

Subclassed by tinyopt::gd::SolverOptions, tinyopt::solvers::Options2

Public Functions

inline Options1()

Public Members

float grad_clipping = 0

Gradient clipping to range [-v, +v], disabled if 0.

bool enable = true
struct tinyopt::solvers::Options1 log
struct Options2 : public tinyopt::Options1

Subclassed by tinyopt::nlls::gn::Options, tinyopt::nlls::lm::Options

Export Options

tinyopt::Options2::Export save

Public Functions

inline Options2(const Options1 &options = {})
struct Options2 : public tinyopt::solvers::Options1

Subclassed by tinyopt::nlls::lm::SolverOptions

Error scaling options (mostly for NLLS solvers really)

Todo:

move it to a tinyopt::solvers::nlls::Options?

struct tinyopt::solvers::Options2::ErrorScaling err

H Properties

bool use_ldlt = true

If not, will use H.inverse() without any checks on invertibility except for Dims==1

bool H_is_full = true

Specify if H is only Upper triangularly or fully filled.

Checks

float check_min_H_diag = 0

Check the the hessian’s diagonal are not all below the threshold. Use 0 to disable the check.

Public Functions

inline Options2(const Options1 &options = {})
template<typename _H_t = std::nullptr_t>
struct Output

Statistics

uint16_t num_residuals = 0

Final number of residuals.

uint16_t num_iters = 0

Final number of iterations.

uint8_t num_failures = 0

Final number of failures to decrease the error.

uint8_t num_consec_failures = 0

Final number of the last consecutive failures to decrease the error.

TimePoint start_time = TimePoint::min()

Starting time of the optimization or std::chrono::system_clock::time_point::min() if not started

float duration_ms = 0

Cumulated optimization duration.

H_t final_H

Final H, excluding any damping (only saved if options.save.H = true)

Per iteration results

std::vector<Scalar> errs

Mean squared accumulated errors of all iterations.

std::vector<Scalar> deltas2

Squared Step sizes of all iterations.

std::vector<bool> successes

Step acceptation status for all iterations.

Public Types

using H_t = _H_t
using Scalar = std::conditional_t<traits::is_nullptr_v<H_t>, double, typename traits::params_trait<H_t>::Scalar>
using Dx_t = std::conditional_t<traits::is_nullptr_v<H_t>, Vector<Scalar, Dynamic>, Vector<Scalar, SQRT(traits::params_trait<H_t>::Dims)>>

Public Functions

inline bool Succeeded() const

Returns true if the stop reason is not a failure to solve or NaNs or missing residuals.

inline bool Converged() const

Returns true if the optimization reached the specified minimal error, delta or gradient.

template<typename H = H_t, std::enable_if_t<!std::is_null_pointer_v<H>, int> = 0>
inline std::optional<H_t> Covariance(bool rescaled = false) const

Computes an approximation of the covariance matrix.

This method calculates the covariance matrix, which is an approximation derived from the inverse of the Hessian matrix (H).

It is only meaningful at convergence, right Johan?

The covariance matrix is calculated as:

  • If rescaled is false: (H)^-1

  • If rescaled is true and num_residuals > H.cols(): (H)^-1 * (Ρ² / (#Ξ΅ - dims))

Where:

  • H is the approximated Hessian matrix.

  • Ξ΅ (final_err) is the sum of residuals.

  • #Ξ΅ (num_residuals) is the number of residuals.

  • dims (H.cols()) is the number of parameters.

The rescaling is applied only when the number of residuals exceeds the number of parameters, indicating an overdetermined system. This rescaling provides a more statistically sound estimate of the covariance in such scenarios, especially when noise modeling was not explicitly included in the observations.

Note

The function relies on the InvCov(H) method to compute the inverse of H. If InvCov returns std::nullopt, indicating that H is not invertible, this method also returns std::nullopt.

Parameters:

rescaled – (optional) If true, the covariance matrix is rescaled by Ρ² / (#Ξ΅ - dims), where Ρ² is the sum of squared residuals (final_errΒ²), #Ξ΅ is the number of residuals (num_residuals), and dims is the number of parameters (H.cols()). This rescaling is useful when observations lack explicit noise modeling and provides a more accurate estimate of the covariance. Defaults to false.

Returns:

std::optional<H_t> Returns the computed covariance matrix if H is invertible, wrapped in a std::optional. Returns std::nullopt if H is not invertible, indicating that the covariance cannot be estimated.

Template Parameters:

H_t – The type of the covariance matrix.

Public Members

Scalar final_err = std::numeric_limits<Scalar>::max()

Last valid error.

Scalar final_rerr_dec = std::numeric_limits<Scalar>::max()
StopReason stop_reason = StopReason::kNone

Stop reason.

template<typename T, typename = void>
struct params_trait

Public Types

using Scalar = typename T::Scalar

Public Static Functions

static inline Index dims(const T &v)
template<typename T2>
static inline auto cast(const T &v)
static inline void PlusEq(T &v, const auto &delta)

Public Static Attributes

static constexpr Index Dims = T::Dims
template<typename _Scalar, int N>
struct params_trait<diff::Jet<_Scalar, N>>

Public Types

using T = diff::Jet<_Scalar, N>
using Scalar = _Scalar

Public Static Functions

static inline auto dims(const T &v)
template<typename T2>
static inline auto cast(const T &v)
static inline void PlusEq(T &v, const auto &delta)

Public Static Attributes

static constexpr Index Dims = params_trait<Scalar>::Dims == Dynamic ? Dynamic : 1
template<typename _Scalar, std::size_t N>
struct params_trait<std::array<_Scalar, N>>

Public Types

using T = typename std::array<_Scalar, N>
using Scalar = _Scalar

Public Static Functions

static inline Index dims(const T &v)
template<typename T2>
static inline auto cast(const T &v)
static inline void PlusEq(T &v, const auto &delta)

Public Static Attributes

static constexpr Index Dims = params_trait<Scalar>::Dims == Dynamic ? Dynamic : N * params_trait<Scalar>::Dims
template<typename T1, typename T2>
struct params_trait<std::pair<T1, T2>>

Public Types

using T = std::pair<T1, T2>
using Scalar = typename params_trait<T1>::Scalar

Public Static Functions

static inline Index dims(const T &v)
template<typename T3>
static inline auto cast(const T &v)
static inline void PlusEq(T &v, const auto &delta)

Public Static Attributes

static constexpr Index Dims = (params_trait<T1>::Dims == Dynamic || params_trait<T2>::Dims == Dynamic) ? Dynamic : params_trait<T1>::Dims + params_trait<T2>::Dims
template<typename _Scalar>
struct params_trait<std::vector<_Scalar>>

Public Types

using T = typename std::vector<_Scalar>
using Scalar = _Scalar

Public Static Functions

static inline Index dims(const T &v)
template<typename T2>
static inline auto cast(const T &v)
static inline void PlusEq(T &v, const auto &delta)

Public Static Attributes

static constexpr Index Dims = Dynamic
template<typename T>
struct params_trait<T, std::enable_if_t<is_matrix_or_array_v<T>>>

Public Types

using Scalar = typename T::Scalar

Public Static Functions

static inline Index dims(const T &m)
template<typename T2>
static inline auto cast(const T &v)
static inline void PlusEq(T &v, const auto &delta)

Public Static Attributes

static constexpr int ColsAtCompileTime = T::ColsAtCompileTime
static constexpr Index Dims = (T::RowsAtCompileTime == Dynamic || T::ColsAtCompileTime == Dynamic) ? Dynamic : T::RowsAtCompileTime * T::ColsAtCompileTime
template<typename T>
struct params_trait<T, std::enable_if_t<is_sparse_matrix_v<T>>>

Public Types

using Scalar = typename T::Scalar

Public Static Functions

static inline Index dims(const T &m)
template<typename T2>
static inline auto cast(const T &v)
static inline void PlusEq(T &v, const auto &delta)

Public Static Attributes

static constexpr Index Dims = Dynamic
template<typename T>
struct params_trait<T, std::enable_if_t<std::is_scalar_v<T>>>

Public Types

using Scalar = T

Public Static Functions

static inline constexpr Index dims(const T&)
template<typename T2>
static inline T2 cast(const T &v)
static inline void PlusEq(T &v, const auto &delta)
static inline void PlusEq(T &v, const Scalar &delta)

Public Static Attributes

static constexpr Index Dims = 1
template<typename _Scalar, int _Dims>
class SolverBase

Subclassed by tinyopt::solvers::SolverGD< Gradient_t >, tinyopt::solvers::SolverGN< Hessian_t >, tinyopt::solvers::SolverLM< Hessian_t >

Public Types

using Scalar = _Scalar

Public Functions

inline SolverBase(const solvers::Options1 &options = {})
template<typename Grad_t>
inline bool Clamp(Grad_t &g, Scalar minmax) const

Clamp the gradient β€˜g’ to within [-minmax, minmax], if minmax is not 0. Returns true if β€˜g’ was clamped.

virtual std::optional<Vector<Scalar, Dims>> Solve() const = 0

Solve the linear system dx = -H^-1 * grad, returns nullopt on failure.

inline virtual void GoodStep(Scalar = 0.0f)
inline virtual void BadStep(Scalar = 0.0f)
inline virtual void FailedStep()
inline virtual void Rebuild(bool)
inline virtual std::string stateAsString() const
inline Scalar Error() const
inline int NumResiduals() const

Public Static Attributes

static constexpr Index Dims = _Dims

Protected Attributes

const solvers::Options1 options_
Scalar err_ = std::numeric_limits<Scalar>::max()
int nerr_ = 0
template<typename Gradient_t = VecX>
class SolverGD : public tinyopt::solvers::SolverBase<Gradient_t::Scalar, traits::params_trait<VecX>::Dims>

Public Types

using Base = SolverBase<typename Gradient_t::Scalar, traits::params_trait<Gradient_t>::Dims>
using Scalar = typename Gradient_t::Scalar
using Grad_t = Gradient_t
using H_t = std::nullptr_t
using Options = gd::SolverOptions

Public Functions

inline explicit SolverGD(const Options &options = {})
inline void InitWith(const Grad_t &g)

Initialize solver with specific gradient and hessian.

inline void reset()

Reset the solver state and clear gradient & hessian.

template<int D = Dims, std::enable_if_t<D == Dynamic, int> = 0>
inline bool resize(int dims)

Resize H and grad if needed, return true if they were resized.

template<int D = Dims, std::enable_if_t<D != Dynamic, int> = 0>
inline bool resize(int dims = Dims)

Resize H and grad if needed, return true if they were resized.

inline void clear()

Set gradient and hessian to 0s.

template<typename X_t>
inline bool ResizeIfNeeded(const X_t &x)

Check whether we need to resize the system (gradient), return true if it did.

template<typename AccFunc>
inline auto GetAccFunc(const AccFunc &acc) const

Ugly function to make sure the returned type is always a pair<scalar, scalar>…

template<typename X_t, typename AccFunc>
inline Scalar Evaluate(const X_t &x, const AccFunc &acc, bool save)

Accumulate residuals and return the final error.

template<typename X_t, typename AccFunc>
inline bool Accumulate(const X_t &x, const AccFunc &acc)

Accumulate residuals and update the gradient, returns true on success.

template<typename X_t, typename AccFunc>
inline bool Build(const X_t &x, const AccFunc &acc, bool resize_and_clear = true)

Build the gradient and hessian by accumulating residuals and their jacobians Returns true on success

inline virtual std::optional<Vector<Scalar, Dims>> Solve() const override

Solve the linear system dx = -lr * grad, returns nullopt on failure.

inline const Grad_t &Gradient() const
inline Grad_t &Gradient()
inline Scalar GradientNorm() const
inline Scalar GradientSquaredNorm() const

Public Static Attributes

static constexpr bool IsNLLS = false
static constexpr bool FirstOrder = true
static constexpr Index Dims = traits::params_trait<Gradient_t>::Dims

Protected Attributes

const Options options_
Grad_t grad_
template<typename Hessian_t = MatX>
class SolverGN : public tinyopt::solvers::SolverBase<Hessian_t::Scalar, SQRT(traits::params_trait<MatX>::Dims)>

Public Types

using Base = SolverBase<typename Hessian_t::Scalar, SQRT(traits::params_trait<Hessian_t>::Dims)>
using Scalar = typename Hessian_t::Scalar
using H_t = Hessian_t
using Grad_t = Vector<Scalar, Dims>
using Options = nlls::gn::SolverOptions

Public Functions

inline explicit SolverGN(const Options &options = {})
inline void InitWith(const Grad_t &g, const H_t &h)

Initialize solver with specific gradient and hessian.

inline void reset()

Reset the solver state and clear gradient & hessian.

template<int D = Dims, std::enable_if_t<D == Dynamic, int> = 0>
inline bool resize(int dims)

Resize H and grad if needed, return true if they were resized.

template<int D = Dims, std::enable_if_t<D != Dynamic, int> = 0>
inline bool resize(int dims = Dims)

Resize H and grad if needed, return true if they were resized.

inline void clear()

Set gradient and hessian to 0s.

template<typename X_t>
inline bool ResizeIfNeeded(const X_t &x)

Check whether we need to resize the system (gradient), return true if it did.

template<typename ResidualsFunc>
inline auto GetAccFunc(const ResidualsFunc &res_func) const

Accumulate residuals and update the gradient, returns true on success.

template<typename X_t, typename ResidualsFunc>
inline Scalar Evaluate(const X_t &x, const ResidualsFunc &res_func, bool save)

Accumulate residuals and return the final error.

template<typename X_t, typename ResidualsFunc>
inline bool Accumulate(const X_t &x, const ResidualsFunc &res_func)

Accumulate residuals and update the gradient, returns true on success.

template<typename X_t, typename ResidualsFunc>
inline bool Build(const X_t &x, const ResidualsFunc &res_func, bool resize_and_clear = true)

Build the gradient and hessian by accumulating residuals and their jacobians Returns true on success

inline virtual std::optional<Vector<Scalar, Dims>> Solve() const override

Solve the linear system dx = -H^-1 * grad, returns nullopt on failure.

inline const H_t &Hessian() const

Latest Hessian approximation (JtJ), un-damped.

inline Scalar MaxStdDev(bool use_damped = true) const

Return the square root of the maximum (co)variance of the H.inv() H being the damped Hessian H_ if use_damped == true (faster) or un-damped Hessian() (accurate)

inline const H_t &H() const

Latest, eventually damped Hessian approximation (JtJ)

inline H_t &H()
inline const Grad_t &Gradient() const
inline Grad_t &Gradient()
inline Scalar GradientNorm() const
inline Scalar GradientSquaredNorm() const

Public Static Attributes

static constexpr bool IsNLLS = true
static constexpr bool FirstOrder = false
static constexpr Index Dims = SQRT(traits::params_trait<Hessian_t>::Dims)

Protected Attributes

const Options options_
H_t H_
Grad_t grad_
template<typename Hessian_t = MatX>
class SolverLM : public tinyopt::solvers::SolverBase<Hessian_t::Scalar, SQRT(traits::params_trait<MatX>::Dims)>

Public Types

using Base = SolverBase<typename Hessian_t::Scalar, SQRT(traits::params_trait<Hessian_t>::Dims)>
using Scalar = typename Hessian_t::Scalar
using H_t = Hessian_t
using Grad_t = Vector<Scalar, Dims>
using Options = nlls::lm::SolverOptions

Public Functions

inline explicit SolverLM(const Options &options = {})
inline void InitWith(const Grad_t &g, const H_t &h)

Initialize solver with specific gradient and hessian.

inline void reset()

Reset the solver state and clear gradient & hessian.

template<int D = Dims, std::enable_if_t<D == Dynamic, int> = 0>
inline bool resize(int dims)

Resize H and grad if needed, return true if they were resized.

template<int D = Dims, std::enable_if_t<D != Dynamic, int> = 0>
inline bool resize(int dims = Dims)

Resize H and grad if needed, return true if they were resized.

inline void clear()

Set gradient and hessian to 0s.

template<typename X_t>
inline bool ResizeIfNeeded(const X_t &x)

Check whether we need to resize the system (gradient), return true if it did.

template<typename ResidualsFunc>
inline auto GetAccFunc(const ResidualsFunc &res_func) const

Accumulate residuals and update the gradient, returns true on success.

template<typename X_t, typename ResidualsFunc>
inline Scalar Evaluate(const X_t &x, const ResidualsFunc &res_func, bool save)

Accumulate residuals and return the final error.

template<typename X_t, typename ResidualsFunc>
inline bool Accumulate(const X_t &x, const ResidualsFunc &res_func)

Accumulate residuals and update the gradient, returns true on success.

template<typename X_t, typename ResidualsFunc>
inline bool Build(const X_t &x, const ResidualsFunc &res_func, bool resize_and_clear = true)

Build the gradient and hessian by accumulating residuals and their jacobians Returns true on success

inline virtual std::optional<Vector<Scalar, Dims>> Solve() const override

Solve the linear system dx = -H^-1 * grad, returns nullopt on failure.

inline void GoodStep(Scalar quality) override
inline void BadStep(Scalar = 0.0f) override
inline virtual void FailedStep() override
inline virtual void Rebuild(bool b) override
inline virtual std::string stateAsString() const override
inline auto Hessian() const

Latest Hessian approximation (JtJ), un-damped.

inline Scalar MaxStdDev(bool use_damped = true) const

Return the square root of the maximum (co)variance of the H.inv() H being the damped Hessian H_ if use_damped == true (faster) or un-damped Hessian() (accurate)

inline const H_t &H() const

Latest, eventually damped Hessian approximation (JtJ)

inline H_t &H()
inline const Grad_t &Gradient() const
inline Grad_t &Gradient()
inline Scalar GradientNorm() const
inline Scalar GradientSquaredNorm() const

Public Static Attributes

static constexpr bool IsNLLS = true
static constexpr bool FirstOrder = false
static constexpr Index Dims = SQRT(traits::params_trait<Hessian_t>::Dims)

Protected Attributes

const Options options_
H_t H_

Hessian Approximate (Jt*J)

Grad_t grad_

Gradient (Jt*residuals)

Scalar lambda_ = 1e-4f
Scalar prev_lambda_ = 0
Scalar bad_factor_ = 2.0f
int steps_count_ = 0
bool rebuild_linear_system_ = true
struct SolverOptions : public tinyopt::solvers::Options1

Public Functions

inline SolverOptions(const solvers::Options1 &options = {})

Public Members

float lr = 1

Learning rate. The step dx will be -lr * gradient.

struct SolverOptions : public tinyopt::solvers::Options2

Damping options

float damping_init = 1e-4f

Min and max damping values (only used when damping_init != 0)

Initial damping factor. If 0, the damping is disable (it will behave like Gauss-Newton)

std::array<float, 2> damping_range = {{1e-9f, 1e9f}}
float good_factor = 1.0f / 3.0f

Scale to apply to the damping for good steps.

float bad_factor = 2.0f

Scale to apply to the damping for bad steps.

Public Functions

inline SolverOptions(const solvers::Options2 options = {})
namespace std

STL namespace.

namespace tinyopt

Typedefs

using Index = Eigen::Index
template<typename Scalar, int Rows = Dynamic, int Cols = Dynamic, int Options = 0, int MaxRows = Rows, int MaxCols = Cols>
using Matrix = std::conditional_t<Rows != 1 || (Rows == 1 && Cols == 1), Eigen::Matrix<Scalar, Rows, Cols, Options, MaxRows, MaxCols>, Eigen::Matrix<Scalar, Rows, Cols, Eigen::RowMajor, MaxRows, MaxCols>>
template<typename Scalar, int Rows = Dynamic>
using Vector = Matrix<Scalar, Rows, 1>
template<typename T, int Cols>
using RowVector = Matrix<T, 1, Cols, Eigen::RowMajor>
template<typename Scalar = double, int Options = 0, typename StorageIndex = int>
using SparseMatrix = Eigen::SparseMatrix<Scalar, Options, StorageIndex>
using SparseMatrixf = SparseMatrix<float>
template<typename T>
using MatrixBase = Eigen::MatrixBase<T>
template<typename T>
using DenseBase = Eigen::DenseBase<T>
template<typename T>
using ArrayBase = Eigen::ArrayBase<T>
using MatX = Matrix<double, Dynamic, Dynamic>
using MatXf = Matrix<float, Dynamic, Dynamic>
using Mat1 = Matrix<double, 1, 1>
using Mat2 = Matrix<double, 2, 2>
using Mat3 = Matrix<double, 3, 3>
using Mat4 = Matrix<double, 4, 4>
using Mat5 = Matrix<double, 5, 5>
using Mat6 = Matrix<double, 6, 6>
using Mat1f = Matrix<float, 1, 1>
using Mat2f = Matrix<float, 2, 2>
using Mat3f = Matrix<float, 3, 3>
using Mat4f = Matrix<float, 4, 4>
using Mat5f = Matrix<float, 5, 5>
using Mat6f = Matrix<float, 6, 6>
using VecX = Vector<double, Dynamic>
using Vec1 = Vector<double, 1>
using Vec2 = Vector<double, 2>
using Vec3 = Vector<double, 3>
using Vec4 = Vector<double, 4>
using Vec5 = Vector<double, 5>
using Vec6 = Vector<double, 6>
using Vec12 = Vector<double, 12>
using Vec21 = Vector<double, 21>
using VecXf = Vector<float, Dynamic>
using Vec1f = Vector<float, 1>
using Vec2f = Vector<float, 2>
using Vec3f = Vector<float, 3>
using Vec4f = Vector<float, 4>
using Vec5f = Vector<float, 5>
using Vec6f = Vector<float, 6>
using Vec12f = Vector<float, 12>
using Vec21f = Vector<float, 21>
using Mat2X = Matrix<double, 2, Dynamic>
using Mat3X = Matrix<double, 3, Dynamic>
using Mat2Xf = Matrix<float, 2, Dynamic>
using Mat3Xf = Matrix<float, 3, Dynamic>
using Mat23 = Matrix<double, 2, 3>
using Mat32 = Matrix<double, 3, 2>
using Mat23f = Matrix<float, 2, 3>
using Mat32f = Matrix<float, 3, 2>
using Clock = std::chrono::high_resolution_clock
using TimePoint = std::chrono::high_resolution_clock::time_point

Enums

enum StopReason

The reason why the optimization terminated.

Values:

enumerator kOutOfMemory

Out of memory when allocating the system (Hessian(s)

enumerator kSolverFailed

Failed to solve the normal equations (H is not definite positive)

enumerator kSystemHasNaNOrInf

Residuals or Jacobians have NaNs or Infinity.

enumerator kSkipped

The system has no residuals or nothing to optimize or H is all 0s.

enumerator kNone

No stop, used by Step() or when no iterations done (success)

enumerator kMinError

Minimal error reached (success)

enumerator kMinRelError

Minimal relative error decrease reached (success)

enumerator kMinDeltaNorm

Minimal delta norm reached (success)

enumerator kMinGradNorm

Minimal gradient reached (success)

enumerator kMaxIters

Maximum number of iterations reached (success)

enumerator kMaxNoDecr

Failed to decrease error too many times (success)

enumerator kMaxConsecNoDecr

Failed to decrease error consecutively too many times (success)

enumerator kTimedOut

Total allocated time reached (success)

enumerator kUserStopped

User stopped the process (success)

Functions

template<bool IsNLLS, typename X_t, typename ResidualsFunc, typename OptimizeFunc, typename OptionsType>
inline auto OptimizeWithAutoDiff(X_t &X, const ResidualsFunc &residuals, const OptimizeFunc &optimize, const OptionsType &options)
std::string format2(const std::string &format_string, const std::vector<std::string> &args)

Dummy function that replaces {*} with the arg. Does not support formatting as such!

template<typename ...Args>
std::string format(const std::string &format_string, Args&&... args)

Dummy function that replaces {*} with the arg. Does not support formatting as such!

template<typename Derived>
std::optional<Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime>> InvCov(const MatrixBase<Derived> &m)

Computes the inverse of a symmetric, semi-definite matrix.

This function calculates the inverse of a symmetric matrix, commonly used for covariance or information matrices. It leverages the LDLT decomposition for efficiency and numerical stability.

Eigen::MatrixXd covarianceMatrix(3, 3);
covarianceMatrix << 1.0, 0.5, 0.2,
0.5, 2.0, 0.8,
0.2, 0.8, 3.0;

auto inverseCovariance = InvCov(covarianceMatrix);
if (inverseCovariance)
  std::cout << "Inverse Covariance Matrix:\n" << inverseCovariance.value() << std::endl;

Note

The input matrix is assumed to be symmetric. If only the upper triangular part is filled, the function implicitly uses the symmetry to construct the full matrix for the LDLT decomposition.

Parameters:

m – [in] The symmetric input matrix. It can be filled either fully or only in the upper triangular part.

Template Parameters:

Scalar – The scalar type of the matrix elements.

Returns:

The inverse of the input matrix or nullopt, with the same dimensions and scalar type as the input.

template<typename Scalar>
std::optional<SparseMatrix<Scalar>> InvCov(const SparseMatrix<Scalar> &m)

Computes the inverse of a symmetric, semi-definite sparse matrix.

This function calculates the inverse of a symmetric matrix, commonly used for covariance or information matrices. It leverages the LDLT decomposition for efficiency and numerical stability.

Note

The input matrix is assumed to be symmetric. If only the upper triangular part is filled, the function implicitly uses the symmetry to construct the full matrix for the LDLT decomposition.

Parameters:

m – [in] The symmetric input matrix. It can be filled either fully or only in the upper triangular part.

Template Parameters:

Scalar – The scalar type of the matrix elements.

Returns:

The inverse of the input matrix or nullopt, with the same dimensions and scalar type as the input.

template<typename Derived, typename Derived2>
std::optional<Vector<typename Derived::Scalar, Derived::RowsAtCompileTime>> SolveLDLT(const MatrixBase<Derived> &A, const MatrixBase<Derived2> &b)

Solves the linear system A * X = B for X, where A is a symmetric positive-definite matrix, return optional otherwise.

This function utilizes the LDLT decomposition (Cholesky decomposition) to efficiently solve the linear system. It assumes that A is a symmetric positive-definite matrix, which is a requirement for the LDLT decomposition.

Eigen::MatrixXd A = Eigen::MatrixXd::Random(3, 3);
A = A * A.transpose(); // Ensure A is symmetric positive-definite (or semi-definite)
Eigen::VectorXd b = Eigen::VectorXd::Random(3);

auto x_opt = Solve(A, b);

if (x_opt) {
 Eigen::VectorXd x = x_opt.value();
 std::cout << "Solution X: " << x.transpose() << std::endl;
 std::cout << "A * X: " << (A * x).transpose() << std::endl;
 std::cout << "Expected -B: " << b.transpose() << std::endl;
} else {
 std::cout << "Matrix A is not positive-definite. System cannot be solved." << std::endl;
}

Note

The function returns chol.solve(b), which implies that the solution returned is actually the solution of A * X = B.

Template Parameters:
  • Derived – The type of the matrix A, which must be a square, symmetric, and positive-definite matrix.

  • Derived2 – The type of the vector B.

Parameters:
  • A – The coefficient matrix A.

  • b – The right-hand side vector B.

Throws:

(Implicitly) – Throws exceptions from Eigen if the LDLT decomposition fails due to reasons other than non-positive definiteness (e.g., memory allocation failure).

Returns:

An std::optional containing the solution vector X if the system is solvable (A is positive-definite), or std::nullopt if the system is not solvable (A is not positive-definite).

template<typename Scalar, int RowsAtCompileTime = Dynamic>
std::optional<Vector<Scalar, RowsAtCompileTime>> SolveLDLT(const SparseMatrix<Scalar> &A, const Vector<Scalar, RowsAtCompileTime> &b)

Solves the linear system A * X = B for X, where A is a symmetric positive-definite sparse matrix, return optional otherwise.

This function utilizes the LDLT decomposition (Cholesky decomposition) to efficiently solve the linear system. It assumes that A is a symmetric positive-definite matrix, which is a requirement for the LDLT decomposition.

Note

The function returns chol.solve(b), which implies that the solution returned is actually the solution of A * X = B.

Template Parameters:
  • Derived – The type of the matrix A, which must be a square, symmetric, and positive-definite matrix.

  • Derived2 – The type of the vector B.

Parameters:
  • A – The sparse coefficient matrix A.

  • b – The right-hand side vector B.

Throws:

(Implicitly) – Throws exceptions from Eigen if the LDLT decomposition fails due to reasons other than non-positive definiteness (e.g., memory allocation failure).

Returns:

An std::optional containing the solution vector X if the system is solvable (A is positive-definite), or std::nullopt if the system is not solvable (A is not positive-definite).

inline constexpr int SQRT(int N)

Integer square root function.

template<typename Scalar = double>
inline constexpr Scalar FloatEpsilon()
template<typename Scalar = double>
inline constexpr Scalar FloatEpsilon2()
inline std::nullptr_t &SuperNul()
template<typename Optimizer, typename X_t, typename Res_t>
inline auto Optimize(X_t &x, const Res_t &func, const typename Optimizer::Options &options = {})

Simplest interface to optimize x and minimize residuals (loss function). Internally call the optimizer and run the optimization.

inline auto tic()

Returns the current time point using the high-resolution clock.

This function provides a convenient way to capture the current time for benchmarking or timing purposes.

Returns:

The current time point as a std::chrono::time_point<Clock>.

inline double toc_ms(const std::chrono::time_point<Clock> &t0)

Calculates the elapsed time in milliseconds between a given time point and the current time.

This function measures the duration between a previously recorded time point (t0) and the current time, returning the result in milliseconds. It utilizes the high-resolution clock for precise timing.

Parameters:

t0 – The starting time point (std::chrono::time_point<Clock>).

Returns:

The elapsed time in milliseconds (double).

inline double dt_ms(const std::chrono::time_point<Clock> &t0, const std::chrono::time_point<Clock> &t1)

Calculates the elapsed time in milliseconds between two given time points.

This function measures the duration between two specified time points (t0 and t1), returning the result in milliseconds. It relies on the high-resolution clock for accurate timing.

Parameters:
  • t0 – The starting time point (std::chrono::time_point<Clock>).

  • t1 – The ending time point (std::chrono::time_point<Clock>).

Returns:

The elapsed time in milliseconds (double).

Variables

static constexpr int Dynamic = Eigen::Dynamic
static constexpr int Lower = Eigen::Lower
static constexpr int Upper = Eigen::Upper
static constexpr int Infinity = Eigen::Infinity
namespace diff

Typedefs

template<typename T, int N>
using Jet = ceres::Jet<T, N>

The Automatix differentiation Jet struct.

Enums

enum Method

Specifies the method used for numerical differentiation.

This enumeration controls the type of numerical differentiation performed. Numerical differentiation is used to approximate the derivative of a function at a given point by using finite differences.

Values:

enumerator kForward

Forward difference method.

Approximates the derivative using the function values at the current point and a point slightly ahead. It is a first-order approximation.

Formula: (f(x + h) - f(x)) / h

enumerator kCentral

Central difference method.

Approximates the derivative using the function values at points slightly before and after the current point. It is a second-order approximation, generally more accurate than the forward difference method.

Formula: (f(x + h) - f(x - h)) / (2 * h)

enumerator kFastCentral

Fast central difference method.

An optimized variant of the central difference method, potentially offering improved performance in certain scenarios, e.g when using a Manifold. In this case, accuracy will be traded with speed.

Formula: (f(xh) - f(xh - 2 * h)) / (2 * h), with xh = x+h.

Functions

template<typename X_t, typename Func>
auto CalculateJac(const X_t &x, const Func &cost)

Estimate the jacobian of d f(x)/d(x) around x using automatic differentiation.

template<typename X_t, typename Func>
auto EstimateNumJac(const X_t &x, const Func &f, const diff::Method &method = diff::Method::kCentral, typename traits::params_trait<X_t>::Scalar h = FloatEpsilon<typename traits::params_trait<X_t>::Scalar>())

Estimate the jacobian of d f(x)/d(x) around x using numerical differentiation.

template<typename X_t, typename ResidualsFunc>
auto CreateNumDiffFunc1(X_t&, const ResidualsFunc &residuals, const Method &method = Method::kCentral, typename traits::params_trait<X_t>::Scalar h = FloatEpsilon<typename traits::params_trait<X_t>::Scalar>())

Creates a numerical differentiation function for a given residuals function.

This function generates a callable object (std::function) that, when invoked, calculates the residuals and gradient of the provided residuals function at a given input x. It offers different numerical differentiation methods, such as forward, backward, and central differences.

// Example residuals function
auto loss = [](const std::vector<double>& x) {
return x[0] * x[0] + x[1];
};

std::vector<double> x = {1.0, 2.0};
auto acc_loss = CreateNumDiffFunc1(x, loss, Method::kCentral);

Eigen::Vector2d grad;

double norm = acc_loss(x, g, H);

The returned function can be passed to an optimizer, e.g.
auto optimizer = Optimizer<SolverGD<Vec2>>();
optimizer(x, acc_loss);

Note

The Method enum should be defined elsewhere and include values like Method::kForward, Method::kBackward, Method::kCentral and Method::kFastCentral.

Note

The step size used for numerical differentiation is determined internally and may be adapted based on the magnitude of the input x to avoid numerical instability.

Note

The generated std::function does not modify the input x.

Template Parameters:
  • X_t – The type of the input vector x. It should support arithmetic operations and element-wise access.

  • ResidualsFunc – The type of the residuals function. It should be a callable object that takes an X_t and returns a scalar value.

  • Scalar – The scalar type of the residuals, Jacobian, and Gradient.

  • Dims – The dimension of the input vector x.

Parameters:
  • residuals – The residuals function to be differentiated.

  • method – The numerical differentiation method to use. Defaults to Method::kCentral.

  • h – The delta added to β€˜x’ to compute the difference on each dimension

Returns:

A std::function object that takes an input x, a vector for the residuals, and a matrix for the Jacobian as arguments.

template<typename X_t, typename ResidualsFunc>
auto CreateNumDiffFunc2(X_t&, const ResidualsFunc &residuals, const Method &method = Method::kCentral, typename traits::params_trait<X_t>::Scalar h = FloatEpsilon<typename traits::params_trait<X_t>::Scalar>())

Creates a numerical differentiation function for a given residuals function.

This function generates a callable object (std::function) that, when invoked, calculates the residuals, gradient and Hessian of the provided residuals function at a given input x. It offers different numerical differentiation methods, such as forward, backward, and central differences.

// Example residuals function
auto loss = [](const std::vector<double>& x) {
return x[0] * x[0] + x[1];
};

std::vector<double> x = {1.0, 2.0};
auto acc_loss = CreateNumDiffFunc1(x, loss, Method::kCentral);

Eigen::Vector2d grad;
Eigen::Matrix2d H;

double norm = acc_loss(x, g, H);

The returned function can be passed to an optimizer, e.g.
auto optimizer = Optimizer<SolverLM<Mat2>>();
optimizer(x, acc_loss);

Note

The Method enum should be defined elsewhere and include values like Method::kCentral, Method::kBackward, Method::kCentral and Method::kFastCentral.

Note

The step size used for numerical differentiation is determined internally and may be adapted based on the magnitude of the input x to avoid numerical instability.

Note

The generated std::function does not modify the input x.

Template Parameters:
  • X_t – The type of the input vector x. It should support arithmetic operations and element-wise access.

  • ResidualsFunc – The type of the residuals function. It should be a callable object that takes an X_t and returns a scalar value.

  • Scalar – The scalar type of the residuals, Jacobian, and Hessian.

  • Dims – The dimension of the input vector x.

Parameters:
  • residuals – The residuals function to be differentiated.

  • method – The numerical differentiation method to use. Defaults to Method::kCentral.

  • h – The delta added to β€˜x’ to compute the difference on each dimension.

Returns:

A std::function object that takes an input x, a vector for the residuals, and a matrix for the Jacobian as arguments. It also takes an optional matrix for the Hessian as an argument.

namespace distances

Functions

template<typename TA, typename TB, typename ExportJ = std::nullptr_t>
auto Euclidean(const TA &a, const TB &b, const ExportJ &export_jac = nullptr)

Compute the Euclidean L2 distance (a.k.a L2) between a and b: d(a,b) = ||a-b||.

template<typename TA, typename TB, typename ExportJ = std::nullptr_t>
auto L2(const TA &a, const TB &b, const ExportJ &export_jac = nullptr)
template<typename TA, typename TB, typename ExportJ = std::nullptr_t>
auto Manhattan(const TA &a, const TB &b, const ExportJ &export_jac = nullptr)

Compute the Manhattan distance (a.k.a L1) between a and b: d(a,b) = |a-b|.

template<typename TA, typename TB, typename ExportJ = std::nullptr_t>
auto L1(const TA &a, const TB &b, const ExportJ &export_jac = nullptr)
template<typename TA, typename TB, typename ExportJ = std::nullptr_t>
auto Linf(const TA &a, const TB &b, const ExportJ &export_jac = nullptr)

Compute the L-infinity distance (a.k.a max(x)) between a and b: d(a,b) = max(a-b)

template<typename TA, typename TB, typename ExportJ = std::nullptr_t>
auto Cosine(const TA &a, const TB &b, const ExportJ& = nullptr)

Compute the cosine distance between a and b: d(a,b) = a ∠ b.

template<typename TA, typename TB, typename Cov_t, typename ExportJ = std::nullptr_t>
auto MahaNorm(const TA &a, const TB &b, const Cov_t &cov_or_var, const ExportJ &export_jac = nullptr)

Compute the Mahalanobis distance between a and b with a covariance cov: d(a,b) = ||a-b||Ξ£

namespace gd

Gradient Descent specific solver, optimizer and their options.

Typedefs

template<typename Gradient_t>
using Solver = solvers::SolverGD<Gradient_t>

Gradient Descent Solver.

template<typename Gradient_t>
using Optimizer = optimizers::Optimizer<Solver<Gradient_t>, Options>

Gradient Descent Optimizater type.

Functions

template<typename X_t, typename Res_t>
inline auto Optimize(X_t &x, const Res_t &func, const Options &options = Options())

Gradient Descent Optimize function.

namespace guc
namespace losses

If you land here, put your shades one because those fat functions are still hurting my eyes as of today.

Activation losses

DefineLoss (Sigmoid, 1.0/(1.0+exp(-x)), l.unaryExpr([](auto x) -> Scalar { return x *(Scalar(1.0) - x);}).reshaped().asDiagonal())

Sigmoid = 1/(1+e^-x), derivative = Sigmoid(x) * (1 - Sigmoid(x))

DefineLoss (Tanh,(exp(x) - exp(-x))/(exp(x)+exp(-x)), l.unaryExpr([](auto x) -> Scalar { return Scalar(1.0) - x *x;}).reshaped().asDiagonal())

Tanh = (e^x-e^-x)/(e^x+e^-x), derivative = 1 - Tanh(x)^2.

DefineLoss (ReLU, x > 0 ? x :Scalar(0.0), x.unaryExpr([](auto x) { return x > 0 ? Scalar(1.0) :Scalar(0.0);}).reshaped().asDiagonal())

ReLU = max(0, x), derivative = {x>0:1, x<=0:0}.

DefineLoss2 (LeakyReLU, x > 0 ? x :a *x, x.unaryExpr([a](auto x) { return x > 0 ? 1 :a;}).reshaped().asDiagonal())

LeakyReLU = {x>0:x, x<=0:a*x}, derivative = {x>0:1, x<=0:a}.

Classification losses

template<typename T, typename ExportJ = std::nullptr_t>
auto Softmax(const T &x, const ExportJ &Jx_or_bool = nullptr)

Softmax = e^xi / sum(e^x), jacobian = {i=j: si(x)*(1-si(x)) , i!=j: -si(x)*sj(x)}.

template<typename T, typename ExportJ = std::nullptr_t>
auto SafeSoftmax(const T &x, const ExportJ &Jx_or_bool = nullptr)

Safe Softmax = e^(xi-max(xi)) / sum(e^(x-max(xi))

Mahalanobis Distances

template<typename T, typename Cov_t, typename ExportJ = std::nullptr_t>
auto MahaSquaredNorm(const T &x, const Cov_t &cov_or_var, const ExportJ &Jx_or_bool = nullptr, bool add_scale = true)

Compute the Squared Mahalanobis distance of ´x´ with a covariance cov: n(x) = ||x||Σ²

template<typename T, typename Cov_t, typename ExportJ = std::nullptr_t>
auto MahaNorm(const T &x, const Cov_t &cov_or_var, const ExportJ &Jx_or_bool = nullptr)

Compute the Mahalanobis distance of Β΄xΒ΄ with a covariance cov: n(x) = ||x||Ξ£

template<typename Derived, typename Cov_t, typename ExportJ = std::nullptr_t>
auto MahaWhitened(const MatrixBase<Derived> &res, const Cov_t &cov_stevs, const ExportJ &Jx_or_bool = nullptr)

Return the Whitened/Sphered/Decorrelated residuals (and its jacobian as a option) when applying a Mahalanobis norm when given residuals and a covariance matrix (Upper filled at least)

template<typename Derived, typename DerivedC, typename ExportJ = std::nullptr_t>
auto MahaWhitenedInfoU(const MatrixBase<Derived> &res, const MatrixBase<DerivedC> &U, const ExportJ &Jx_or_bool = nullptr)

Return the Whitened/Sphered/Decorrelated residuals (and its jacobian as a option) when applying a Mahalanobis norm when given residuals and a Upper Information matrix

Classical Norms

template<typename T>
auto SquaredL2(const T &x)

Return the squared L2 norm of a vector or scalar (a.k.a Sum of Squares)

template<typename T, typename ExportJ>
auto SquaredL2(const T &x, const ExportJ &Jx_or_bool, bool add_scale = true)

Return the squared L2 norm of a vector or scalar and its jacobian (x.t())

template<typename T>
auto L2(const T &x)

Return L2 norm of a vector or scalar.

template<typename T, typename ExportJ>
auto L2(const T &x, const ExportJ &Jx_or_bool)

Return L2 norm of a vector or scalar and its jacobian (x.t() / ||x||)

template<typename T>
auto L1(const T &x)

Return L1 norm of a vector or scalar.

template<typename T, typename ExportJ>
auto L1(const T &x, const ExportJ &Jx_or_bool)

Return L1 norm of a vector or scalar and its jacobian.

template<typename T>
auto Linf(const T &x)

Return L infinity norm of a vector or scalar.

template<typename T, typename ExportJ>
auto Linf(const T &x, const ExportJ &Jx_or_bool)

Return L infinity norm of a vector or scalar.

M-Estimators / Robust Losses

These functions return a loss and a scale to apply on the jacobian or the scaled jacobian

The loss’s gradient is then s*res.t(), with β€˜s’ the scale returned by the robust function (e.g. Truncated, huber, Artcan, etc.), typically 1 if the ||res|| < threshold. The scale can then be used to solve the lossal equations: JtJ * dx = Jt*res*s Ex: Huber : s = {1, th/||x||} so that we thereby clip the gradient to a max of theshold β€˜th’ while keeping the direction. grdaient = {inlier:Jt * res; outlier:Jt * res * th / ||res||}

template<typename T, typename ExportJ = std::nullptr_t>
auto Truncated(const T &n2, typename traits::params_trait<T>::Scalar th2, const ExportJ &Jx_or_bool = nullptr)

Hard Truncation: Clip the given loss n to a max of th, also return the scale {0 1} or scaled jacobian if given as input Jx_or_bool

template<typename T, typename ExportJ = std::nullptr_t>
auto TruncatedLoss(const T &x, typename traits::params_trait<T>::Scalar th2, const ExportJ &Jx_or_bool = nullptr)

Hard Truncated L2 Loss: Clip the L2 loss of x to a max of th, also return the scale {0, 1} or scaled jacobian if given as input Jx_or_bool

template<typename T, typename ExportJ = std::nullptr_t>
auto Huber(const T &n2, typename traits::params_trait<T>::Scalar th2, const ExportJ &Jx_or_bool = nullptr)

Huber: Return a scaled loss \( n'= n if n < th, else n'= \sqrt(2.0 * th * n - thΒ²) \), also return the scale {1, th/n} or scaled jacobian if given as input Jx_or_bool

template<typename T, typename ExportJ = std::nullptr_t>
auto HuberLoss(const T &x, typename traits::params_trait<T>::Scalar th2, const ExportJ &Jx_or_bool = nullptr)

Huber: Return a scaled loss of β€˜x’ such that \( n'Β² = ||x||Β² if ||x||Β² < th2, else n'Β² = 2.0 * th * n - thΒ² \), also return the scale {1, th/n} or scaled jacobian if given as input Jx_or_bool

Note

The gradient is then {inlier:x.t(), outlier:x.t() * s}, with \(s = th / ||x||\).

template<typename T, typename ExportJ = std::nullptr_t>
auto Tukey(const T &n2, typename traits::params_trait<T>::Scalar th2, const ExportJ &Jx_or_bool = nullptr)

Tukey: Return a scaled loss \( n'Β² = nΒ² if n < th, else n'Β² = 2.0 * h * n - thΒ² \), also return the scale or scaled jacobian if given as input Jx_or_bool

template<typename T, typename ExportJ = std::nullptr_t>
auto TukeyLoss(const T &x, typename traits::params_trait<T>::Scalar th2, const ExportJ &Jx_or_bool = nullptr)

Tukey: Return a scaled loss of β€˜x’ such that \( n'Β² = ||x||Β² if ||x|| < th, else n'Β² = 2.0 * h * n - thΒ² \), also return the scale or scaled jacobian if given as input Jx_or_bool

template<typename T, typename ExportJ = std::nullptr_t>
auto Arctan(const T &n2, typename traits::params_trait<T>::Scalar th2, const ExportJ &Jx_or_bool = nullptr)

Arctan: Return a scaled loss \( n'Β² = th * \atan(nΒ² / th) \), also return the scale or scaled jacobian if given as input Jx_or_bool

template<typename T, typename ExportJ = std::nullptr_t>
auto ArctanLoss(const T &x, typename traits::params_trait<T>::Scalar th2, const ExportJ &Jx_or_bool = nullptr)

Arctan: Return a scaled loss of β€˜x’ such that \( n'Β² = th * \atan(||x||Β² / th) \), also return the scale or scaled jacobian if given as input Jx_or_bool

template<typename T, typename ExportJ = std::nullptr_t>
auto Cauchy(const T &n2, typename traits::params_trait<T>::Scalar th2, const ExportJ &Jx_or_bool = nullptr)

Cauchy: Return a scaled loss \( n'Β² = thΒ² * \log(1 + nΒ² / thΒ²) \), also return the scale or scaled jacobian if given as input Jx_or_bool

template<typename T, typename ExportJ = std::nullptr_t>
auto CauchyLoss(const T &x, typename traits::params_trait<T>::Scalar th2, const ExportJ &Jx_or_bool = nullptr)

Cauchy: Return a scaled loss of β€˜x’ such that \( n'Β² = thΒ² * \log(1 + ||x||Β² / thΒ²) \), also return the scale or scaled jacobian if given as input Jx_or_bool

template<typename T, typename ExportJ = std::nullptr_t>
auto GemanMcClure(const T &n2, typename traits::params_trait<T>::Scalar th2, const ExportJ &Jx_or_bool = nullptr)

GemanMcClure: Return a scaled loss \( n'Β² = nΒ² / (||x||Β² + thΒ²) \), also return the scale or scaled jacobian if given as input Jx_or_bool

template<typename T, typename ExportJ = std::nullptr_t>
auto GemanMcClureLoss(const T &x, typename traits::params_trait<T>::Scalar th2, const ExportJ &Jx_or_bool = nullptr)

GemanMcClure: Return a scaled loss of β€˜x’ such that \( n' = ||x|| / \sqrt(||x||Β² + thΒ²) \), also return the scale or scaled jacobian if given as input Jx_or_bool

template<typename T, typename ExportJ = std::nullptr_t>
auto BlakeZisserman(const T &n2, typename traits::params_trait<T>::Scalar th2, const ExportJ &Jx_or_bool = nullptr)

BlakeZisserman: Return a scaled loss nsuch that @f$ n'Β² = -\log(\exp(-nΒ²) + \exp(-thΒ²)) @f$, also return the scale or scaled jacobian if given as inputJx_or_bool`.

template<typename T, typename ExportJ = std::nullptr_t>
auto BlakeZissermanLoss(const T &x, typename traits::params_trait<T>::Scalar th2, const ExportJ &Jx_or_bool = nullptr)

BlakeZisserman: Return a scaled loss of β€˜x’ such that \( n'Β² = -\log(\exp(-||x||Β²) + \exp(-thΒ²)) \), also return the scale or scaled jacobian if given as input Jx_or_bool

namespace nlls
namespace gn

Gauss-Newton specific solver, optimizer and their options.

Typedefs

template<typename Hessian_t>
using Solver = solvers::SolverGN<Hessian_t>

Gauss-Newton Solver.

template<typename Hessian_t = SparseMatrix<double>>
using SparseSolver = solvers::SolverGN<Hessian_t>

Gauss-Newton Sparse Solver.

template<typename Hessian_t>
using Optimizer = optimizers::Optimizer<Solver<Hessian_t>, Options>

Gauss-Newton Optimizater type.

using SolverOptions = solvers::Options2

Functions

template<typename X_t, typename Res_t>
inline auto Optimize(X_t &x, const Res_t &func, const Options &options = Options())

Gauss-Newton Optimize function.

namespace lm

Levenberg-Marquardt specific solver, optimizer and their options.

Typedefs

template<typename Hessian_t>
using Solver = solvers::SolverLM<Hessian_t>

Levenberg-Marquardt Solver.

template<typename Hessian_t = SparseMatrix<double>>
using SparseSolver = solvers::SolverLM<Hessian_t>

Levenberg-Marquardt Sparse Solver.

template<typename Hessian_t>
using Optimizer = optimizers::Optimizer<Solver<Hessian_t>, Options>

Levenberg-Marquardt Optimizater type.

Functions

template<typename X_t, typename Res_t>
inline auto Optimize(X_t &x, const Res_t &func, const Options &options = Options())

Levenberg-Marquardt Optimize function.

namespace optimizers
namespace solvers
namespace traits

Variables

template<typename T>
constexpr bool is_jet_type_v = is_jet_type<T>::value
template<typename T>
constexpr bool is_nullptr_v = is_nullptr_t<T>::value
template<typename T>
constexpr bool is_bool_v = is_bool<T>::value
template<typename T>
constexpr bool is_scalar_v = is_scalar<T>::value
template<typename T>
constexpr bool is_pair_v = is_pair<std::decay_t<T>>::value
template<typename T>
constexpr bool is_matrix_or_array_v = is_matrix_or_array<std::decay_t<T>>::value
template<typename T>
constexpr bool is_sparse_matrix_v = is_sparse_matrix<std::decay_t<T>>::value
template<typename T>
constexpr bool is_matrix_or_scalar_v = (std::is_scalar_v<T> && !std::is_same_v<T, bool>) || is_sparse_matrix_v<T> || is_matrix_or_array_v<T>
template<typename T>
constexpr bool is_streamable_v = is_streamable<T>::value
file auto_diff.h
#include <tinyopt/diff/jet.h>
#include <tinyopt/math.h>
#include <type_traits>
file jet.h
#include <tinyopt/3rdparty/ceres/jet.h>
#include <tinyopt/math.h>
file jet_traits.h
#include <tinyopt/diff/jet.h>
#include <tinyopt/traits.h>
file num_diff.h
#include <tinyopt/math.h>
#include <tinyopt/traits.h>
file optimize_autodiff.h
#include <cassert>
#include <cstddef>
#include <type_traits>
#include <utility>
#include <tinyopt/math.h>
#include <tinyopt/traits.h>
#include <tinyopt/diff/jet.h>
file distances.h
#include <tinyopt/math.h>
#include <tinyopt/traits.h>
#include <type_traits>
file formatters.h
#include <format>
#include <sstream>
#include <tinyopt/math.h>
file log.h
#include <iostream>
#include <sstream>
#include <stdexcept>
#include <string>
#include <vector>
#include <tinyopt/traits.h>
#include β€œtinyopt/formatters.h”

Defines

TINYOPT_LOG(...)
TINYOPT_FORMAT_NS
TINYOPT_LOG_MAT(m)
file activations.h
file classif.h
#include <tinyopt/math.h>
#include <tinyopt/traits.h>
file helpers.h
#include <tinyopt/math.h>
#include <tinyopt/traits.h>

Defines

DefineLoss(name, code, jac_code)
DefineLoss2(name, code, jac_code)
file losses.h
file mahalanobis.h
#include <tinyopt/math.h>
#include <tinyopt/traits.h>
file norms.h
#include <tinyopt/math.h>
#include <tinyopt/traits.h>
file robust_norms.h
#include <tinyopt/math.h>
#include <tinyopt/traits.h>
file math.h
#include <Eigen/Core>
#include <Eigen/Dense>
#include <Eigen/Sparse>
#include <Eigen/src/SparseCore/SparseSelfAdjointView.h>
#include <Eigen/SparseCholesky>
#include <cstddef>
#include <optional>
file optimize.h
file nlls.h
file optimizer.h
#include <cmath>
#include <cstddef>
#include <cstdint>
#include <limits>
#include <optional>
#include <variant>
#include β€œtinyopt/math.h”
#include <tinyopt/log.h>
#include <tinyopt/output.h>
#include <tinyopt/time.h>
#include <tinyopt/traits.h>
file optimizers.h
file unconstrained.h
file output.h
#include <cstddef>
#include <cstdint>
#include <optional>
#include <vector>
#include <Eigen/Core>
#include <tinyopt/math.h>
#include <tinyopt/time.h>
#include <tinyopt/traits.h>
file base.h
#include <algorithm>
#include <cassert>
#include <limits>
#include <type_traits>
#include <tinyopt/log.h>
#include <tinyopt/math.h>
#include <tinyopt/output.h>
file gd.h
#include <tinyopt/math.h>
#include <tinyopt/optimize.h>
file gd.h
#include β€œtinyopt/traits.h”
file gn.h
#include <tinyopt/math.h>
#include <tinyopt/optimize.h>
file gn.h
#include <type_traits>
#include β€œtinyopt/traits.h”
file lm.h
#include <tinyopt/math.h>
#include <tinyopt/optimize.h>
#include <type_traits>
file lm.h
#include <cstddef>
#include <stdexcept>
#include β€œtinyopt/math.h”
#include <tinyopt/log.h>
file options.h
#include <cstdint>
#include <functional>
#include β€œtinyopt/math.h”
#include <tinyopt/log.h>
file options.h
file solvers.h
file stop_reasons.h
#include <ostream>
#include <sstream>
#include <tinyopt/traits.h>
file time.h
#include <chrono>
file tinyopt.h
#include <tinyopt/math.h>
#include <tinyopt/traits.h>
file traits.h
#include <type_traits>
#include <tinyopt/math.h>
dir include/tinyopt/diff
dir include
dir include/tinyopt/losses
dir include/tinyopt/optimizers
dir include/tinyopt/solvers
dir include/tinyopt

Indices and tables