QMCPACK
DescentEngine Class Reference
+ Collaboration diagram for DescentEngine:

Public Member Functions

 DescentEngine (Communicate *comm, const xmlNodePtr cur)
 Constructor for engine. More...
 
bool processXML (const xmlNodePtr cur)
 process xml node More...
 
void prepareStorage (const int num_replicas, const int num_optimizables)
 Prepare for taking samples. More...
 
void takeSample (const int replica_id, const std::vector< FullPrecValueType > &der_rat_samp, const std::vector< FullPrecValueType > &le_der_samp, const std::vector< FullPrecValueType > &ls_der_samp, ValueType vgs_samp, ValueType weight_samp)
 Function that Take Sample Data from the Host Code. More...
 
void sample_finish ()
 Function that reduces all vector information from all processors to the root processor. More...
 
void mpi_unbiased_ratio_of_means (int numSamples, std::vector< ValueType > &weights, std::vector< ValueType > &numerSamples, std::vector< ValueType > &denomSamples, ValueType &mean, ValueType &variance, ValueType &stdErr)
 Function for computing ratios of the form <f>/<g> as well as the associated variance and standard error. More...
 
const std::vector< ValueType > & getAveragedDerivatives () const
 Returns the derivatives of the cost function we are minimizing. More...
 
void updateParameters ()
 helper method for updating parameter values with descent More...
 
ValueType setStepSize (int i)
 helper method for seting step sizes for different parameter types in descent optimization More...
 
void storeDerivRecord ()
 stores derivatives so they can be used in accelerated descent algorithm on later iterations More...
 
void setupUpdate (const optimize::VariableSet &my_vars)
 helper method for transferring information on parameter names and types to the engine More...
 
void storeVectors (std::vector< ValueType > &current_params)
 Store a vector of parameter differences to be used by the BLM in a hybrid optimization. More...
 
void computeFinalizationUncertainties (std::vector< ValueType > &weights, std::vector< ValueType > &numerSamples, std::vector< ValueType > &denomSamples)
 Compute uncertainties for energy/target function and variance over a history of samples from a set of iterations. More...
 
int retrieveStoreFrequency () const
 Returns number of times a parameter difference vector will be stored in the optimization. More...
 
const std::vector< std::vector< ValueType > > & retrieveHybridBLM_Input () const
 Returns the set of stored parameter difference vectors that will be given to the BLM. More...
 
const std::vector< ValueType > & retrieveNewParams () const
 Returns the current set of parameter values. More...
 
int getDescentNum () const
 Returns number of optimization steps that have been taken with descent. More...
 
ValueType getOmega () const
 Returns current value of omega. More...
 
ValueType getEnergy () const
 Returns value of average energy. More...
 
ValueType getVariance () const
 Returns variance of the energy. More...
 
ValueType getSD () const
 Returns standard deviation of energy. More...
 
bool targetingExcited () const
 Returns whether an excited state is being targeted. More...
 
int getFinalDescentNum () const
 Returns the descent iteration number when on a finalizing descent section. More...
 
void resetStorageCount ()
 Resets the number of vectors stored to 0 for next hybrid method macro-iteration. More...
 
void setDerivs (std::vector< ValueType > &test_derivs)
 Function for setting averaged derivatives, currently only used as part of a unit test of the engine's parameter update. More...
 
void setParamVal (int index, ValueType value)
 Function for setting parameter value, used to keep descent parameter values up to date with changes that occur on BLM steps of hybrid method. More...
 

Private Types

using FullPrecValueType = qmcplusplus::QMCTraits::FullPrecValueType
 
using ValueType = qmcplusplus::QMCTraits::ValueType
 
using RealType = qmcplusplus::QMCTraits::RealType
 
using FullPrecRealType = qmcplusplus::QMCTraits::FullPrecRealType
 

Private Attributes

std::vector< FullPrecValueTypeavg_le_der_samp_
 Vector for local energy parameter derivatives. More...
 
std::vector< std::vector< FullPrecValueType > > replica_le_der_samp_
 Vector for local energy parameter derivatives on one thread. More...
 
std::vector< FullPrecValueTypeavg_der_rat_samp_
 Vector for WF parameter derivatives. More...
 
std::vector< std::vector< FullPrecValueType > > replica_der_rat_samp_
 Vector for WF parameter derivatives on one thread. More...
 
std::vector< FullPrecValueTypeavg_numer_der_samp_
 Vector for target function numerator parameter derivatives. More...
 
std::vector< std::vector< FullPrecValueType > > replica_numer_der_samp_
 Vector for target function numerator parameter derivatives on one thread. More...
 
std::vector< FullPrecValueTypeavg_denom_der_samp_
 Vector for target function denominator parameter derivatives. More...
 
std::vector< std::vector< FullPrecValueType > > replica_denom_der_samp_
 Vector for target function denominator parameter derivatives on one thread. More...
 
ValueType w_sum_
 Total sum of weights. More...
 
ValueType e_avg_
 Average energy on a descent step. More...
 
ValueType e_var_
 Variance of the energy. More...
 
ValueType e_sd_
 Standard deviation of the energy. More...
 
ValueType e_err_
 Standard error of the energy. More...
 
ValueType numer_avg_
 Average target function numerator on a descent step. More...
 
ValueType numer_var_
 Variance of the target function numerator. More...
 
ValueType numer_err_
 Standard error of the target function numerator. More...
 
ValueType denom_avg_
 Average target function denominator on a descent step. More...
 
ValueType denom_var_
 Variance of the target function denominator. More...
 
ValueType denom_err_
 Standard error of the target function denominator. More...
 
ValueType target_avg_
 Average target function value on a descent step. More...
 
ValueType target_var_
 Variance of the target function. More...
 
ValueType target_err_
 Standard error of the target function. More...
 
std::vector< ValueTypevg_history_
 history of sampled |value/guiding|^2 ratios for one iteration More...
 
std::vector< std::vector< ValueType > > replica_vg_history_
 
std::vector< ValueTypefinal_vg_history_
 history of sampled |value/guiding|^2 ratios during the descent finalization phase More...
 
std::vector< std::vector< ValueType > > replica_final_vg_history_
 
std::vector< ValueTypew_history_
 history of sampled configuration weights for one iteration More...
 
std::vector< std::vector< ValueType > > replica_w_history_
 
std::vector< ValueTypefinal_w_history_
 history of sampled configuration weights during descent finalization phase More...
 
std::vector< std::vector< ValueType > > replica_final_w_history_
 
std::vector< ValueTypelev_history_
 a history of sampled local energies times the |value/guiding|^2 raitos for one iteration More...
 
std::vector< std::vector< ValueType > > replica_lev_history_
 
std::vector< ValueTypefinal_lev_history_
 history of sampled local energies times the |value/guiding|^2 raitos during the descent finalization phase More...
 
std::vector< std::vector< ValueType > > replica_final_lev_history_
 
std::vector< ValueTypefinal_le_avg_history_
 a vector to store the averages of the energy during the descent finalization phase More...
 
std::vector< ValueTypefinal_var_avg_history_
 a vector to store the variances of the energy during the descent finalization phase More...
 
std::vector< ValueTypetnv_history_
 a history of target function numerator times the |value/guiding|^2 ratios for one iteration More...
 
std::vector< std::vector< ValueType > > replica_tnv_history_
 
std::vector< ValueTypefinal_tnv_history_
 a history of target function numerator times the |value/guiding|^2 ratios during the descent finalization phase More...
 
std::vector< std::vector< ValueType > > replica_final_tnv_history_
 
std::vector< ValueTypetdv_history_
 a history of target function denominator times the |value/guiding|^2 ratios for one iteration More...
 
std::vector< std::vector< ValueType > > replica_tdv_history_
 
std::vector< ValueTypefinal_tdv_history_
 a history of target function denomerator times the |value/guiding|^2 ratios during the descent finalization phase More...
 
std::vector< std::vector< ValueType > > replica_final_tdv_history_
 
std::vector< ValueTypefinal_tar_avg_history_
 a vector to store the averages of the target function during the descent finalization phase More...
 
std::vector< ValueTypefinal_tar_var_history_
 a vector to store the variances of the target function during the descent finalization phase More...
 
std::vector< ValueTypelderivs_
 Vector that stores the final averaged derivatives of the cost function. More...
 
Communicatemy_comm_
 Communicator handles MPI reduction. More...
 
bool engine_target_excited_
 Whether to target excited state. More...
 
int num_params_
 Number of optimizable parameters. More...
 
std::vector< ValueTypeparams_copy_
 Vector for storing parameter values from previous optimization step. More...
 
std::vector< ValueTypecurrent_params_
 Vector for storing parameter values for current optimization step. More...
 
std::vector< std::vector< ValueType > > deriv_records_
 Vector for storing Lagrangian derivatives from previous optimization steps. More...
 
std::vector< ValueTypedenom_records_
 Vector for storing step size denominator values from previous optimization step. More...
 
std::vector< ValueTypenumer_records_
 Vector for storing step size numerator values from previous optimization step. More...
 
ValueType lambda_ = 0.0
 Parameter for accelerated descent recursion relation. More...
 
std::vector< ValueTypetaus_
 Vector for storing step sizes from previous optimization step. More...
 
std::vector< ValueTypederivs_squared_
 Vector for storing running average of squares of the derivatives. More...
 
int descent_num_
 Integer for keeping track of only number of descent steps taken. More...
 
std::string flavor_
 What variety of gradient descent will be used. More...
 
ValueType tjf_2body_eta_
 Step sizes for different types of parameters. More...
 
ValueType tjf_1body_eta_
 
ValueType f_eta_
 
ValueType gauss_eta_
 
ValueType ci_eta_
 
ValueType orb_eta_
 
bool ramp_eta_
 Whether to gradually ramp up step sizes in descent. More...
 
int ramp_num_
 Number of steps over which to ramp up step size. More...
 
int store_num_
 Number of parameter difference vectors stored when descent is used in a hybrid optimization. More...
 
int store_count_
 Counter of how many vectors have been stored so far. More...
 
std::vector< std::string > engine_param_names_
 Vectors of parameter names and types, used in the assignment of step sizes. More...
 
std::vector< int > engine_param_types_
 
std::vector< ValueTypeparams_for_diff_
 Vector for storing parameter values for calculating differences to be given to hybrid method. More...
 
std::vector< std::vector< ValueType > > hybrid_blm_input_
 Vector for storing the input vectors to the BLM steps of hybrid method. More...
 
ValueType omega_
 Value of omega in excited state functional. More...
 
int collection_step_
 Iteration to start collecting samples for final average and error blocking analysis. More...
 
int compute_step_
 Iteration to start computing averages and errors from the stored values during the finalization phase. More...
 
bool collect_count_ = false
 Whether to start collecting samples for the histories in the finalization phase. More...
 
int final_descent_num_ = 0
 Counter for the number of descent steps taken in the finalization phase. More...
 
std::string print_deriv_
 Whether to print out derivative terms for each parameter. More...
 

Detailed Description

Definition at line 27 of file DescentEngine.h.

Member Typedef Documentation

◆ FullPrecRealType

Definition at line 32 of file DescentEngine.h.

◆ FullPrecValueType

◆ RealType

Definition at line 31 of file DescentEngine.h.

◆ ValueType

Definition at line 30 of file DescentEngine.h.

Constructor & Destructor Documentation

◆ DescentEngine()

DescentEngine ( Communicate comm,
const xmlNodePtr  cur 
)

Constructor for engine.

Definition at line 29 of file DescentEngine.cpp.

References DescentEngine::descent_num_, DescentEngine::processXML(), and DescentEngine::store_count_.

30  : my_comm_(comm),
32  num_params_(0),
33  flavor_("RMSprop"),
34  tjf_2body_eta_(.01),
35  tjf_1body_eta_(.01),
36  f_eta_(.001),
37  gauss_eta_(.001),
38  ci_eta_(.01),
39  orb_eta_(.001),
40  ramp_eta_(false),
41  ramp_num_(30),
42  store_num_(5),
43  omega_(0.0),
44  collection_step_(-1),
45  compute_step_(-1),
46  print_deriv_("no")
47 {
48  descent_num_ = 0;
49  store_count_ = 0;
50  processXML(cur);
51 }
int num_params_
Number of optimizable parameters.
Communicate * my_comm_
Communicator handles MPI reduction.
int store_count_
Counter of how many vectors have been stored so far.
int ramp_num_
Number of steps over which to ramp up step size.
ValueType tjf_2body_eta_
Step sizes for different types of parameters.
int store_num_
Number of parameter difference vectors stored when descent is used in a hybrid optimization.
int collection_step_
Iteration to start collecting samples for final average and error blocking analysis.
bool engine_target_excited_
Whether to target excited state.
int compute_step_
Iteration to start computing averages and errors from the stored values during the finalization phase...
bool processXML(const xmlNodePtr cur)
process xml node
ValueType omega_
Value of omega in excited state functional.
bool ramp_eta_
Whether to gradually ramp up step sizes in descent.
int descent_num_
Integer for keeping track of only number of descent steps taken.
std::string flavor_
What variety of gradient descent will be used.
std::string print_deriv_
Whether to print out derivative terms for each parameter.

