25 #ifndef SRC_NUMERICS_INCLUDE_CORBO_NUMERICS_FINITE_DIFFERENCES_INTERFACE_H_ 26 #define SRC_NUMERICS_INCLUDE_CORBO_NUMERICS_FINITE_DIFFERENCES_INTERFACE_H_ 32 #ifdef MESSAGE_SUPPORT 33 #include <corbo-communication/messages/numerics/finite_differences.pb.h> 63 using Ptr = std::shared_ptr<FiniteDifferencesInterface>;
64 using UPtr = std::unique_ptr<FiniteDifferencesInterface>;
145 virtual void computeJacobian2(std::function<
void(
int,
const double&)> inc_fun, std::function<
void(Eigen::VectorXd&)> eval_fun,
196 virtual void computeHessian2(std::function<
void(
int,
const double&)> inc_fun, std::function<
void(Eigen::VectorXd&)> eval_fun,
int dim_f,
230 const double* multipliers =
nullptr) = 0;
249 virtual void computeJacobianAndHessian2(std::function<
void(
int,
const double&)> inc_fun, std::function<
void(Eigen::VectorXd&)> eval_fun,
251 const double* multipliers =
nullptr) = 0;
253 #ifdef MESSAGE_SUPPORT 254 virtual void toMessage(messages::FiniteDifferences& message)
const {}
257 virtual void fromMessage(
const messages::FiniteDifferences& message, std::stringstream* issues =
nullptr) {}
262 #define FACTORY_REGISTER_FINITE_DIFFERENCES(type) FACTORY_REGISTER_OBJECT(type, FiniteDifferencesInterface) 266 #endif // SRC_NUMERICS_INCLUDE_CORBO_NUMERICS_FINITE_DIFFERENCES_INTERFACE_H_
Eigen::VectorXd StateVector
virtual FiniteDifferencesInterface::Ptr getInstance() const =0
Return a newly allocated instances of the inherited class.
virtual void computeJacobianAndHessian2(std::function< void(int, const double &)> inc_fun, std::function< void(Eigen::VectorXd &)> eval_fun, Eigen::Ref< Eigen::MatrixXd > jacobian, Eigen::Ref< Eigen::MatrixXd > hessian, const double *multipliers=nullptr)=0
Compute Jacobian and Hessian of a desired function (overload which accepts a slightly different callb...
std::shared_ptr< FiniteDifferencesInterface > Ptr
virtual ~FiniteDifferencesInterface()
Virtual destructor.
virtual void computeJacobian2(std::function< void(int, const double &)> inc_fun, std::function< void(Eigen::VectorXd &)> eval_fun, Eigen::Ref< Eigen::MatrixXd > jacobian)=0
Compute Jacobian of a desired function (overload which accepts a slightly different callback function...
std::unique_ptr< FiniteDifferencesInterface > UPtr
A matrix or vector expression mapping an existing expression.
virtual void computeJacobian(std::function< void(int, const double &)> inc_fun, std::function< void(Eigen::Ref< Eigen::VectorXd >)> eval_fun, Eigen::Ref< Eigen::MatrixXd > jacobian)=0
Compute Jacobian of a desired function.
virtual void computeHessian2(std::function< void(int, const double &)> inc_fun, std::function< void(Eigen::VectorXd &)> eval_fun, int dim_f, Eigen::Ref< Eigen::MatrixXd > hessian, const double *multipliers=nullptr)=0
Compute Hessian of a desired function (overload which accepts a slightly different callback function)...
virtual void computeHessian(std::function< void(int, const double &)> inc_fun, std::function< void(Eigen::Ref< Eigen::VectorXd >)> eval_fun, int dim_f, Eigen::Ref< Eigen::MatrixXd > hessian, const double *multipliers=nullptr)=0
Compute Hessian of a desired function.
virtual void computeJacobianAndHessian(std::function< void(int, const double &)> inc_fun, std::function< void(Eigen::Ref< Eigen::VectorXd >)> eval_fun, Eigen::Ref< Eigen::MatrixXd > jacobian, Eigen::Ref< Eigen::MatrixXd > hessian, const double *multipliers=nullptr)=0
Compute Jacobian and Hessian of a desired function.
Eigen::VectorXd InputVector
Interface class for finite difference approaches.