mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat > Class Template Reference
Collaboration diagram for mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >:

Documentation

template<typename Tile, typename Policy, typename Gradient, typename Mat>
class mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >

Implements non-linear conjugate gradient for a given gradient expression.

Public Member Functions

virtual ~NonLinearConjGradient ()=default
 
 NonLinearConjGradient (const KeyVal &kv)
 
double inner_product (Mat const &bra_mat, Mat const &ket_mat) const
 
std::vector< Mat > rotate_W (const int nkpts, double mu, std::vector< Mat > const &Hn, std::vector< Mat > const &Wn) const
 
int find_max_idx (Eigen::VectorXd const &list, double const max_value) const
 
Eigen::VectorXd fit_data (Eigen::VectorXd const &mu_list, Eigen::VectorXd const &P_list) const
 
void evaluate_P_on_range (Eigen::VectorXd &mu_list, Eigen::VectorXd &P_list, const double lower_bound, const double upper_bound, const int p, const int nkpts, Eigen::Vector3i const &nk, std::vector< Mat > const &Hn, std::vector< Mat > const &Wn, std::vector< std::vector< Mat >> const &Q_list) const
 
double get_upper_bound (double const &lower_bound, double const &upper_bound, Eigen::VectorXd &mu_list, Eigen::VectorXd &P_list, const int p, const double P_initial, Eigen::Vector3i const &nk, std::vector< Mat > const &Hn, std::vector< Mat > const &Wn, std::vector< std::vector< Mat >> const &Q_list, const int nkpts) const
 
double solve_root (const double &lower_bound, double &upper_bound, const int p, const double P_initial, const int nkpts, Eigen::Vector3i const &nk, std::vector< Mat > const &Hn, std::vector< Mat > const &Wn, std::vector< std::vector< Mat >> const &Q_list) const
 
std::vector< Mat > update_W_fit_functional (const int nkpts, const int p, double const &lower_bound, double &upper_bound, double &P_initial, Eigen::Vector3i const &nk, std::vector< Mat > const &Hn, std::vector< Mat > const &Wn, std::vector< std::vector< Mat >> const &Q_list) const
 
std::vector< Mat > compute (std::vector< Mat > const &C, std::vector< Mat > &C_rotated, Eigen::Vector3i const &nk, int64_t nocc, size_t ncols_of_C_to_skip=0) const final
 
- Public Member Functions inherited from mpqc::lcao::zOrbitalLocalizer< Tile, Policy, Mat >
virtual ~zOrbitalLocalizer ()=default
 

Protected Member Functions

void init (std::shared_ptr< detail::zOrbitalLocalizationGradient< Mat >> gradient)
 
 NonLinearConjGradient (int max_iter, double target_precision)
 

Protected Attributes

std::shared_ptr< WavefunctionWorldwfn_world_
 
std::shared_ptr< lcao::pbc::gaussian::AOFactory< TileD, Policy > > ao_factory_
 
double convergence_threshold_
 
size_t max_iter_
 
std::string update_factor_
 
std::string initial_W_
 
std::shared_ptr< detail::zOrbitalLocalizationGradient< Mat > > gradient_
 

Constructor & Destructor Documentation

◆ ~NonLinearConjGradient()

template<typename Tile , typename Policy , typename Gradient , typename Mat >
virtual mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::~NonLinearConjGradient ( )
virtualdefault

◆ NonLinearConjGradient() [1/2]

template<typename Tile , typename Policy , typename Gradient , typename Mat >
mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::NonLinearConjGradient ( const KeyVal kv)
inlineexplicit

KeyVal constructor for NonLinearConjGradient

Parameters
kvthe KeyVal object; it will be queried for the following keywords:
Keyword Type Default Description
wfn_world class none the WavefunctionWorld
convergence_threshold double 1e-8 the non-linear conjugate gradient solver is converged when the norm of the gradient falls below this value
max_iter int 50 the maximum number of iterations
update_factor string sdsa the means of forming gamma

◆ NonLinearConjGradient() [2/2]

template<typename Tile , typename Policy , typename Gradient , typename Mat >
mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::NonLinearConjGradient ( int  max_iter,
double  target_precision 
)
inlineprotected

Member Function Documentation

◆ compute()

template<typename Tile , typename Policy , typename Gradient , typename Mat >
std::vector<Mat> mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::compute ( std::vector< Mat > const &  C,
std::vector< Mat > &  C_rotated,
Eigen::Vector3i const &  nk,
int64_t  nocc,
size_t  ncols_of_C_to_skip = 0 
) const
inlinefinalvirtual
Parameters
Cinput LCAOs
[in]ncols_of_C_to_skipthe number of columns of C to keep non-localized, presumably because they are already localized
Returns
transformation matrix U that converts C to localized LCAOs, i.e.
Cao("mu,k") * U("k,i")
computes the AO coefficients of localized MOs from the AO coefficients of input MOs";

Step 1: Get the maximum of all omega_k

Step 2: Determine the order q of the cost function

Step 3: Determine the value of Tmu = 2 pi / (q * max_omega)

Step 1: Get the maximum of all omega_k

Step 2: Determine the order q of the cost function

Step 3: Determine the value of Tmu = 2 pi / (q * max_omega)