Member Function Documentation

◆ computeFinalizationUncertainties()

void computeFinalizationUncertainties ( std::vector< ValueType > &  weights,
std::vector< ValueType > &  numerSamples,
std::vector< ValueType > &  denomSamples 
)

Compute uncertainties for energy/target function and variance over a history of samples from a set of iterations.

Definition at line 966 of file DescentEngine.cpp.

References qmcplusplus::abs(), qmcplusplus::app_log(), copy(), DescentEngine::mpi_unbiased_ratio_of_means(), qmcplusplus::n, and qmcplusplus::sqrt().

Referenced by DescentEngine::updateParameters().

969 {
970  // Make copies of input vectors to do recursive blocking for error estimates
971  std::vector<ValueType> wtv(weights);
972  std::vector<ValueType> nmv(numerSamples);
973  std::vector<ValueType> dnv(denomSamples);
974 
975  // Reproduce LM engine's crude way of estimating uncertainty in the variance
976  int n = nmv.size();
977  int blocks = 100;
978  int section_len = n / 100;
979 
980  std::vector<ValueType> bbvars(blocks);
981  ValueType bbv = 0.0;
982  ValueType bbv2 = 0.0;
983 
984  for (int i = 0; i < blocks; i++)
985  {
986  std::vector<ValueType> sub_wtv(section_len);
987  std::vector<ValueType> sub_nmv(section_len);
988  std::vector<ValueType> sub_dnv(section_len);
989 
990  ValueType temp_e, temp_v;
991  std::copy(wtv.begin() + i * section_len, wtv.begin() + (i + 1) * section_len, sub_wtv.begin());
992  std::copy(nmv.begin() + i * section_len, nmv.begin() + (i + 1) * section_len, sub_nmv.begin());
993  std::copy(dnv.begin() + i * section_len, dnv.begin() + (i + 1) * section_len, sub_dnv.begin());
994 
995  ValueType blockVar;
996  ValueType blockMean;
997  ValueType blockErr;
998 
999  this->mpi_unbiased_ratio_of_means(section_len, sub_wtv, sub_nmv, sub_dnv, blockMean, blockVar, blockErr);
1000 
1001  bbvars[i] = blockVar;
1002 
1003  bbv += bbvars.at(i) / static_cast<ValueType>(blocks);
1004  bbv2 += bbvars.at(i) * bbvars.at(i) / static_cast<ValueType>(blocks);
1005  sub_wtv.clear();
1006  sub_nmv.clear();
1007  sub_dnv.clear();
1008  }
1009 
1010  const ValueType bbvov = bbv2 - bbv * bbv;
1011  ValueType var_uncertainty = std::sqrt(std::abs(bbvov) / blocks);
1012  // Depending on when this function is called, this will be the uncertainty in
1013  // the variance
1014  // of either the energy or the target function.
1015  // Which one should be clear from the preceding print statements in the
1016  // output file.
1017  app_log() << "Uncertainty in variance of averaged quantity: " << var_uncertainty << std::endl;
1018 
1019  // To compute the standard error on the energy/target function itself,
1020  // blocking is performed on the whole set of samples in the history of many
1021  // iterations.
1022  while (n > 128)
1023  {
1024  ValueType currentVar;
1025  ValueType currentMean;
1026  ValueType currentErr;
1027  this->mpi_unbiased_ratio_of_means(n, wtv, nmv, dnv, currentMean, currentVar, currentErr);
1028 
1029  app_log() << "Blocking analysis error for per process length " << n << " is: " << currentErr << std::endl;
1030 
1031  std::vector<ValueType> tmp1, tmp2, tmp3;
1032 
1033  for (int i = 0; i < n; i += 2)
1034  {
1035  ValueType avgW = (wtv[i] + wtv[i + 1]) / static_cast<ValueType>(2.0);
1036  ValueType avgNumer = (nmv[i] + nmv[i + 1]) / static_cast<ValueType>(2.0);
1037  ValueType avgDenom = (dnv[i] + dnv[i + 1]) / static_cast<ValueType>(2.0);
1038 
1039  tmp1.push_back(avgW);
1040  tmp2.push_back(avgNumer);
1041  tmp3.push_back(avgDenom);
1042  }
1043  wtv = tmp1;
1044  nmv = tmp2;
1045  dnv = tmp3;
1046  n = nmv.size();
1047  }
1048 }
MakeReturn< UnaryNode< FnFabs, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t abs(const Vector< T1, C1 > &l)
std::ostream & app_log()
Definition: OutputManager.h:65
void copy(const Array< T1, 3 > &src, Array< T2, 3 > &dest)
Definition: Blitz.h:639
MakeReturn< UnaryNode< FnSqrt, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t sqrt(const Vector< T1, C1 > &l)
qmcplusplus::QMCTraits::ValueType ValueType
Definition: DescentEngine.h:30
void mpi_unbiased_ratio_of_means(int numSamples, std::vector< ValueType > &weights, std::vector< ValueType > &numerSamples, std::vector< ValueType > &denomSamples, ValueType &mean, ValueType &variance, ValueType &stdErr)
Function for computing ratios of the form <f>/<g> as well as the associated variance and standard err...

◆ getAveragedDerivatives()

const std::vector<ValueType>& getAveragedDerivatives ( ) const
inline

Returns the derivatives of the cost function we are minimizing.

Definition at line 293 of file DescentEngine.h.

References DescentEngine::lderivs_.

293 { return lderivs_; }
std::vector< ValueType > lderivs_
Vector that stores the final averaged derivatives of the cost function.

◆ getDescentNum()

int getDescentNum ( ) const
inline

Returns number of optimization steps that have been taken with descent.

Definition at line 332 of file DescentEngine.h.

References DescentEngine::descent_num_.

332 { return descent_num_; }
int descent_num_
Integer for keeping track of only number of descent steps taken.

◆ getEnergy()

ValueType getEnergy ( ) const
inline

Returns value of average energy.

Definition at line 338 of file DescentEngine.h.

References DescentEngine::e_avg_.

338 { return e_avg_; }
ValueType e_avg_
Average energy on a descent step.
Definition: DescentEngine.h:59

◆ getFinalDescentNum()

int getFinalDescentNum ( ) const
inline

Returns the descent iteration number when on a finalizing descent section.

Definition at line 350 of file DescentEngine.h.

References DescentEngine::final_descent_num_.

350 { return final_descent_num_; }
int final_descent_num_
Counter for the number of descent steps taken in the finalization phase.

◆ getOmega()

ValueType getOmega ( ) const
inline

Returns current value of omega.

Definition at line 335 of file DescentEngine.h.

References DescentEngine::omega_.

335 { return omega_; }
ValueType omega_
Value of omega in excited state functional.

◆ getSD()

ValueType getSD ( ) const
inline

Returns standard deviation of energy.

Definition at line 344 of file DescentEngine.h.

References DescentEngine::e_sd_.

344 { return e_sd_; }
ValueType e_sd_
Standard deviation of the energy.
Definition: DescentEngine.h:63

◆ getVariance()

ValueType getVariance ( ) const
inline

Returns variance of the energy.

Definition at line 341 of file DescentEngine.h.

References DescentEngine::e_var_.

341 { return e_var_; }
ValueType e_var_
Variance of the energy.
Definition: DescentEngine.h:61

◆ mpi_unbiased_ratio_of_means()

void mpi_unbiased_ratio_of_means ( int  numSamples,
std::vector< ValueType > &  weights,
std::vector< ValueType > &  numerSamples,
std::vector< ValueType > &  denomSamples,
ValueType mean,
ValueType variance,
ValueType stdErr 
)

Function for computing ratios of the form <f>/<g> as well as the associated variance and standard error.

Definition at line 481 of file DescentEngine.cpp.

References Communicate::allreduce(), mean(), qmcplusplus::Units::mass::mp, DescentEngine::my_comm_, qmcplusplus::n, qmcplusplus::Units::time::ns, qmcplusplus::sqrt(), and DescentEngine::w_sum_.

Referenced by DescentEngine::computeFinalizationUncertainties(), and DescentEngine::sample_finish().

488 {
489  std::vector<ValueType> y(7);
490  y[0] = 0.0; // normalization constant
491  y[1] = 0.0; // mean of numerator
492  y[2] = 0.0; // mean of denominator
493  y[3] = 0.0; // mean of the square of the numerator terms
494  y[4] = 0.0; // mean of the square of the denominator terms
495  y[5] = 0.0; // mean of the product of numerator times denominator
496  y[6] = ValueType(numSamples); // number of samples
497 
498  for (int i = 0; i < numSamples; i++)
499  {
500  ValueType n = numerSamples[i];
501  ValueType d = denomSamples[i];
502  ValueType weight = weights[i];
503 
504  y[0] += weight;
505  y[1] += weight * n;
506  y[2] += d;
507  y[3] += weight * n * n;
508  y[4] += weight * d * d;
509  y[5] += weight * n * d;
510  }
511 
512  my_comm_->allreduce(y);
513 
514  ValueType mf = y[1] / y[0]; // mean of numerator
515  ValueType mg = y[2] / y[0]; // mean of denominator
516  ValueType sf = y[3] / y[0]; // mean of the square of the numerator terms
517  ValueType sg = y[4] / y[0]; // mean of the square of the denominator terms
518  ValueType mp = y[5] / y[0]; // mean of the product of numerator times denominator
519  ValueType ns = y[6]; // number of samples
520 
521  ValueType vf = (sf - mf * mf) * ns / (ns - static_cast<ValueType>(1.0));
522  ValueType vg = (sg - mg * mg) * ns / (ns - static_cast<ValueType>(1.0));
523  ValueType cv = (mp - mf * mg) * ns / (ns - static_cast<ValueType>(1.0));
524 
525  w_sum_ = y[0];
526  mean = (mf / mg) / (static_cast<ValueType>(1.0) + (vg / mg / mg - cv / mf / mg) / ns);
527  variance = (mf * mf / mg / mg) * (vf / mf / mf + vg / mg / mg - static_cast<ValueType>(2.0) * cv / mf / mg);
528  stdErr = std::sqrt(variance / ns);
529 }
Communicate * my_comm_
Communicator handles MPI reduction.
void allreduce(T &)
ValueType w_sum_
Total sum of weights.
Definition: DescentEngine.h:56
MakeReturn< UnaryNode< FnSqrt, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t sqrt(const Vector< T1, C1 > &l)
qmcplusplus::QMCTraits::ValueType ValueType
Definition: DescentEngine.h:30
ACC::value_type mean(const ACC &ac)
Definition: accumulators.h:147

◆ prepareStorage()

void prepareStorage ( const int  num_replicas,
const int  num_optimizables 
)

Prepare for taking samples.

Definition at line 103 of file DescentEngine.cpp.

References DescentEngine::avg_denom_der_samp_, DescentEngine::avg_der_rat_samp_, DescentEngine::avg_le_der_samp_, DescentEngine::avg_numer_der_samp_, DescentEngine::collect_count_, DescentEngine::collection_step_, DescentEngine::denom_avg_, DescentEngine::denom_err_, DescentEngine::denom_var_, DescentEngine::e_avg_, DescentEngine::e_err_, DescentEngine::e_sd_, DescentEngine::e_var_, DescentEngine::engine_target_excited_, DescentEngine::final_descent_num_, DescentEngine::lderivs_, DescentEngine::num_params_, DescentEngine::numer_avg_, DescentEngine::numer_err_, DescentEngine::numer_var_, DescentEngine::replica_denom_der_samp_, DescentEngine::replica_der_rat_samp_, DescentEngine::replica_final_lev_history_, DescentEngine::replica_final_tdv_history_, DescentEngine::replica_final_tnv_history_, DescentEngine::replica_final_vg_history_, DescentEngine::replica_final_w_history_, DescentEngine::replica_le_der_samp_, DescentEngine::replica_lev_history_, DescentEngine::replica_numer_der_samp_, DescentEngine::replica_tdv_history_, DescentEngine::replica_tnv_history_, DescentEngine::replica_vg_history_, DescentEngine::replica_w_history_, DescentEngine::target_avg_, DescentEngine::target_err_, DescentEngine::target_var_, and DescentEngine::w_sum_.

