Stan  2.10.0
probability, sampling & optimization
Public Types | Public Member Functions | List of all members
stan::optimization::BFGSUpdate_HInv< Scalar, DimAtCompile > Class Template Reference

#include <bfgs_update.hpp>

Public Types

typedef Eigen::Matrix< Scalar, DimAtCompile, 1 > VectorT
 
typedef Eigen::Matrix< Scalar, DimAtCompile, DimAtCompile > HessianT
 

Public Member Functions

Scalar update (const VectorT &yk, const VectorT &sk, bool reset=false)
 Update the inverse Hessian approximation. More...
 
void search_direction (VectorT &pk, const VectorT &gk) const
 Compute the search direction based on the current (inverse) Hessian approximation and given gradient. More...
 

Detailed Description

template<typename Scalar = double, int DimAtCompile = Eigen::Dynamic>
class stan::optimization::BFGSUpdate_HInv< Scalar, DimAtCompile >

Definition at line 10 of file bfgs_update.hpp.

Member Typedef Documentation

template<typename Scalar = double, int DimAtCompile = Eigen::Dynamic>
typedef Eigen::Matrix<Scalar, DimAtCompile, DimAtCompile> stan::optimization::BFGSUpdate_HInv< Scalar, DimAtCompile >::HessianT

Definition at line 13 of file bfgs_update.hpp.

template<typename Scalar = double, int DimAtCompile = Eigen::Dynamic>
typedef Eigen::Matrix<Scalar, DimAtCompile, 1> stan::optimization::BFGSUpdate_HInv< Scalar, DimAtCompile >::VectorT

Definition at line 12 of file bfgs_update.hpp.

Member Function Documentation

template<typename Scalar = double, int DimAtCompile = Eigen::Dynamic>
void stan::optimization::BFGSUpdate_HInv< Scalar, DimAtCompile >::search_direction ( VectorT pk,
const VectorT gk 
) const
inline

Compute the search direction based on the current (inverse) Hessian approximation and given gradient.

Parameters
[out]pkThe negative product of the inverse Hessian and gradient direction gk.
[in]gkGradient direction.

Definition at line 55 of file bfgs_update.hpp.

template<typename Scalar = double, int DimAtCompile = Eigen::Dynamic>
Scalar stan::optimization::BFGSUpdate_HInv< Scalar, DimAtCompile >::update ( const VectorT yk,
const VectorT sk,
bool  reset = false 
)
inline

Update the inverse Hessian approximation.

Parameters
ykDifference between the current and previous gradient vector.
skDifference between the current and previous state vector.
resetWhether to reset the approximation, forgetting about previous values.
Returns
In the case of a reset, returns the optimal scaling of the initial Hessian approximation which is useful for predicting step-sizes.

Definition at line 26 of file bfgs_update.hpp.


The documentation for this class was generated from the following file:

     [ Stan Home Page ] © 2011–2016, Stan Development Team.