Janus 2.0.0
High-performance C++20 dual-mode numerical framework
Loading...
Searching...
No Matches
OptiSol.hpp
Go to the documentation of this file.
1
5
6#pragma once
7
10#include <casadi/casadi.hpp>
11#include <map>
12#include <optional>
13#include <string>
14#include <vector>
15
16namespace janus {
17
32class OptiSol {
33 public:
38 explicit OptiSol(casadi::OptiSol cas_sol) : cas_sol_(std::move(cas_sol)) {}
39
45 double value(const SymbolicScalar &var) const {
46 casadi::DM result = cas_sol_.value(var);
47 return static_cast<double>(result);
48 }
49
56 // Convert SymbolicVector (Eigen<MX>) to MX for CasADi
57 casadi::MX mx_var = janus::to_mx(var);
58 casadi::DM result = cas_sol_.value(mx_var);
59
60 // Convert DM to NumericVector
61 std::vector<double> elements = static_cast<std::vector<double>>(result);
62 NumericVector vec(elements.size());
63 for (size_t i = 0; i < elements.size(); ++i) {
64 vec(static_cast<Eigen::Index>(i)) = elements[i];
65 }
66 return vec;
67 }
68
75 casadi::MX mx_var = janus::to_mx(var);
76 casadi::DM result = cas_sol_.value(mx_var);
77
78 std::vector<double> elements = static_cast<std::vector<double>>(result);
79 NumericMatrix mat(result.size1(), result.size2());
80
81 // CasADi uses column-major (Fortran) order
82 for (casadi_int j = 0; j < result.size2(); ++j) {
83 for (casadi_int i = 0; i < result.size1(); ++i) {
84 mat(i, j) = elements[j * result.size1() + i];
85 }
86 }
87 return mat;
88 }
89
94 casadi::Dict stats() const { return cas_sol_.stats(); }
95
100 std::optional<int> num_function_evals() const {
101 auto s = stats();
102 if (s.count("n_call_nlp_f")) {
103 return static_cast<int>(s.at("n_call_nlp_f"));
104 }
105 return std::nullopt;
106 }
107
112 std::optional<int> num_iterations() const {
113 auto s = stats();
114 if (s.count("iter_count")) {
115 return static_cast<int>(s.at("iter_count"));
116 }
117 return std::nullopt;
118 }
119
124 const casadi::OptiSol &casadi_sol() const { return cas_sol_; }
125
132 void save(const std::string &filename,
133 const std::map<std::string, SymbolicScalar> &named_vars) const {
134 std::map<std::string, std::vector<double>> data;
135
136 for (const auto &[name, var] : named_vars) {
137 // Evaluate variable using helper
138 double val = value(var);
139 data[name] = {val};
140 }
141
142 janus::utils::write_json(filename, data);
143 }
144
150 void save(const std::string &filename,
151 const std::map<std::string, SymbolicVector> &named_vars) const {
152 std::map<std::string, std::vector<double>> data;
153
154 for (const auto &[name, var] : named_vars) {
155 // Evaluate variable using helper
156 NumericVector result = value(var);
157
158 // Convert Eigen Vector to std::vector
159 std::vector<double> elements(result.data(), result.data() + result.size());
160 data[name] = elements;
161 }
162
163 janus::utils::write_json(filename, data);
164 }
165
173 static std::map<std::string, std::vector<double>> load(const std::string &filename) {
174 return janus::utils::read_json(filename);
175 }
176
177 private:
178 casadi::OptiSol cas_sol_;
179};
180
181} // namespace janus
Core type aliases for numeric and symbolic Eigen/CasADi interop.
Simple JSON read/write utilities for flat string-to-vector maps.
casadi::Dict stats() const
Get solver statistics.
Definition OptiSol.hpp:94
NumericVector value(const SymbolicVector &var) const
Extract vector value at optimum.
Definition OptiSol.hpp:55
static std::map< std::string, std::vector< double > > load(const std::string &filename)
Load solution data from JSON file.
Definition OptiSol.hpp:173
double value(const SymbolicScalar &var) const
Extract scalar value at optimum.
Definition OptiSol.hpp:45
std::optional< int > num_iterations() const
Get number of iterations.
Definition OptiSol.hpp:112
void save(const std::string &filename, const std::map< std::string, SymbolicVector > &named_vars) const
Save solution map of vectors to JSON file.
Definition OptiSol.hpp:150
OptiSol(casadi::OptiSol cas_sol)
Construct from CasADi solution.
Definition OptiSol.hpp:38
void save(const std::string &filename, const std::map< std::string, SymbolicScalar > &named_vars) const
Save solution to JSON file.
Definition OptiSol.hpp:132
std::optional< int > num_function_evals() const
Get number of function evaluations.
Definition OptiSol.hpp:100
const casadi::OptiSol & casadi_sol() const
Access underlying CasADi solution.
Definition OptiSol.hpp:124
NumericMatrix value(const SymbolicMatrix &var) const
Extract matrix value at optimum.
Definition OptiSol.hpp:74
std::map< std::string, std::vector< double > > read_json(const std::string &filename)
Parse a flat JSON file into a string-to-vector<double> map.
Definition JsonUtils.hpp:55
void write_json(const std::string &filename, const std::map< std::string, std::vector< double > > &data)
Write a map of string-to-vector<double> as JSON.
Definition JsonUtils.hpp:22
Definition Diagnostics.hpp:19
JanusVector< SymbolicScalar > SymbolicVector
Eigen vector of MX elements.
Definition JanusTypes.hpp:72
JanusMatrix< NumericScalar > NumericMatrix
Eigen::MatrixXd equivalent.
Definition JanusTypes.hpp:66
JanusVector< NumericScalar > NumericVector
Eigen::VectorXd equivalent.
Definition JanusTypes.hpp:67
casadi::MX to_mx(const Eigen::MatrixBase< Derived > &e)
Convert Eigen matrix of MX (or numeric) to CasADi MX.
Definition JanusTypes.hpp:189
JanusMatrix< SymbolicScalar > SymbolicMatrix
Eigen matrix of MX elements.
Definition JanusTypes.hpp:71
casadi::MX SymbolicScalar
CasADi MX symbolic scalar.
Definition JanusTypes.hpp:70