Referenced by DescentEngineHandle::prepareSampling().

104 {
105  lderivs_.resize(num_optimizables);
106 
107  //Resize the history vectors for the current iteration for the number of threads present
108  replica_vg_history_.resize(num_replicas);
109  replica_w_history_.resize(num_replicas);
110  replica_lev_history_.resize(num_replicas);
111 
113  {
114  replica_tnv_history_.resize(num_replicas);
115  replica_tdv_history_.resize(num_replicas);
116 
117  }
118 
119 
120  //Also resize the history vectors for descent finalization if necessary
122  {
123  replica_final_vg_history_.resize(num_replicas);
124  replica_final_w_history_.resize(num_replicas);
125  replica_final_lev_history_.resize(num_replicas);
126 
127 
128 
130  {
131  replica_final_tnv_history_.resize(num_replicas);
132  replica_final_tdv_history_.resize(num_replicas);
133 
134  }
135  }
136 
137 
138  // Ground state case
140  {
141  avg_le_der_samp_.resize(num_optimizables);
142  avg_der_rat_samp_.resize(num_optimizables);
143 
144  num_params_ = num_optimizables;
145 
146  std::fill(avg_le_der_samp_.begin(), avg_le_der_samp_.end(), 0.0);
147  std::fill(avg_der_rat_samp_.begin(), avg_der_rat_samp_.end(), 0.0);
148 
149  replica_le_der_samp_.resize(num_replicas);
150  replica_der_rat_samp_.resize(num_replicas);
151 
152  for (int i = 0; i < num_replicas; i++)
153  {
154  replica_le_der_samp_[i].resize(num_optimizables);
155  std::fill(replica_le_der_samp_[i].begin(), replica_le_der_samp_[i].end(), 0.0);
156 
157  replica_der_rat_samp_[i].resize(num_optimizables);
158  std::fill(replica_der_rat_samp_[i].begin(), replica_der_rat_samp_[i].end(), 0.0);
159  }
160  }
161  // Excited state case
162  else
163  {
164  replica_numer_der_samp_.resize(num_replicas);
165  replica_denom_der_samp_.resize(num_replicas);
166 
167  avg_numer_der_samp_.resize(num_optimizables);
168  avg_denom_der_samp_.resize(num_optimizables);
169 
170  std::fill(avg_numer_der_samp_.begin(), avg_numer_der_samp_.end(), 0.0);
171  std::fill(avg_denom_der_samp_.begin(), avg_denom_der_samp_.end(), 0.0);
172 
173  for (int i = 0; i < num_replicas; i++)
174  {
175  replica_numer_der_samp_[i].resize(num_optimizables);
176  replica_denom_der_samp_[i].resize(num_optimizables);
177 
178  std::fill(replica_numer_der_samp_[i].begin(), replica_numer_der_samp_[i].end(), 0.0);
179  std::fill(replica_denom_der_samp_[i].begin(), replica_denom_der_samp_[i].end(), 0.0);
180  }
181  }
182 
183  w_sum_ = 0;
184  e_avg_ = 0;
185  e_var_ = 0;
186  e_sd_ = 0;
187  e_err_ = 0;
188 
189  numer_avg_ = 0;
190  numer_var_ = 0;
191  numer_err_ = 0;
192 
193  denom_avg_ = 0;
194  denom_var_ = 0;
195  denom_err_ = 0;
196 
197  target_avg_ = 0;
198  target_var_ = 0;
199  target_err_ = 0;
200 }
int num_params_
Number of optimizable parameters.
std::vector< std::vector< ValueType > > replica_final_tnv_history_
ValueType numer_err_
Standard error of the target function numerator.
Definition: DescentEngine.h:72
int final_descent_num_
Counter for the number of descent steps taken in the finalization phase.
std::vector< std::vector< FullPrecValueType > > replica_denom_der_samp_
Vector for target function denominator parameter derivatives on one thread.
Definition: DescentEngine.h:53
ValueType target_avg_
Average target function value on a descent step.
Definition: DescentEngine.h:82
std::vector< std::vector< ValueType > > replica_final_w_history_
ValueType denom_avg_
Average target function denominator on a descent step.
Definition: DescentEngine.h:75
std::vector< std::vector< FullPrecValueType > > replica_numer_der_samp_
Vector for target function numerator parameter derivatives on one thread.
Definition: DescentEngine.h:48
int collection_step_
Iteration to start collecting samples for final average and error blocking analysis.
ValueType denom_err_
Standard error of the target function denominator.
Definition: DescentEngine.h:79
std::vector< FullPrecValueType > avg_le_der_samp_
Vector for local energy parameter derivatives.
Definition: DescentEngine.h:36
ValueType e_err_
Standard error of the energy.
Definition: DescentEngine.h:65
std::vector< FullPrecValueType > avg_denom_der_samp_
Vector for target function denominator parameter derivatives.
Definition: DescentEngine.h:51
std::vector< FullPrecValueType > avg_numer_der_samp_
Vector for target function numerator parameter derivatives.
Definition: DescentEngine.h:46
std::vector< FullPrecValueType > avg_der_rat_samp_
Vector for WF parameter derivatives.
Definition: DescentEngine.h:41
bool engine_target_excited_
Whether to target excited state.
std::vector< std::vector< ValueType > > replica_final_lev_history_
ValueType denom_var_
Variance of the target function denominator.
Definition: DescentEngine.h:77
std::vector< std::vector< ValueType > > replica_tnv_history_
std::vector< ValueType > lderivs_
Vector that stores the final averaged derivatives of the cost function.
std::vector< std::vector< ValueType > > replica_lev_history_
std::vector< std::vector< ValueType > > replica_tdv_history_
ValueType e_sd_
Standard deviation of the energy.
Definition: DescentEngine.h:63
ValueType target_var_
Variance of the target function.
Definition: DescentEngine.h:84
ValueType w_sum_
Total sum of weights.
Definition: DescentEngine.h:56
ValueType numer_avg_
Average target function numerator on a descent step.
Definition: DescentEngine.h:68
ValueType e_avg_
Average energy on a descent step.
Definition: DescentEngine.h:59
std::vector< std::vector< FullPrecValueType > > replica_le_der_samp_
Vector for local energy parameter derivatives on one thread.
Definition: DescentEngine.h:38
bool collect_count_
Whether to start collecting samples for the histories in the finalization phase.
std::vector< std::vector< FullPrecValueType > > replica_der_rat_samp_
Vector for WF parameter derivatives on one thread.
Definition: DescentEngine.h:43
std::vector< std::vector< ValueType > > replica_final_vg_history_
Definition: DescentEngine.h:95
std::vector< std::vector< ValueType > > replica_final_tdv_history_
ValueType e_var_
Variance of the energy.
Definition: DescentEngine.h:61
std::vector< std::vector< ValueType > > replica_w_history_
ValueType target_err_
Standard error of the target function.
Definition: DescentEngine.h:86
std::vector< std::vector< ValueType > > replica_vg_history_
Definition: DescentEngine.h:90
ValueType numer_var_
Variance of the target function numerator.
Definition: DescentEngine.h:70

◆ processXML()

bool processXML ( const xmlNodePtr  cur)

process xml node

Definition at line 53 of file DescentEngine.cpp.

References ParameterSet::add(), qmcplusplus::app_log(), DescentEngine::ci_eta_, DescentEngine::collect_count_, DescentEngine::collection_step_, DescentEngine::compute_step_, DescentEngine::engine_target_excited_, DescentEngine::f_eta_, DescentEngine::flavor_, DescentEngine::gauss_eta_, DescentEngine::omega_, DescentEngine::orb_eta_, DescentEngine::print_deriv_, ParameterSet::put(), DescentEngine::ramp_eta_, DescentEngine::ramp_num_, DescentEngine::store_num_, DescentEngine::tjf_1body_eta_, and DescentEngine::tjf_2body_eta_.

Referenced by DescentEngine::DescentEngine(), and qmcplusplus::TEST_CASE().

54 {
55  std::string excited("no");
56  std::string ramp_eta_str("no");
57 
58  ParameterSet m_param;
59  m_param.add(excited, "targetExcited");
60  m_param.add(omega_, "omega");
61  // Type of descent method being used
62  m_param.add(flavor_, "flavor");
63  // Step size inputs for different parameter types
64  m_param.add(tjf_2body_eta_, "TJF_2Body_eta");
65  m_param.add(tjf_1body_eta_, "TJF_1Body_eta");
66  m_param.add(f_eta_, "F_eta");
67  m_param.add(ci_eta_, "CI_eta");
68  m_param.add(gauss_eta_, "Gauss_eta");
69  m_param.add(orb_eta_, "Orb_eta");
70  // Whether to gradually ramp up step sizes and over how many steps
71  m_param.add(ramp_eta_str, "Ramp_eta");
72  m_param.add(ramp_num_, "Ramp_num");
73  // If using descent as part of hybrid method, how many parameter difference
74  // vectors to store
75  m_param.add(store_num_, "Stored_Vectors");
76  m_param.add(print_deriv_, "print_derivs");
77 
78  // When to start storing samples for a final average and when to start
79  // computing it
80  m_param.add(collection_step_, "collection_step");
81  m_param.add(compute_step_, "compute_step");
82 
83  //app_log() << "Omega from input file: " << omega_ << std::endl;
84  //app_log() << "Current collection step: " << collection_step_ << std::endl;
85  // Use -1 as a default value when you don't collect history. Would want to
86  // collect only during the descent finalization section
87  if (collection_step_ != -1)
88  {
89  app_log() << "On descent finalization, have collect_count as true" << std::endl;
90  collect_count_ = true;
91  }
92 
93  m_param.put(cur);
94 
95  engine_target_excited_ = (excited == "yes");
96 
97  ramp_eta_ = (ramp_eta_str == "yes");
98 
99  return true;
100 }
int ramp_num_
Number of steps over which to ramp up step size.
ValueType tjf_2body_eta_
Step sizes for different types of parameters.
int store_num_
Number of parameter difference vectors stored when descent is used in a hybrid optimization.
std::ostream & app_log()
Definition: OutputManager.h:65
int collection_step_
Iteration to start collecting samples for final average and error blocking analysis.
bool put(std::istream &is) override
read from std::istream
Definition: ParameterSet.h:42
bool engine_target_excited_
Whether to target excited state.
int compute_step_
Iteration to start computing averages and errors from the stored values during the finalization phase...
class to handle a set of parameters
Definition: ParameterSet.h:27
void add(PDT &aparam, const std::string &aname_in, std::vector< PDT > candidate_values={}, TagStatus status=TagStatus::OPTIONAL)
add a new parameter corresponding to an xmlNode <parameter>
ValueType omega_
Value of omega in excited state functional.
bool ramp_eta_
Whether to gradually ramp up step sizes in descent.
bool collect_count_
Whether to start collecting samples for the histories in the finalization phase.
std::string flavor_
What variety of gradient descent will be used.
std::string print_deriv_
Whether to print out derivative terms for each parameter.

◆ resetStorageCount()

void resetStorageCount ( )
inline

Resets the number of vectors stored to 0 for next hybrid method macro-iteration.

Definition at line 354 of file DescentEngine.h.

References DescentEngine::store_count_.

354 { store_count_ = 0; }
int store_count_
Counter of how many vectors have been stored so far.

◆ retrieveHybridBLM_Input()

const std::vector<std::vector<ValueType> >& retrieveHybridBLM_Input ( ) const
inline

Returns the set of stored parameter difference vectors that will be given to the BLM.

Definition at line 326 of file DescentEngine.h.

References DescentEngine::hybrid_blm_input_.

326 { return hybrid_blm_input_; }
std::vector< std::vector< ValueType > > hybrid_blm_input_
Vector for storing the input vectors to the BLM steps of hybrid method.

◆ retrieveNewParams()

const std::vector<ValueType>& retrieveNewParams ( ) const
inline

Returns the current set of parameter values.

Definition at line 329 of file DescentEngine.h.

References DescentEngine::current_params_.

329 { return current_params_; }
std::vector< ValueType > current_params_
Vector for storing parameter values for current optimization step.

◆ retrieveStoreFrequency()

int retrieveStoreFrequency ( ) const
inline

Returns number of times a parameter difference vector will be stored in the optimization.

