123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185 |
- #include "ceres/c_api.h"
- #include <iostream>
- #include <memory>
- #include <string>
- #include <vector>
- #include "ceres/cost_function.h"
- #include "ceres/loss_function.h"
- #include "ceres/problem.h"
- #include "ceres/solver.h"
- #include "ceres/types.h"
- #include "glog/logging.h"
- using ceres::Problem;
- void ceres_init() {
-
-
- char message[] = "<unknown>";
- google::InitGoogleLogging(message);
- }
- ceres_problem_t* ceres_create_problem() {
- return reinterpret_cast<ceres_problem_t*>(new Problem);
- }
- void ceres_free_problem(ceres_problem_t* problem) {
- delete reinterpret_cast<Problem*>(problem);
- }
- class CERES_NO_EXPORT CallbackCostFunction final : public ceres::CostFunction {
- public:
- CallbackCostFunction(ceres_cost_function_t cost_function,
- void* user_data,
- int num_residuals,
- int num_parameter_blocks,
- int* parameter_block_sizes)
- : cost_function_(cost_function), user_data_(user_data) {
- set_num_residuals(num_residuals);
- for (int i = 0; i < num_parameter_blocks; ++i) {
- mutable_parameter_block_sizes()->push_back(parameter_block_sizes[i]);
- }
- }
- bool Evaluate(double const* const* parameters,
- double* residuals,
- double** jacobians) const final {
- return (*cost_function_)(
- user_data_, const_cast<double**>(parameters), residuals, jacobians);
- }
- private:
- ceres_cost_function_t cost_function_;
- void* user_data_;
- };
- class CallbackLossFunction final : public ceres::LossFunction {
- public:
- explicit CallbackLossFunction(ceres_loss_function_t loss_function,
- void* user_data)
- : loss_function_(loss_function), user_data_(user_data) {}
- void Evaluate(double sq_norm, double* rho) const final {
- (*loss_function_)(user_data_, sq_norm, rho);
- }
- private:
- ceres_loss_function_t loss_function_;
- void* user_data_;
- };
- void* ceres_create_huber_loss_function_data(double a) {
- return new ceres::HuberLoss(a);
- }
- void* ceres_create_softl1_loss_function_data(double a) {
- return new ceres::SoftLOneLoss(a);
- }
- void* ceres_create_cauchy_loss_function_data(double a) {
- return new ceres::CauchyLoss(a);
- }
- void* ceres_create_arctan_loss_function_data(double a) {
- return new ceres::ArctanLoss(a);
- }
- void* ceres_create_tolerant_loss_function_data(double a, double b) {
- return new ceres::TolerantLoss(a, b);
- }
- void ceres_free_stock_loss_function_data(void* loss_function_data) {
- delete reinterpret_cast<ceres::LossFunction*>(loss_function_data);
- }
- void ceres_stock_loss_function(void* user_data,
- double squared_norm,
- double out[3]) {
- reinterpret_cast<ceres::LossFunction*>(user_data)->Evaluate(squared_norm,
- out);
- }
- ceres_residual_block_id_t* ceres_problem_add_residual_block(
- ceres_problem_t* problem,
- ceres_cost_function_t cost_function,
- void* cost_function_data,
- ceres_loss_function_t loss_function,
- void* loss_function_data,
- int num_residuals,
- int num_parameter_blocks,
- int* parameter_block_sizes,
- double** parameters) {
- auto* ceres_problem = reinterpret_cast<Problem*>(problem);
- auto callback_cost_function =
- std::make_unique<CallbackCostFunction>(cost_function,
- cost_function_data,
- num_residuals,
- num_parameter_blocks,
- parameter_block_sizes);
- std::unique_ptr<ceres::LossFunction> callback_loss_function;
- if (loss_function != nullptr) {
- callback_loss_function = std::make_unique<CallbackLossFunction>(
- loss_function, loss_function_data);
- }
- std::vector<double*> parameter_blocks(parameters,
- parameters + num_parameter_blocks);
- return reinterpret_cast<ceres_residual_block_id_t*>(
- ceres_problem->AddResidualBlock(callback_cost_function.release(),
- callback_loss_function.release(),
- parameter_blocks));
- }
- void ceres_solve(ceres_problem_t* c_problem) {
- auto* problem = reinterpret_cast<Problem*>(c_problem);
-
-
-
- ceres::Solver::Options options;
- options.max_num_iterations = 100;
- options.linear_solver_type = ceres::DENSE_QR;
- options.minimizer_progress_to_stdout = true;
- ceres::Solver::Summary summary;
- ceres::Solve(options, problem, &summary);
- std::cout << summary.FullReport() << "\n";
- }
|