Implements mpqc::lcao::zOrbitalLocalizer< Tile, Policy, Mat >.

◆ evaluate_P_on_range()

template<typename Tile , typename Policy , typename Gradient , typename Mat >
void mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::evaluate_P_on_range ( Eigen::VectorXd &  mu_list,
Eigen::VectorXd &  P_list,
const double  lower_bound,
const double  upper_bound,
const int  p,
const int  nkpts,
Eigen::Vector3i const &  nk,
std::vector< Mat > const &  Hn,
std::vector< Mat > const &  Wn,
std::vector< std::vector< Mat >> const &  Q_list 
) const
inline

Step 5: Evaluate R(mu) at each mu_i

Step 6: For each R(mu), evaluate the PM functional

◆ find_max_idx()

template<typename Tile , typename Policy , typename Gradient , typename Mat >
int mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::find_max_idx ( Eigen::VectorXd const &  list,
double const  max_value 
) const
inline

Given an Eigen::VectorXd of items and the max item in that vector, find the index of the max item

◆ fit_data()

template<typename Tile , typename Policy , typename Gradient , typename Mat >
Eigen::VectorXd mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::fit_data ( Eigen::VectorXd const &  mu_list,
Eigen::VectorXd const &  P_list 
) const
inline

◆ get_upper_bound()

template<typename Tile , typename Policy , typename Gradient , typename Mat >
double mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::get_upper_bound ( double const &  lower_bound,
double const &  upper_bound,
Eigen::VectorXd &  mu_list,
Eigen::VectorXd &  P_list,
const int  p,
const double  P_initial,
Eigen::Vector3i const &  nk,
std::vector< Mat > const &  Hn,
std::vector< Mat > const &  Wn,
std::vector< std::vector< Mat >> const &  Q_list,
const int  nkpts 
) const
inline

◆ init()

template<typename Tile , typename Policy , typename Gradient , typename Mat >
void mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::init ( std::shared_ptr< detail::zOrbitalLocalizationGradient< Mat >>  gradient)
inlineprotected

◆ inner_product()

template<typename Tile , typename Policy , typename Gradient , typename Mat >
double mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::inner_product ( Mat const &  bra_mat,
Mat const &  ket_mat 
) const
inline

◆ rotate_W()

template<typename Tile , typename Policy , typename Gradient , typename Mat >
std::vector<Mat> mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::rotate_W ( const int  nkpts,
double  mu,
std::vector< Mat > const &  Hn,
std::vector< Mat > const &  Wn 
) const
inline

I believe that this is the correct one to use Given an appropriate double mu and vector of matrices Hn, generate the rotated W matrices Rotate each W with the same mu but a different H

◆ solve_root()

template<typename Tile , typename Policy , typename Gradient , typename Mat >
double mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::solve_root ( const double &  lower_bound,
double &  upper_bound,
const int  p,
const double  P_initial,
const int  nkpts,
Eigen::Vector3i const &  nk,
std::vector< Mat > const &  Hn,
std::vector< Mat > const &  Wn,
std::vector< std::vector< Mat >> const &  Q_list 
) const
inline

◆ update_W_fit_functional()

template<typename Tile , typename Policy , typename Gradient , typename Mat >
std::vector<Mat> mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::update_W_fit_functional ( const int  nkpts,
const int  p,
double const &  lower_bound,
double &  upper_bound,
double &  P_initial,
Eigen::Vector3i const &  nk,
std::vector< Mat > const &  Hn,
std::vector< Mat > const &  Wn,
std::vector< std::vector< Mat >> const &  Q_list 
) const
inline

Compute the optimal step size, \mu_{opt} and update W by fitting the P(mu) values and then finding the roots of the first derivative of this function

Parameters
[in]nkptsthe number of k points
[in]pthe order of the polynomial
[in]Hnthe current iteration's search descent direction: one matrix for each k point
[in]Wnthe current iteration's estimate of the orbital transformation matrix: one matrix for each k point
Returns
mu_n a vector containing the current optimal step size for each k point

Member Data Documentation

◆ ao_factory_

template<typename Tile , typename Policy , typename Gradient , typename Mat >
std::shared_ptr<lcao::pbc::gaussian::AOFactory<TileD, Policy> > mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::ao_factory_
protected

◆ convergence_threshold_

template<typename Tile , typename Policy , typename Gradient , typename Mat >
double mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::convergence_threshold_
protected

◆ gradient_

template<typename Tile , typename Policy , typename Gradient , typename Mat >
std::shared_ptr<detail::zOrbitalLocalizationGradient<Mat> > mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::gradient_
protected

◆ initial_W_

template<typename Tile , typename Policy , typename Gradient , typename Mat >
std::string mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::initial_W_
protected

◆ max_iter_

template<typename Tile , typename Policy , typename Gradient , typename Mat >
size_t mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::max_iter_
protected

◆ update_factor_

template<typename Tile , typename Policy , typename Gradient , typename Mat >
std::string mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::update_factor_
protected

◆ wfn_world_

template<typename Tile , typename Policy , typename Gradient , typename Mat >
std::shared_ptr<WavefunctionWorld> mpqc::lcao::pbc::NonLinearConjGradient< Tile, Policy, Gradient, Mat >::wfn_world_
protected

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