Definition at line 322 of file DescentEngine.h.

References DescentEngine::store_num_.

322 { return store_num_; }
int store_num_
Number of parameter difference vectors stored when descent is used in a hybrid optimization.

◆ sample_finish()

void sample_finish ( )

Function that reduces all vector information from all processors to the root processor.

Definition at line 292 of file DescentEngine.cpp.

References Communicate::allreduce(), qmcplusplus::app_log(), DescentEngine::avg_denom_der_samp_, DescentEngine::avg_der_rat_samp_, DescentEngine::avg_le_der_samp_, DescentEngine::avg_numer_der_samp_, DescentEngine::collect_count_, DescentEngine::collection_step_, DescentEngine::denom_avg_, DescentEngine::denom_err_, DescentEngine::denom_var_, DescentEngine::e_avg_, DescentEngine::e_err_, DescentEngine::e_sd_, DescentEngine::e_var_, DescentEngine::engine_target_excited_, DescentEngine::final_descent_num_, DescentEngine::final_le_avg_history_, DescentEngine::final_lev_history_, DescentEngine::final_tar_avg_history_, DescentEngine::final_tar_var_history_, DescentEngine::final_tdv_history_, DescentEngine::final_tnv_history_, DescentEngine::final_var_avg_history_, DescentEngine::final_vg_history_, DescentEngine::final_w_history_, DescentEngine::lderivs_, DescentEngine::lev_history_, DescentEngine::mpi_unbiased_ratio_of_means(), DescentEngine::my_comm_, DescentEngine::numer_avg_, DescentEngine::numer_err_, DescentEngine::numer_var_, DescentEngine::print_deriv_, DescentEngine::replica_denom_der_samp_, DescentEngine::replica_der_rat_samp_, DescentEngine::replica_final_lev_history_, DescentEngine::replica_final_tdv_history_, DescentEngine::replica_final_tnv_history_, DescentEngine::replica_final_vg_history_, DescentEngine::replica_final_w_history_, DescentEngine::replica_le_der_samp_, DescentEngine::replica_lev_history_, DescentEngine::replica_numer_der_samp_, DescentEngine::replica_tdv_history_, DescentEngine::replica_tnv_history_, DescentEngine::replica_vg_history_, DescentEngine::replica_w_history_, qmcplusplus::sqrt(), DescentEngine::target_avg_, DescentEngine::target_err_, DescentEngine::target_var_, DescentEngine::tdv_history_, DescentEngine::tnv_history_, DescentEngine::vg_history_, DescentEngine::w_history_, and DescentEngine::w_sum_.

Referenced by DescentEngineHandle::finishSampling().

293 {
294 
295  //After a potentially multithreaded section, cocatenate the replica history vectors
296  int num_threads = replica_vg_history_.size();
297  for(int i =0; i< num_threads; i++)
298  {
299  vg_history_.insert(vg_history_.end(), replica_vg_history_[i].begin(),replica_vg_history_[i].end());
300  w_history_.insert(w_history_.end(), replica_w_history_[i].begin(),replica_w_history_[i].end());
301  lev_history_.insert(lev_history_.end(),replica_lev_history_[i].begin(),replica_lev_history_[i].end());
302 
303  //Clear the individual thread history vectors for the next iteration after their values have been collected
304  replica_vg_history_[i].clear();
305  replica_w_history_[i].clear();
306  replica_lev_history_[i].clear();
307 
308 
310  {
311 
312  tnv_history_.insert(tnv_history_.end(), replica_tnv_history_[i].begin(),replica_tnv_history_[i].end());
313  tdv_history_.insert(tdv_history_.end(), replica_tdv_history_[i].begin(),replica_tdv_history_[i].end());
314 
315  replica_tnv_history_[i].clear();
316  replica_tdv_history_[i].clear();
317  }
318  }
319 
320  //Do the same for finalization histories if necessary
322  {
323  for(int i =0; i< num_threads; i++)
324  {
328 
329  replica_final_vg_history_[i].clear();
330  replica_final_w_history_[i].clear();
331  replica_final_lev_history_[i].clear();
332 
333 
335  {
336 
339 
340  replica_final_tnv_history_[i].clear();
341  replica_final_tdv_history_[i].clear();
342  }
343  }
344  }
345 
346 
347  int numSamples = lev_history_.size();
348 
349  //Compute average energy and variance for this iteration
351 
352  app_log() << "Energy Average: " << std::setprecision(9) << e_avg_ << std::endl;
353  app_log() << "Energy Variance: " << std::setprecision(9) << e_var_ << std::endl;
354  app_log() << "Weight Total: " << std::setprecision(9) << w_sum_ << std::endl;
356  app_log() << "Energy Standard Deviation: " << std::setprecision(9) << e_sd_ << std::endl;
357  app_log() << "Energy Standard Error: " << std::setprecision(9) << e_err_ << std::endl;
358 
359  //Store average values during descent finalization
361  {
362  final_le_avg_history_.push_back(e_avg_);
363  final_var_avg_history_.push_back(e_var_);
364  }
365 
367  {
368  for (int i = 0; i < replica_le_der_samp_.size(); i++)
369  {
370  for (int j = 0; j < lderivs_.size(); j++)
371  {
372  avg_le_der_samp_[j] += replica_le_der_samp_[i].at(j);
374  }
375  }
376 
379  }
380  else
381  {
383  numer_err_);
384 
386  denom_err_);
387 
389  target_err_);
390 
391  app_log() << "Target Function Average: " << std::setprecision(9) << target_avg_ << std::endl;
392  app_log() << "Target Function Variance: " << std::setprecision(9) << target_var_ << std::endl;
393  app_log() << "Target Function Error: " << std::setprecision(9) << target_err_ << std::endl;
394 
396  {
399  }
400  for (int i = 0; i < replica_numer_der_samp_.size(); i++)
401  {
402  for (int j = 0; j < lderivs_.size(); j++)
403  {
406  }
407  }
408 
411  }
412 
413  int num_optimizables = lderivs_.size();
414 
415  // Vectors for parts of excited state functional derivatives
416  std::vector<ValueType> numer_term1(num_optimizables, 0.0);
417  std::vector<ValueType> numer_term2(num_optimizables, 0.0);
418  std::vector<ValueType> denom(num_optimizables, 0.0);
419 
420  ValueType gradNorm = 0.0;
421 
422  // Compute contribution to derivatives
423  for (int i = 0; i < num_optimizables; i++)
424  {
425  // Ground state case
427  {
428  avg_le_der_samp_.at(i) = avg_le_der_samp_.at(i) / static_cast<FullPrecValueType>(w_sum_);
429  avg_der_rat_samp_.at(i) = avg_der_rat_samp_.at(i) / static_cast<FullPrecValueType>(w_sum_);
430 
431  lderivs_.at(i) = static_cast<FullPrecValueType>(2.0) *
432  (avg_le_der_samp_.at(i) - static_cast<FullPrecValueType>(e_avg_) * avg_der_rat_samp_.at(i));
433  if (print_deriv_ == "yes")
434  {
435  app_log() << "Parameter # " << i << " Hamiltonian term: " << avg_le_der_samp_.at(i) << std::endl;
436  app_log() << "Parameter # " << i << " Overlap term: " << avg_der_rat_samp_.at(i) << std::endl;
437  app_log() << "Derivative for param # " << i << " : " << lderivs_.at(i) << std::endl;
438  }
439  }
440  // Excited state case
441  else
442  {
443  avg_numer_der_samp_.at(i) = avg_numer_der_samp_.at(i) / static_cast<FullPrecValueType>(w_sum_);
444  avg_denom_der_samp_.at(i) = avg_denom_der_samp_.at(i) / static_cast<FullPrecValueType>(w_sum_);
445 
446  if (print_deriv_ == "yes")
447  {
448  app_log() << "Parameter # " << i << " Numer Deriv: " << avg_numer_der_samp_.at(i) << std::endl;
449  app_log() << "Parameter # " << i << " Denom Deriv: " << avg_denom_der_samp_.at(i) << std::endl;
450  }
451 
452  numer_term1.at(i) = avg_numer_der_samp_.at(i) * static_cast<FullPrecValueType>(denom_avg_);
453 
454  numer_term2.at(i) = avg_denom_der_samp_.at(i) * static_cast<FullPrecValueType>(numer_avg_);
455 
456  denom.at(i) = denom_avg_ * denom_avg_;
457 
458  lderivs_.at(i) = (numer_term1.at(i) - numer_term2.at(i)) / denom.at(i);
459  }
460 
461  gradNorm += lderivs_.at(i) * lderivs_.at(i);
462  }
463 
464  gradNorm = std::sqrt(gradNorm);
465  app_log() << "Norm of gradient vector is: " << gradNorm << std::endl;
466 
467  // Clear the history vectors for next iteration once sample_finish is done
468  lev_history_.clear();
469  vg_history_.clear();
470  w_history_.clear();
471 
473  {
474  tnv_history_.clear();
475  tdv_history_.clear();
476  }
477 }
Communicate * my_comm_
Communicator handles MPI reduction.
std::vector< std::vector< ValueType > > replica_final_tnv_history_
std::vector< ValueType > final_tar_var_history_
a vector to store the variances of the target function during the descent finalization phase ...
std::vector< ValueType > final_tdv_history_
a history of target function denomerator times the |value/guiding|^2 ratios during the descent finali...
std::vector< ValueType > final_vg_history_
history of sampled |value/guiding|^2 ratios during the descent finalization phase ...
Definition: DescentEngine.h:94
ValueType numer_err_
Standard error of the target function numerator.
Definition: DescentEngine.h:72
int final_descent_num_
Counter for the number of descent steps taken in the finalization phase.
std::vector< ValueType > tnv_history_
a history of target function numerator times the |value/guiding|^2 ratios for one iteration ...
std::vector< std::vector< FullPrecValueType > > replica_denom_der_samp_
Vector for target function denominator parameter derivatives on one thread.
Definition: DescentEngine.h:53
std::vector< ValueType > vg_history_
history of sampled |value/guiding|^2 ratios for one iteration
Definition: DescentEngine.h:89
ValueType target_avg_
Average target function value on a descent step.
Definition: DescentEngine.h:82
std::ostream & app_log()
Definition: OutputManager.h:65
std::vector< std::vector< ValueType > > replica_final_w_history_
ValueType denom_avg_
Average target function denominator on a descent step.
Definition: DescentEngine.h:75
std::vector< ValueType > final_le_avg_history_
a vector to store the averages of the energy during the descent finalization phase ...
std::vector< std::vector< FullPrecValueType > > replica_numer_der_samp_
Vector for target function numerator parameter derivatives on one thread.
Definition: DescentEngine.h:48
int collection_step_
Iteration to start collecting samples for final average and error blocking analysis.
ValueType denom_err_
Standard error of the target function denominator.
Definition: DescentEngine.h:79
std::vector< ValueType > tdv_history_
a history of target function denominator times the |value/guiding|^2 ratios for one iteration ...
std::vector< FullPrecValueType > avg_le_der_samp_
Vector for local energy parameter derivatives.
Definition: DescentEngine.h:36
ValueType e_err_
Standard error of the energy.
Definition: DescentEngine.h:65
std::vector< FullPrecValueType > avg_denom_der_samp_
Vector for target function denominator parameter derivatives.
Definition: DescentEngine.h:51
std::vector< FullPrecValueType > avg_numer_der_samp_
Vector for target function numerator parameter derivatives.
Definition: DescentEngine.h:46
std::vector< FullPrecValueType > avg_der_rat_samp_
Vector for WF parameter derivatives.
Definition: DescentEngine.h:41
bool engine_target_excited_
Whether to target excited state.
std::vector< std::vector< ValueType > > replica_final_lev_history_
void allreduce(T &)
ValueType denom_var_
Variance of the target function denominator.
Definition: DescentEngine.h:77
std::vector< std::vector< ValueType > > replica_tnv_history_
std::vector< ValueType > lderivs_
Vector that stores the final averaged derivatives of the cost function.
std::vector< ValueType > w_history_
history of sampled configuration weights for one iteration
Definition: DescentEngine.h:99
std::vector< std::vector< ValueType > > replica_lev_history_
std::vector< ValueType > final_w_history_
history of sampled configuration weights during descent finalization phase
std::vector< std::vector< ValueType > > replica_tdv_history_
ValueType e_sd_
Standard deviation of the energy.
Definition: DescentEngine.h:63
ValueType target_var_
Variance of the target function.
Definition: DescentEngine.h:84
std::vector< ValueType > final_var_avg_history_
a vector to store the variances of the energy during the descent finalization phase ...
std::vector< ValueType > final_tar_avg_history_
a vector to store the averages of the target function during the descent finalization phase ...
ValueType w_sum_
Total sum of weights.
Definition: DescentEngine.h:56
ValueType numer_avg_
Average target function numerator on a descent step.
Definition: DescentEngine.h:68
ValueType e_avg_
Average energy on a descent step.
Definition: DescentEngine.h:59
MakeReturn< UnaryNode< FnSqrt, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t sqrt(const Vector< T1, C1 > &l)
qmcplusplus::QMCTraits::ValueType ValueType
Definition: DescentEngine.h:30
std::vector< ValueType > final_lev_history_
history of sampled local energies times the |value/guiding|^2 raitos during the descent finalization ...
std::vector< std::vector< FullPrecValueType > > replica_le_der_samp_
Vector for local energy parameter derivatives on one thread.
Definition: DescentEngine.h:38
bool collect_count_
Whether to start collecting samples for the histories in the finalization phase.
void mpi_unbiased_ratio_of_means(int numSamples, std::vector< ValueType > &weights, std::vector< ValueType > &numerSamples, std::vector< ValueType > &denomSamples, ValueType &mean, ValueType &variance, ValueType &stdErr)
Function for computing ratios of the form <f>/<g> as well as the associated variance and standard err...
std::vector< std::vector< FullPrecValueType > > replica_der_rat_samp_
Vector for WF parameter derivatives on one thread.
Definition: DescentEngine.h:43
std::vector< std::vector< ValueType > > replica_final_vg_history_
Definition: DescentEngine.h:95
std::vector< std::vector< ValueType > > replica_final_tdv_history_
qmcplusplus::QMCTraits::FullPrecValueType FullPrecValueType
Definition: DescentEngine.h:29
std::vector< ValueType > lev_history_
a history of sampled local energies times the |value/guiding|^2 raitos for one iteration ...
ValueType e_var_
Variance of the energy.
Definition: DescentEngine.h:61
std::string print_deriv_
Whether to print out derivative terms for each parameter.
std::vector< std::vector< ValueType > > replica_w_history_
ValueType target_err_
Standard error of the target function.
Definition: DescentEngine.h:86
std::vector< std::vector< ValueType > > replica_vg_history_
Definition: DescentEngine.h:90
ValueType numer_var_
Variance of the target function numerator.
Definition: DescentEngine.h:70
std::vector< ValueType > final_tnv_history_
a history of target function numerator times the |value/guiding|^2 ratios during the descent finaliza...

◆ setDerivs()

void setDerivs ( std::vector< ValueType > &  test_derivs)
inline

Function for setting averaged derivatives, currently only used as part of a unit test of the engine's parameter update.

Definition at line 358 of file DescentEngine.h.

References DescentEngine::lderivs_.

358 { lderivs_ = test_derivs; }
std::vector< ValueType > lderivs_
Vector that stores the final averaged derivatives of the cost function.

◆ setParamVal()

void setParamVal ( int  index,
ValueType  value 
)
inline

Function for setting parameter value, used to keep descent parameter values up to date with changes that occur on BLM steps of hybrid method.

Definition at line 362 of file DescentEngine.h.

References DescentEngine::current_params_.

362 { current_params_[index] = value; }
std::vector< ValueType > current_params_
Vector for storing parameter values for current optimization step.

◆ setStepSize()

DescentEngine::ValueType setStepSize ( int  i)

helper method for seting step sizes for different parameter types in descent optimization

Definition at line 844 of file DescentEngine.cpp.

References DescentEngine::ci_eta_, DescentEngine::descent_num_, DescentEngine::engine_param_names_, DescentEngine::engine_param_types_, DescentEngine::f_eta_, DescentEngine::gauss_eta_, DescentEngine::orb_eta_, DescentEngine::ramp_eta_, DescentEngine::ramp_num_, DescentEngine::tjf_1body_eta_, and DescentEngine::tjf_2body_eta_.

Referenced by DescentEngine::updateParameters().

845 {
846  ValueType type_eta;
847 
848  std::string name = engine_param_names_[i];
849 
850  int type = engine_param_types_[i];
851 
852  // Step sizes are assigned according to parameter type identified from the
853  // variable name.
854  // Other parameter types could be added to this section as other wave function
855  // ansatzes are developed.
856  if ((name.find("uu") != std::string::npos) || (name.find("ud") != std::string::npos))
857  {
858  type_eta = tjf_2body_eta_;
859  }
860  // If parameter name doesn't have "uu" or "ud" in it and is of type 1, assume
861  // it is a 1 body Jastrow parameter.
862  else if (type == 1)
863  {
864  type_eta = tjf_1body_eta_;
865  }
866  else if (name.find("F_") != std::string::npos)
867  {
868  type_eta = f_eta_;
869  }
870  else if (name.find("CIcoeff_") != std::string::npos || name.find("CSFcoeff_") != std::string::npos)
871  {
872  type_eta = ci_eta_;
873  }
874  else if (name.find("orb_rot_") != std::string::npos)
875  {
876  type_eta = orb_eta_;
877  }
878  else if (name.find("g") != std::string::npos)
879  {
880  // Gaussian parameters are rarely optimized in practice but the descent code
881  // allows for it.
882  type_eta = gauss_eta_;
883  }
884  else
885  {
886  // If there is some other parameter type that isn't in one of the categories
887  // with a default/input, use a conservative default step size.
888  type_eta = .001;
889  }
890 
892  {
893  type_eta = type_eta * static_cast<ValueType>((descent_num_ + 1) / (double)ramp_num_);
894  }
895 
896  return type_eta;
897 }
int ramp_num_
Number of steps over which to ramp up step size.
ValueType tjf_2body_eta_
Step sizes for different types of parameters.
std::vector< int > engine_param_types_
std::vector< std::string > engine_param_names_
Vectors of parameter names and types, used in the assignment of step sizes.
qmcplusplus::QMCTraits::ValueType ValueType
Definition: DescentEngine.h:30
bool ramp_eta_
Whether to gradually ramp up step sizes in descent.
int descent_num_
Integer for keeping track of only number of descent steps taken.

◆ setupUpdate()

void setupUpdate ( const optimize::VariableSet my_vars)

helper method for transferring information on parameter names and types to the engine

Definition at line 901 of file DescentEngine.cpp.

References qmcplusplus::app_log(), DescentEngine::current_params_, DescentEngine::engine_param_names_, DescentEngine::engine_param_types_, VariableSet::getType(), VariableSet::name(), DescentEngine::num_params_, DescentEngine::params_copy_, DescentEngine::params_for_diff_, VariableSet::size(), and VariableSet::where().

902 {
903  // omega_ = omega_input;
904 
905  num_params_ = my_vars.size();
906  app_log() << "This is num_params_: " << num_params_ << std::endl;
907  for (int i = 0; i < num_params_; i++)
908  {
909  // app_log() << "Variable #" << i << ": " << my_vars[i] << " with index val:
910  // " << my_vars.where(i) << std::endl;
911  if (my_vars.where(i) != -1)
912  {
913  engine_param_names_.push_back(my_vars.name(i));
914  engine_param_types_.push_back(my_vars.getType(i));
915  params_copy_.push_back(my_vars[i]);
916  current_params_.push_back(my_vars[i]);
917  params_for_diff_.push_back(my_vars[i]);
918  }
919  }
920 }
int num_params_
Number of optimizable parameters.
std::ostream & app_log()
Definition: OutputManager.h:65
std::vector< int > engine_param_types_
std::vector< std::string > engine_param_names_
Vectors of parameter names and types, used in the assignment of step sizes.
std::vector< ValueType > params_copy_
Vector for storing parameter values from previous optimization step.
std::vector< ValueType > current_params_
Vector for storing parameter values for current optimization step.
int where(int i) const
return the locator of the i-th Index
Definition: VariableSet.h:90
int getType(int i) const
get the i-th parameter&#39;s type
Definition: VariableSet.h:204
size_type size() const
return the size
Definition: VariableSet.h:88
const std::string & name(int i) const
return the name of i-th variable
Definition: VariableSet.h:189
std::vector< ValueType > params_for_diff_
Vector for storing parameter values for calculating differences to be given to hybrid method...

◆ storeDerivRecord()

void storeDerivRecord ( )
inline

stores derivatives so they can be used in accelerated descent algorithm on later iterations

Definition at line 304 of file DescentEngine.h.

References DescentEngine::deriv_records_, and DescentEngine::lderivs_.

304 { deriv_records_.push_back(lderivs_); }
std::vector< std::vector< ValueType > > deriv_records_
Vector for storing Lagrangian derivatives from previous optimization steps.
std::vector< ValueType > lderivs_
Vector that stores the final averaged derivatives of the cost function.

◆ storeVectors()

void storeVectors ( std::vector< ValueType > &  current_params)

Store a vector of parameter differences to be used by the BLM in a hybrid optimization.

Definition at line 924 of file DescentEngine.cpp.

References qmcplusplus::app_log(), DescentEngine::hybrid_blm_input_, DescentEngine::params_for_diff_, and DescentEngine::store_count_.

925 {
926  std::vector<ValueType> row_vec(current_params.size(), 0.0);
927 
928  // Take difference between current parameter values and the values from some
929  // number of
930  // iterations before to be stored as input to BLM.
931  // The current parameter values are then copied to params_for_diff_ to be used
932  // if storeVectors is called again later in the optimization.
933  for (int i = 0; i < current_params.size(); i++)
934  {
935  row_vec[i] = current_params[i] - params_for_diff_[i];
936  params_for_diff_[i] = current_params[i];
937  }
938 
939  // If on first store of descent section, clear anything that was in vector
940  if (store_count_ == 0)
941  {
942  hybrid_blm_input_.clear();
943  hybrid_blm_input_.push_back(row_vec);
944  }
945  else
946  {
947  hybrid_blm_input_.push_back(row_vec);
948  }
949 
950 #if !defined(QMC_COMPLEX)
951  for (int i = 0; i < hybrid_blm_input_.size(); i++)
952  {
953  std::string entry = "";
954  for (int j = 0; j < hybrid_blm_input_.at(i).size(); j++)
955  {
956  entry = entry + std::to_string(hybrid_blm_input_.at(i).at(j)) + ",";
957  }
958  app_log() << "Stored Vector: " << entry << std::endl;
959  }
960 #endif
961  store_count_++;
962 }
int store_count_
Counter of how many vectors have been stored so far.
std::ostream & app_log()
Definition: OutputManager.h:65
std::vector< std::vector< ValueType > > hybrid_blm_input_
Vector for storing the input vectors to the BLM steps of hybrid method.
std::vector< ValueType > params_for_diff_
Vector for storing parameter values for calculating differences to be given to hybrid method...

◆ takeSample()

void takeSample ( const int  replica_id,
const std::vector< FullPrecValueType > &  der_rat_samp,
const std::vector< FullPrecValueType > &  le_der_samp,
const std::vector< FullPrecValueType > &  ls_der_samp,
ValueType  vgs_samp,
ValueType  weight_samp 
)

Function that Take Sample Data from the Host Code.

Parameters
[in]der_rat_samp<n|Psi_i>/<n|Psi> (i = 0 (|Psi>), 1, ... N_var )
[in]le_der_samp<n|H|Psi_i>/<n|Psi> (i = 0 (|Psi>), 1, ... N_var )
[in]ls_der_samp<|S^2|Psi_i>/<n|Psi> (i = 0 (|Psi>), 1, ... N_var )
[in]vgs_samp|<n|value_fn>/<n|guiding_fn>|^2
[in]weight_sampweight for this sample

Definition at line 214 of file DescentEngine.cpp.

References DescentEngine::collect_count_, DescentEngine::collection_step_, DescentEngine::engine_target_excited_, DescentEngine::final_descent_num_, qmcplusplus::n, DescentEngine::omega_, DescentEngine::replica_denom_der_samp_, DescentEngine::replica_der_rat_samp_, DescentEngine::replica_final_lev_history_, DescentEngine::replica_final_tdv_history_, DescentEngine::replica_final_tnv_history_, DescentEngine::replica_final_vg_history_, DescentEngine::replica_final_w_history_, DescentEngine::replica_le_der_samp_, DescentEngine::replica_lev_history_, DescentEngine::replica_numer_der_samp_, DescentEngine::replica_tdv_history_, DescentEngine::replica_tnv_history_, DescentEngine::replica_vg_history_, and DescentEngine::replica_w_history_.

Referenced by DescentEngineHandle::takeSample().

220 {
221 
222  const size_t num_optimizables = der_rat_samp.size() - 1;
223 
224  ValueType etmp = static_cast<ValueType>(le_der_samp.at(0));
225 
226 
227  // Store a history of samples for the current iteration
228  replica_lev_history_[replica_id].push_back(etmp * vgs_samp);
229  replica_vg_history_[replica_id].push_back(vgs_samp);
230  replica_w_history_[replica_id].push_back(static_cast<ValueType>(1.0));
231 
232  // If on descent finalizing section and past the collection step, store each
233  // local energy for this iteration
235  {
236  replica_final_lev_history_[replica_id].push_back(etmp * vgs_samp);
237  replica_final_vg_history_[replica_id].push_back(vgs_samp);
238  replica_final_w_history_[replica_id].push_back(1.0);
239  }
240 
241  // Ground State Case
243  {
244  for (int i = 0; i < num_optimizables; i++)
245  {
246  replica_le_der_samp_[replica_id].at(i) += le_der_samp.at(i + 1) * static_cast<FullPrecValueType>(vgs_samp);
247  replica_der_rat_samp_[replica_id].at(i) += der_rat_samp.at(i + 1) * static_cast<FullPrecValueType>(vgs_samp);
248  }
249  }
250  // Excited State Case
251  else
252  {
253  ValueType n;
254  ValueType d;
255 
256  n = (omega_ - etmp) * vgs_samp;
257  d = (omega_ * omega_ - static_cast<ValueType>(2) * omega_ * etmp + etmp * etmp) * vgs_samp;
258 
259  replica_tnv_history_[replica_id].push_back(n);
260  replica_tdv_history_[replica_id].push_back(d);
261 
263  {
264  replica_final_tnv_history_[replica_id].push_back(n);
265  replica_final_tdv_history_[replica_id].push_back(d);
266  }
267 
268 
269  for (int i = 0; i < num_optimizables; i++)
270  {
271  // Combination of derivative ratios for target function numerator <psi |
272  // omega -H | psi>/<psi | psi>
273  replica_numer_der_samp_[replica_id].at(i) += static_cast<FullPrecValueType>(2) *
274  (static_cast<FullPrecValueType>(omega_) * der_rat_samp.at(i + 1) - le_der_samp.at(i + 1)) *
275  static_cast<FullPrecValueType>(vgs_samp);
276 
277  // Combination of derivative ratios for target function denominator <psi |
278  // (omega -H)^2 | psi>/<psi | psi>
279  replica_denom_der_samp_[replica_id].at(i) += static_cast<FullPrecValueType>(2) *
280  ((static_cast<FullPrecValueType>(omega_) * der_rat_samp.at(i + 1) - le_der_samp.at(i + 1)) *
281  (static_cast<FullPrecValueType>(omega_) * der_rat_samp.at(0) - le_der_samp.at(0))) *
282  static_cast<FullPrecValueType>(vgs_samp);
283  }
284  }
285 }
std::vector< std::vector< ValueType > > replica_final_tnv_history_
int final_descent_num_
Counter for the number of descent steps taken in the finalization phase.
std::vector< std::vector< FullPrecValueType > > replica_denom_der_samp_
Vector for target function denominator parameter derivatives on one thread.
Definition: DescentEngine.h:53
std::vector< std::vector< ValueType > > replica_final_w_history_
std::vector< std::vector< FullPrecValueType > > replica_numer_der_samp_
Vector for target function numerator parameter derivatives on one thread.
Definition: DescentEngine.h:48
int collection_step_
Iteration to start collecting samples for final average and error blocking analysis.
bool engine_target_excited_
Whether to target excited state.
std::vector< std::vector< ValueType > > replica_final_lev_history_
std::vector< std::vector< ValueType > > replica_tnv_history_
std::vector< std::vector< ValueType > > replica_lev_history_
std::vector< std::vector< ValueType > > replica_tdv_history_
qmcplusplus::QMCTraits::ValueType ValueType
Definition: DescentEngine.h:30
ValueType omega_
Value of omega in excited state functional.
std::vector< std::vector< FullPrecValueType > > replica_le_der_samp_
Vector for local energy parameter derivatives on one thread.
Definition: DescentEngine.h:38
bool collect_count_
Whether to start collecting samples for the histories in the finalization phase.
std::vector< std::vector< FullPrecValueType > > replica_der_rat_samp_
Vector for WF parameter derivatives on one thread.
Definition: DescentEngine.h:43
std::vector< std::vector< ValueType > > replica_final_vg_history_
Definition: DescentEngine.h:95
std::vector< std::vector< ValueType > > replica_final_tdv_history_
qmcplusplus::QMCTraits::FullPrecValueType FullPrecValueType
Definition: DescentEngine.h:29
std::vector< std::vector< ValueType > > replica_w_history_
std::vector< std::vector< ValueType > > replica_vg_history_
Definition: DescentEngine.h:90

◆ targetingExcited()

bool targetingExcited ( ) const
inline

Returns whether an excited state is being targeted.

Definition at line 347 of file DescentEngine.h.

References DescentEngine::engine_target_excited_.

347 { return engine_target_excited_; }
bool engine_target_excited_
Whether to target excited state.

◆ updateParameters()

void updateParameters ( )

helper method for updating parameter values with descent

Definition at line 532 of file DescentEngine.cpp.

References qmcplusplus::abs(), qmcplusplus::app_log(), DescentEngine::ci_eta_, DescentEngine::collect_count_, DescentEngine::collection_step_, DescentEngine::compute_step_, DescentEngine::computeFinalizationUncertainties(), DescentEngine::current_params_, DescentEngine::denom_records_, DescentEngine::deriv_records_, DescentEngine::derivs_squared_, DescentEngine::descent_num_, qmcplusplus::Units::charge::e, DescentEngine::engine_target_excited_, qmcplusplus::exp(), DescentEngine::f_eta_, DescentEngine::final_descent_num_, DescentEngine::final_le_avg_history_, DescentEngine::final_lev_history_, DescentEngine::final_tar_avg_history_, DescentEngine::final_tar_var_history_, DescentEngine::final_tdv_history_, DescentEngine::final_tnv_history_, DescentEngine::final_var_avg_history_, DescentEngine::final_vg_history_, DescentEngine::final_w_history_, DescentEngine::flavor_, qmcplusplus::isnan(), DescentEngine::lambda_, DescentEngine::num_params_, DescentEngine::numer_records_, DescentEngine::orb_eta_, DescentEngine::params_copy_, qmcplusplus::pow(), DescentEngine::setStepSize(), sign(), qmcplusplus::sqrt(), DescentEngine::taus_, DescentEngine::tjf_1body_eta_, and DescentEngine::tjf_2body_eta_.

533 {
534  app_log() << "Number of Parameters: " << num_params_ << std::endl;
535  app_log() << "Descent Number: " << descent_num_ << std::endl;
536 
537  app_log() << "Finalization Descent Num (should be zero if not on last section): " << final_descent_num_ << std::endl;
539  {
540  app_log() << "Should be storing history, length of history on one process is: " << final_lev_history_.size()
541  << std::endl;
542  }
543 
544  app_log() << "Parameter Type step sizes: "
545  << " TJF_2Body_eta=" << tjf_2body_eta_ << " TJF_1Body_eta=" << tjf_1body_eta_ << " F_eta=" << f_eta_
546  << " CI_eta=" << ci_eta_ << " Orb_eta=" << orb_eta_ << std::endl;
547 
548  // Get set of derivatives for current (kth) optimization step
549  std::vector<ValueType> cur_deriv_set = deriv_records_.at(deriv_records_.size() - 1);
550  std::vector<ValueType> prev_deriv_set;
551  if (!taus_.empty())
552  {
553  // Get set of derivatives for previous (k-1th) optimization step
554  prev_deriv_set = deriv_records_.at(deriv_records_.size() - 2);
555  }
556 
557  ValueType denom;
558  ValueType numer;
559  ValueType v;
560  ValueType cor_numer;
561  ValueType cor_v;
562 
563  ValueType epsilon = 1e-8;
564  ValueType type_eta = 0;
565 
566  ValueType tau = 0;
567  // Update parameters according to specified flavor of gradient descent method
568 
569  // RMSprop corresponds to the method used by Booth and co-workers
570  if (flavor_.compare("RMSprop") == 0)
571  {
572  app_log() << "Using RMSprop" << std::endl;
573 
574  // To match up with Booth group paper notation, prevLambda is lambda_k-1,
575  // curLambda is lambda_k, nextLambda is lambda_k+1
576  ValueType cur_lambda = static_cast<ValueType>(.5) +
577  static_cast<ValueType>(.5) *
578  std::sqrt(static_cast<ValueType>(1.0) + static_cast<ValueType>(4.0) * lambda_ * lambda_);
579  ValueType next_lambda = static_cast<ValueType>(.5) +
580  static_cast<ValueType>(.5) *
581  std::sqrt(static_cast<ValueType>(1.0) + static_cast<ValueType>(4.0) * cur_lambda * cur_lambda);
582  ValueType gamma = (static_cast<ValueType>(1.0) - cur_lambda) / next_lambda;
583 
584  // Define damping factor that turns off acceleration of the algorithm
585  // small value of d corresponds to quick damping and effectively using
586  // steepest descent
587  ValueType d = 100;
588  ValueType decay_factor = std::exp(-(static_cast<ValueType>(1.0) / d) * (static_cast<ValueType>(descent_num_)));
589  gamma = gamma * decay_factor;
590 
591  ValueType rho = .9;
592 
593  for (int i = 0; i < num_params_; i++)
594  {
595  ValueType cur_square = std::pow(cur_deriv_set.at(i), static_cast<RealType>(2));
596 
597  // Need to calculate step size tau for each parameter inside loop
598  // In RMSprop, the denominator of the step size depends on a a running
599  // average of past squares of the parameter derivative
600  if (derivs_squared_.size() < num_params_)
601  {
602  cur_square = std::pow(cur_deriv_set.at(i), static_cast<RealType>(2));
603  }
604  else if (derivs_squared_.size() >= num_params_)
605  {
606  cur_square = rho * derivs_squared_.at(i) +
607  (static_cast<ValueType>(1.0) - rho) * std::pow(cur_deriv_set.at(i), static_cast<RealType>(2));
608  }
609 
610  denom = std::sqrt(cur_square + epsilon);
611 
612  // The numerator of the step size is set according to parameter type based
613  // on input choices
614  type_eta = this->setStepSize(i);
615  tau = type_eta / denom;
616 
617  // Include an additional factor to cause step size to eventually decrease
618  // to 0 as number of steps taken increases
619  ValueType step_lambda = .1;
620 
621  ValueType step_decay_denom = static_cast<ValueType>(1.0) + step_lambda * static_cast<ValueType>(descent_num_);
622  tau = tau / step_decay_denom;
623 
624  // Update parameter values
625  // If case corresponds to being after the first descent step
626  if (descent_num_ > 0)
627  {
628  ValueType old_tau = taus_.at(i);
629 
630  current_params_[i] = (static_cast<ValueType>(1.0) - gamma) * (current_params_[i] - tau * cur_deriv_set[i]) +
631  gamma * (params_copy_[i] - old_tau * prev_deriv_set[i]);
632  }
633  else
634  {
635  tau = type_eta;
636 
637  current_params_[i] = current_params_[i] - tau * cur_deriv_set[i];
638  }
639 
640  if (taus_.size() < num_params_)
641  {
642  // For the first optimization step, need to add to the vectors
643  taus_.push_back(tau);
644  derivs_squared_.push_back(cur_square);
645  }
646  else
647  {
648  // When not on the first step, can overwrite the previous stored values
649  taus_[i] = tau;
650  derivs_squared_[i] = cur_square;
651  }
652 
654  }
655 
656  // Store current (kth) lambda value for next optimization step
657  lambda_ = cur_lambda;
658  }
659  // Random uses only the sign of the parameter derivatives and takes a step of
660  // random size within a range.
661  else if (flavor_.compare("Random") == 0)
662  {
663  app_log() << "Using Random" << std::endl;
664 
665  for (int i = 0; i < num_params_; i++)
666  {
667  denom = 1;
668  ValueType alpha = (static_cast<ValueType>(rand() / RAND_MAX));
669  ValueType sign = std::abs(cur_deriv_set[i]) / cur_deriv_set[i];
671  {
672  app_log() << "Got a nan, choosing sign randomly with 50-50 probability" << std::endl;
673 
674  ValueType t = (static_cast<ValueType>(rand() / RAND_MAX));
675  if (std::real(t) > std::real(.5))
676  {
677  sign = 1;
678  }
679  else
680  {
681  sign = -1;
682  }
683  }
684  app_log() << "This is random alpha: " << alpha << " with sign: " << sign << std::endl;
685 
686  current_params_.at(i) = current_params_.at(i) - tau * alpha * sign;
687  }
688  }
689 
690  else
691  {
692  // ADAM method
693  if (flavor_.compare("ADAM") == 0)
694  {
695  app_log() << "Using ADAM" << std::endl;
696 
697  for (int i = 0; i < num_params_; i++)
698  {
699  ValueType cur_square = std::pow(cur_deriv_set.at(i), static_cast<RealType>(2));
700  ValueType beta1 = .9;
701  ValueType beta2 = .99;
702  if (descent_num_ == 0)
703  {
704  numer_records_.push_back(0);
705  denom_records_.push_back(0);
706  }
707  numer = beta1 * numer_records_[i] + (static_cast<ValueType>(1.0) - beta1) * cur_deriv_set[i];
708  v = beta2 * denom_records_[i] + (static_cast<ValueType>(1.0) - beta2) * cur_square;
709 
710  cor_numer = numer / (static_cast<ValueType>(1.0) - std::pow(beta1, static_cast<RealType>(descent_num_ + 1)));
711  cor_v = v / (static_cast<ValueType>(1.0) - std::pow(beta2, static_cast<RealType>(descent_num_ + 1)));
712 
713  denom = std::sqrt(cor_v) + epsilon;
714 
715  type_eta = this->setStepSize(i);
716  tau = type_eta / denom;
717 
718  current_params_.at(i) = current_params_.at(i) - tau * cor_numer;
719 
720  if (taus_.size() < num_params_)
721  {
722  // For the first optimization step, need to add to the vectors
723  taus_.push_back(tau);
724  derivs_squared_.push_back(cur_square);
725  denom_records_[i] = v;
726  numer_records_[i] = numer;
727  }
728  else
729  {
730  // When not on the first step, can overwrite the previous stored
731  // values
732  taus_[i] = tau;
733  derivs_squared_[i] = cur_square;
734  denom_records_[i] = v;
735  numer_records_[i] = numer;
736  }
737 
738  params_copy_[i] = current_params_.at(i);
739  }
740  }
741  // AMSGrad method, similar to ADAM except for form of the step size
742  // denominator
743  else if (flavor_.compare("AMSGrad") == 0)
744  {
745  app_log() << "Using AMSGrad" << std::endl;
746 
747  for (int i = 0; i < num_params_; i++)
748  {
749  ValueType cur_square = std::pow(cur_deriv_set.at(i), static_cast<RealType>(2));
750  ValueType beta1 = .9;
751  ValueType beta2 = .99;
752  if (descent_num_ == 0)
753  {
754  numer_records_.push_back(0);
755  denom_records_.push_back(0);
756  }
757 
758  numer = beta1 * numer_records_[i] + (static_cast<ValueType>(1.0) - beta1) * cur_deriv_set[i];
759  v = beta2 * denom_records_[i] + (static_cast<ValueType>(1.0) - beta2) * cur_square;
760  v = std::max(std::real(denom_records_[i]), std::real(v));
761 
762  denom = std::sqrt(v) + epsilon;
763  type_eta = this->setStepSize(i);
764  tau = type_eta / denom;
765 
766  current_params_.at(i) = current_params_.at(i) - tau * numer;
767 
768  if (taus_.size() < num_params_)
769  {
770  // For the first optimization step, need to add to the vectors
771  taus_.push_back(tau);
772  derivs_squared_.push_back(cur_square);
773  denom_records_[i] = v;
774  numer_records_[i] = numer;
775  }
776  else
777  {
778  // When not on the first step, can overwrite the previous stored
779  // values
780  taus_[i] = tau;
781  derivs_squared_[i] = cur_square;
782  denom_records_[i] = v;
783  numer_records_[i] = numer;
784  }
785 
786  params_copy_[i] = current_params_.at(i);
787  }
788  }
789  }
790 
791  descent_num_++;
792  if (collect_count_)
793  {
795 
796  // Start computing averages and uncertainties from stored history.
797  // This should be done only near the end of the descent finalization section
798  // as it is unnecessary earlier.
800  {
801  app_log() << "Computing average energy and its variance over stored "
802  "steps and its standard error"
803  << std::endl;
804 
805  ValueType collected_steps = final_le_avg_history_.size();
806  ValueType final_e_sum =
807  std::accumulate(final_le_avg_history_.begin(), final_le_avg_history_.end(), static_cast<ValueType>(0.0));
808  ValueType final_e_avg = final_e_sum / collected_steps;
809 
810  ValueType final_var_sum =
811  std::accumulate(final_var_avg_history_.begin(), final_var_avg_history_.end(), static_cast<ValueType>(0.0));
812  ValueType final_var_avg = final_var_sum / collected_steps;
813 
814  app_log() << "Final average energy: " << final_e_avg << std::endl;
815  app_log() << "Final varaince: " << final_var_avg << std::endl;
816 
818 
819  // Do the same for the target function during excited state optimizations.
821  {
822  app_log() << "Computing average target function over stored steps and "
823  "its standard error"
824  << std::endl;
825 
826  ValueType final_tar_sum =
827  std::accumulate(final_tar_avg_history_.begin(), final_tar_avg_history_.end(), static_cast<ValueType>(0.0));
828 
829  ValueType final_tar_avg = final_tar_sum / collected_steps;
830 
831  ValueType final_tar_var_sum =
832  std::accumulate(final_tar_var_history_.begin(), final_tar_var_history_.end(), static_cast<ValueType>(0.0));
833  ValueType final_tar_var_avg = final_tar_var_sum / collected_steps;
834 
835  app_log() << "Final average target function: " << final_tar_avg << std::endl;
836  app_log() << "Final target function varaince: " << final_tar_var_avg << std::endl;
838  }
839  }
840  }
841 }
int num_params_
Number of optimizable parameters.
std::vector< std::vector< ValueType > > deriv_records_
Vector for storing Lagrangian derivatives from previous optimization steps.
std::vector< ValueType > final_tar_var_history_
a vector to store the variances of the target function during the descent finalization phase ...
std::vector< ValueType > final_tdv_history_
a history of target function denomerator times the |value/guiding|^2 ratios during the descent finali...
std::vector< ValueType > final_vg_history_
history of sampled |value/guiding|^2 ratios during the descent finalization phase ...
Definition: DescentEngine.h:94
QMCTraits::RealType real
int final_descent_num_
Counter for the number of descent steps taken in the finalization phase.
std::vector< ValueType > denom_records_
Vector for storing step size denominator values from previous optimization step.
std::vector< ValueType > taus_
Vector for storing step sizes from previous optimization step.
ValueType tjf_2body_eta_
Step sizes for different types of parameters.
MakeReturn< UnaryNode< FnFabs, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t abs(const Vector< T1, C1 > &l)
std::ostream & app_log()
Definition: OutputManager.h:65
std::vector< ValueType > final_le_avg_history_
a vector to store the averages of the energy during the descent finalization phase ...
int collection_step_
Iteration to start collecting samples for final average and error blocking analysis.
bool engine_target_excited_
Whether to target excited state.
MakeReturn< BinaryNode< FnPow, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t, typename CreateLeaf< Vector< T2, C2 > >::Leaf_t > >::Expression_t pow(const Vector< T1, C1 > &l, const Vector< T2, C2 > &r)
int compute_step_
Iteration to start computing averages and errors from the stored values during the finalization phase...
std::vector< ValueType > final_w_history_
history of sampled configuration weights during descent finalization phase
ValueType setStepSize(int i)
helper method for seting step sizes for different parameter types in descent optimization ...
std::vector< ValueType > params_copy_
Vector for storing parameter values from previous optimization step.
std::vector< ValueType > current_params_
Vector for storing parameter values for current optimization step.
void computeFinalizationUncertainties(std::vector< ValueType > &weights, std::vector< ValueType > &numerSamples, std::vector< ValueType > &denomSamples)
Compute uncertainties for energy/target function and variance over a history of samples from a set of...
std::vector< ValueType > final_var_avg_history_
a vector to store the variances of the energy during the descent finalization phase ...
std::vector< ValueType > final_tar_avg_history_
a vector to store the averages of the target function during the descent finalization phase ...
MakeReturn< UnaryNode< FnExp, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t exp(const Vector< T1, C1 > &l)
double sign(double x)
Definition: Standard.h:73
MakeReturn< UnaryNode< FnSqrt, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t sqrt(const Vector< T1, C1 > &l)
QMCTraits::RealType RealType
qmcplusplus::QMCTraits::ValueType ValueType
Definition: DescentEngine.h:30
std::vector< ValueType > final_lev_history_
history of sampled local energies times the |value/guiding|^2 raitos during the descent finalization ...
std::vector< ValueType > derivs_squared_
Vector for storing running average of squares of the derivatives.
int descent_num_
Integer for keeping track of only number of descent steps taken.
bool collect_count_
Whether to start collecting samples for the histories in the finalization phase.
ValueType lambda_
Parameter for accelerated descent recursion relation.
std::string flavor_
What variety of gradient descent will be used.
std::vector< ValueType > numer_records_
Vector for storing step size numerator values from previous optimization step.
bool isnan(float a)
return true if the value is NaN.
Definition: math.cpp:18
std::vector< ValueType > final_tnv_history_
a history of target function numerator times the |value/guiding|^2 ratios during the descent finaliza...

Member Data Documentation

◆ avg_denom_der_samp_

std::vector<FullPrecValueType> avg_denom_der_samp_
private

Vector for target function denominator parameter derivatives.

Definition at line 51 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ avg_der_rat_samp_

std::vector<FullPrecValueType> avg_der_rat_samp_
private

Vector for WF parameter derivatives.

Definition at line 41 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ avg_le_der_samp_

std::vector<FullPrecValueType> avg_le_der_samp_
private

Vector for local energy parameter derivatives.

Definition at line 36 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ avg_numer_der_samp_

std::vector<FullPrecValueType> avg_numer_der_samp_
private

Vector for target function numerator parameter derivatives.

Definition at line 46 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ ci_eta_

◆ collect_count_

bool collect_count_ = false
private

Whether to start collecting samples for the histories in the finalization phase.

Definition at line 236 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), DescentEngine::processXML(), DescentEngine::sample_finish(), DescentEngine::takeSample(), and DescentEngine::updateParameters().

◆ collection_step_

int collection_step_
private

Iteration to start collecting samples for final average and error blocking analysis.

Definition at line 228 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), DescentEngine::processXML(), DescentEngine::sample_finish(), DescentEngine::takeSample(), and DescentEngine::updateParameters().

◆ compute_step_

int compute_step_
private

Iteration to start computing averages and errors from the stored values during the finalization phase.

Definition at line 232 of file DescentEngine.h.

Referenced by DescentEngine::processXML(), and DescentEngine::updateParameters().

◆ current_params_

std::vector<ValueType> current_params_
private

Vector for storing parameter values for current optimization step.

Definition at line 165 of file DescentEngine.h.

Referenced by DescentEngine::retrieveNewParams(), DescentEngine::setParamVal(), DescentEngine::setupUpdate(), and DescentEngine::updateParameters().

◆ denom_avg_

ValueType denom_avg_
private

Average target function denominator on a descent step.

Definition at line 75 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ denom_err_

ValueType denom_err_
private

Standard error of the target function denominator.

Definition at line 79 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ denom_records_

std::vector<ValueType> denom_records_
private

Vector for storing step size denominator values from previous optimization step.

Definition at line 172 of file DescentEngine.h.

Referenced by DescentEngine::updateParameters().

◆ denom_var_

ValueType denom_var_
private

Variance of the target function denominator.

Definition at line 77 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ deriv_records_

std::vector<std::vector<ValueType> > deriv_records_
private

Vector for storing Lagrangian derivatives from previous optimization steps.

Definition at line 168 of file DescentEngine.h.

Referenced by DescentEngine::storeDerivRecord(), and DescentEngine::updateParameters().

◆ derivs_squared_

std::vector<ValueType> derivs_squared_
private

Vector for storing running average of squares of the derivatives.

Definition at line 183 of file DescentEngine.h.

Referenced by DescentEngine::updateParameters().

◆ descent_num_

int descent_num_
private

Integer for keeping track of only number of descent steps taken.

Definition at line 186 of file DescentEngine.h.

Referenced by DescentEngine::DescentEngine(), DescentEngine::getDescentNum(), DescentEngine::setStepSize(), and DescentEngine::updateParameters().

◆ e_avg_

ValueType e_avg_
private

Average energy on a descent step.

Definition at line 59 of file DescentEngine.h.

Referenced by DescentEngine::getEnergy(), DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ e_err_

ValueType e_err_
private

Standard error of the energy.

Definition at line 65 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ e_sd_

ValueType e_sd_
private

Standard deviation of the energy.

Definition at line 63 of file DescentEngine.h.

Referenced by DescentEngine::getSD(), DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ e_var_

ValueType e_var_
private

Variance of the energy.

Definition at line 61 of file DescentEngine.h.

Referenced by DescentEngine::getVariance(), DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ engine_param_names_

std::vector<std::string> engine_param_names_
private

Vectors of parameter names and types, used in the assignment of step sizes.

Definition at line 213 of file DescentEngine.h.

Referenced by DescentEngine::setStepSize(), and DescentEngine::setupUpdate().

◆ engine_param_types_

std::vector<int> engine_param_types_
private

Definition at line 214 of file DescentEngine.h.

Referenced by DescentEngine::setStepSize(), and DescentEngine::setupUpdate().

◆ engine_target_excited_

◆ f_eta_

◆ final_descent_num_

int final_descent_num_ = 0
private

Counter for the number of descent steps taken in the finalization phase.

Definition at line 239 of file DescentEngine.h.

Referenced by DescentEngine::getFinalDescentNum(), DescentEngine::prepareStorage(), DescentEngine::sample_finish(), DescentEngine::takeSample(), and DescentEngine::updateParameters().

◆ final_le_avg_history_

std::vector<ValueType> final_le_avg_history_
private

a vector to store the averages of the energy during the descent finalization phase

Definition at line 117 of file DescentEngine.h.

Referenced by DescentEngine::sample_finish(), and DescentEngine::updateParameters().

◆ final_lev_history_

std::vector<ValueType> final_lev_history_
private

history of sampled local energies times the |value/guiding|^2 raitos during the descent finalization phase

Definition at line 112 of file DescentEngine.h.

Referenced by DescentEngine::sample_finish(), and DescentEngine::updateParameters().

◆ final_tar_avg_history_

std::vector<ValueType> final_tar_avg_history_
private

a vector to store the averages of the target function during the descent finalization phase

Definition at line 144 of file DescentEngine.h.

Referenced by DescentEngine::sample_finish(), and DescentEngine::updateParameters().

◆ final_tar_var_history_

std::vector<ValueType> final_tar_var_history_
private

a vector to store the variances of the target function during the descent finalization phase

Definition at line 147 of file DescentEngine.h.

Referenced by DescentEngine::sample_finish(), and DescentEngine::updateParameters().

◆ final_tdv_history_

std::vector<ValueType> final_tdv_history_
private

a history of target function denomerator times the |value/guiding|^2 ratios during the descent finalization phase

Definition at line 139 of file DescentEngine.h.

Referenced by DescentEngine::sample_finish(), and DescentEngine::updateParameters().

◆ final_tnv_history_

std::vector<ValueType> final_tnv_history_
private

a history of target function numerator times the |value/guiding|^2 ratios during the descent finalization phase

Definition at line 130 of file DescentEngine.h.

Referenced by DescentEngine::sample_finish(), and DescentEngine::updateParameters().

◆ final_var_avg_history_

std::vector<ValueType> final_var_avg_history_
private

a vector to store the variances of the energy during the descent finalization phase

Definition at line 121 of file DescentEngine.h.

Referenced by DescentEngine::sample_finish(), and DescentEngine::updateParameters().

◆ final_vg_history_

std::vector<ValueType> final_vg_history_
private

history of sampled |value/guiding|^2 ratios during the descent finalization phase

Definition at line 94 of file DescentEngine.h.

Referenced by DescentEngine::sample_finish(), and DescentEngine::updateParameters().

◆ final_w_history_

std::vector<ValueType> final_w_history_
private

history of sampled configuration weights during descent finalization phase

Definition at line 102 of file DescentEngine.h.

Referenced by DescentEngine::sample_finish(), and DescentEngine::updateParameters().

◆ flavor_

std::string flavor_
private

What variety of gradient descent will be used.

Definition at line 189 of file DescentEngine.h.

Referenced by DescentEngine::processXML(), and DescentEngine::updateParameters().

◆ gauss_eta_

ValueType gauss_eta_
private

Definition at line 195 of file DescentEngine.h.

Referenced by DescentEngine::processXML(), and DescentEngine::setStepSize().

◆ hybrid_blm_input_

std::vector<std::vector<ValueType> > hybrid_blm_input_
private

Vector for storing the input vectors to the BLM steps of hybrid method.

Definition at line 221 of file DescentEngine.h.

Referenced by DescentEngine::retrieveHybridBLM_Input(), and DescentEngine::storeVectors().

◆ lambda_

ValueType lambda_ = 0.0
private

Parameter for accelerated descent recursion relation.

Definition at line 179 of file DescentEngine.h.

Referenced by DescentEngine::updateParameters().

◆ lderivs_

std::vector<ValueType> lderivs_
private

Vector that stores the final averaged derivatives of the cost function.

Definition at line 150 of file DescentEngine.h.

Referenced by DescentEngine::getAveragedDerivatives(), DescentEngine::prepareStorage(), DescentEngine::sample_finish(), DescentEngine::setDerivs(), and DescentEngine::storeDerivRecord().

◆ lev_history_

std::vector<ValueType> lev_history_
private

a history of sampled local energies times the |value/guiding|^2 raitos for one iteration

Definition at line 107 of file DescentEngine.h.

Referenced by DescentEngine::sample_finish().

◆ my_comm_

Communicate* my_comm_
private

Communicator handles MPI reduction.

Definition at line 153 of file DescentEngine.h.

Referenced by DescentEngine::mpi_unbiased_ratio_of_means(), and DescentEngine::sample_finish().

◆ num_params_

int num_params_
private

Number of optimizable parameters.

Definition at line 159 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), DescentEngine::setupUpdate(), and DescentEngine::updateParameters().

◆ numer_avg_

ValueType numer_avg_
private

Average target function numerator on a descent step.

Definition at line 68 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ numer_err_

ValueType numer_err_
private

Standard error of the target function numerator.

Definition at line 72 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ numer_records_

std::vector<ValueType> numer_records_
private

Vector for storing step size numerator values from previous optimization step.

Definition at line 176 of file DescentEngine.h.

Referenced by DescentEngine::updateParameters().

◆ numer_var_

ValueType numer_var_
private

Variance of the target function numerator.

Definition at line 70 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ omega_

ValueType omega_
private

Value of omega in excited state functional.

Definition at line 224 of file DescentEngine.h.

Referenced by DescentEngine::getOmega(), DescentEngine::processXML(), and DescentEngine::takeSample().

◆ orb_eta_

◆ params_copy_

std::vector<ValueType> params_copy_
private

Vector for storing parameter values from previous optimization step.

Definition at line 162 of file DescentEngine.h.

Referenced by DescentEngine::setupUpdate(), and DescentEngine::updateParameters().

◆ params_for_diff_

std::vector<ValueType> params_for_diff_
private

Vector for storing parameter values for calculating differences to be given to hybrid method.

Definition at line 218 of file DescentEngine.h.

Referenced by DescentEngine::setupUpdate(), and DescentEngine::storeVectors().

◆ print_deriv_

std::string print_deriv_
private

Whether to print out derivative terms for each parameter.

Definition at line 242 of file DescentEngine.h.

Referenced by DescentEngine::processXML(), and DescentEngine::sample_finish().

◆ ramp_eta_

bool ramp_eta_
private

Whether to gradually ramp up step sizes in descent.

Definition at line 200 of file DescentEngine.h.

Referenced by DescentEngine::processXML(), and DescentEngine::setStepSize().

◆ ramp_num_

int ramp_num_
private

Number of steps over which to ramp up step size.

Definition at line 203 of file DescentEngine.h.

Referenced by DescentEngine::processXML(), and DescentEngine::setStepSize().

◆ replica_denom_der_samp_

std::vector<std::vector<FullPrecValueType> > replica_denom_der_samp_
private

Vector for target function denominator parameter derivatives on one thread.

Definition at line 53 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), DescentEngine::sample_finish(), and DescentEngine::takeSample().

◆ replica_der_rat_samp_

std::vector<std::vector<FullPrecValueType> > replica_der_rat_samp_
private

Vector for WF parameter derivatives on one thread.

Definition at line 43 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), DescentEngine::sample_finish(), and DescentEngine::takeSample().

◆ replica_final_lev_history_

std::vector<std::vector<ValueType> > replica_final_lev_history_
private

◆ replica_final_tdv_history_

std::vector<std::vector<ValueType> > replica_final_tdv_history_
private

◆ replica_final_tnv_history_

std::vector<std::vector<ValueType> > replica_final_tnv_history_
private

◆ replica_final_vg_history_

std::vector<std::vector<ValueType> > replica_final_vg_history_
private

◆ replica_final_w_history_

std::vector<std::vector<ValueType> > replica_final_w_history_
private

◆ replica_le_der_samp_

std::vector<std::vector<FullPrecValueType> > replica_le_der_samp_
private

Vector for local energy parameter derivatives on one thread.

Definition at line 38 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), DescentEngine::sample_finish(), and DescentEngine::takeSample().

◆ replica_lev_history_

std::vector<std::vector<ValueType> > replica_lev_history_
private

◆ replica_numer_der_samp_

std::vector<std::vector<FullPrecValueType> > replica_numer_der_samp_
private

Vector for target function numerator parameter derivatives on one thread.

Definition at line 48 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), DescentEngine::sample_finish(), and DescentEngine::takeSample().

◆ replica_tdv_history_

std::vector<std::vector<ValueType> > replica_tdv_history_
private

◆ replica_tnv_history_

std::vector<std::vector<ValueType> > replica_tnv_history_
private

◆ replica_vg_history_

std::vector<std::vector<ValueType> > replica_vg_history_
private

◆ replica_w_history_

std::vector<std::vector<ValueType> > replica_w_history_
private

◆ store_count_

int store_count_
private

Counter of how many vectors have been stored so far.

Definition at line 210 of file DescentEngine.h.

Referenced by DescentEngine::DescentEngine(), DescentEngine::resetStorageCount(), and DescentEngine::storeVectors().

◆ store_num_

int store_num_
private

Number of parameter difference vectors stored when descent is used in a hybrid optimization.

Definition at line 207 of file DescentEngine.h.

Referenced by DescentEngine::processXML(), and DescentEngine::retrieveStoreFrequency().

◆ target_avg_

ValueType target_avg_
private

Average target function value on a descent step.

Definition at line 82 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ target_err_

ValueType target_err_
private

Standard error of the target function.

Definition at line 86 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ target_var_

ValueType target_var_
private

Variance of the target function.

Definition at line 84 of file DescentEngine.h.

Referenced by DescentEngine::prepareStorage(), and DescentEngine::sample_finish().

◆ taus_

std::vector<ValueType> taus_
private

Vector for storing step sizes from previous optimization step.

Definition at line 181 of file DescentEngine.h.

Referenced by DescentEngine::updateParameters().

◆ tdv_history_

std::vector<ValueType> tdv_history_
private

a history of target function denominator times the |value/guiding|^2 ratios for one iteration

Definition at line 135 of file DescentEngine.h.

Referenced by DescentEngine::sample_finish().

◆ tjf_1body_eta_

ValueType tjf_1body_eta_
private

◆ tjf_2body_eta_

ValueType tjf_2body_eta_
private

Step sizes for different types of parameters.

Definition at line 192 of file DescentEngine.h.

Referenced by DescentEngine::processXML(), DescentEngine::setStepSize(), and DescentEngine::updateParameters().

◆ tnv_history_

std::vector<ValueType> tnv_history_
private

a history of target function numerator times the |value/guiding|^2 ratios for one iteration

Definition at line 125 of file DescentEngine.h.

Referenced by DescentEngine::sample_finish().

◆ vg_history_

std::vector<ValueType> vg_history_
private

history of sampled |value/guiding|^2 ratios for one iteration

Definition at line 89 of file DescentEngine.h.

Referenced by DescentEngine::sample_finish().

◆ w_history_

std::vector<ValueType> w_history_
private

history of sampled configuration weights for one iteration

Definition at line 99 of file DescentEngine.h.

Referenced by DescentEngine::sample_finish().

◆ w_sum_

ValueType w_sum_
private

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