QMCPACK
QMCFixedSampleLinearOptimizeBatched.h
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////////////////
2 // This file is distributed under the University of Illinois/NCSA Open Source License.
3 // See LICENSE file in top directory for details.
4 //
5 // Copyright (c) 2020 QMCPACK developers.
6 //
7 // File developed by: Jeremy McMinnis, jmcminis@gmail.com, University of Illinois at Urbana-Champaign
8 // Mark A. Berrill, berrillma@ornl.gov, Oak Ridge National Laboratory
9 // Mark Dewing, mdewing@anl.gov, Argonne National Laboratory
10 //
11 // File created by: Jeongnim Kim, jeongnim.kim@gmail.com, University of Illinois at Urbana-Champaign
12 //////////////////////////////////////////////////////////////////////////////////////
13 
14 
15 /** @file QMCFixedSampleLinearOptimizeBatched.h
16  * @brief Definition of QMCDriver which performs VMC and optimization.
17  */
18 #ifndef QMCPLUSPLUS_QMCFSLINEAROPTIMIZATION_BATCHED_H
19 #define QMCPLUSPLUS_QMCFSLINEAROPTIMIZATION_BATCHED_H
20 
25 #include "Optimize/NRCOptimization.h"
26 #include "Optimize/NRCOptimizationFunctionWrapper.h"
27 #ifdef HAVE_LMY_ENGINE
28 #include "formic/utils/matrix.h"
29 #include "formic/utils/lmyengine/engine.h"
30 #endif
33 #include "OutputMatrix.h"
34 #include "LinearMethod.h"
35 
36 namespace qmcplusplus
37 {
38 /** @ingroup QMCDrivers
39  * @brief Implements wave-function optimization
40  *
41  * Optimization by correlated sampling method with configurations
42  * generated from VMC.
43  */
44 
45 
46 ///forward declaration of a cost function
47 class QMCCostFunctionBase;
48 
49 class GradientTest;
50 
51 
53 {
54 public:
55  ///Constructor.
57  QMCDriverInput&& qmcdriver_input,
58  const std::optional<EstimatorManagerInput>& global_emi,
59  VMCDriverInput&& vmcdriver_input,
61  MCPopulation&& population,
62  SampleStack& samples,
63  Communicate* comm);
64 
65  ///Destructor
67 
69 
70  void setWaveFunctionNode(xmlNodePtr cur) { wfNode = cur; }
71 
72  ///Run the Optimization algorithm.
73  bool run() override;
74  ///preprocess xml node
75  void process(xmlNodePtr cur) override;
76  ///process xml node value (parameters for both VMC and OPT) for the actual optimization
77  bool processOptXML(xmlNodePtr cur, const std::string& vmcMove, bool reportH5, bool useGPU);
78 
80 
81  ///common operation to start optimization
82  void start();
83 
84 #ifdef HAVE_LMY_ENGINE
86  void engine_start(cqmc::engine::LMYEngine<ValueType>* EngineObj,
88  std::string MinMethod);
89 #endif
90 
91 
92  ///common operation to finish optimization, used by the derived classes
93  void finish();
94 
95  void generateSamples();
96 
97 
98 private:
99  NRCOptimizationFunctionWrapper<QMCFixedSampleLinearOptimizeBatched> objFuncWrapper_;
100 
101  inline bool ValidCostFunction(bool valid)
102  {
103  if (!valid)
104  app_log() << " Cost Function is Invalid. If this frequently, try reducing the step size of the line minimization "
105  "or reduce the number of cycles. "
106  << std::endl;
107  return valid;
108  }
109 
110  // check if the proposed new cost function value is the best available
111  bool is_best_cost(const int ii,
112  const std::vector<RealType>& cv,
113  const std::vector<double>& sh,
114  const RealType ic) const;
115 
116  // perform the adaptive three-shift update
118 
119  // perform the single-shift update, no sample regeneration
120  bool one_shift_run();
121 
122  // perform optimization using a gradient descent algorithm
123  bool descent_run();
124 
125  // Previous linear optimizers ("quartic" and "rescale")
127 
128 
129 #ifdef HAVE_LMY_ENGINE
130  // use hybrid approach of descent and blocked linear method for optimization
131  bool hybrid_run();
132 #endif
133 
134  // Perform test of gradients
135  bool test_run();
136 
137  std::unique_ptr<GradientTest> testEngineObj;
138 
139 
140  void solveShiftsWithoutLMYEngine(const std::vector<double>& shifts_i,
141  const std::vector<double>& shiffts_s,
142  std::vector<std::vector<RealType>>& parameterDirections);
143 
144 #ifdef HAVE_LMY_ENGINE
145  formic::VarDeps vdeps;
146  cqmc::engine::LMYEngine<ValueType>* EngineObj;
147 #endif
148 
149  //engine for running various gradient descent based algorithms for optimization
150  std::unique_ptr<DescentEngine> descentEngineObj;
151 
152  //engine for controlling a optimization using a hybrid combination of linear method and descent
153  std::unique_ptr<HybridEngine> hybridEngineObj;
154 
155  // prepare a vector of shifts to try
156  std::vector<double> prepare_shifts(const double central_shift) const;
157 
158  // previous update
159 #ifdef HAVE_LMY_ENGINE
160  std::vector<formic::ColVec<double>> previous_update;
161 #endif
162 
164  void print_cost_summary(const double si,
165  const double ss,
166  const RealType mc,
167  const RealType cv,
168  const int ind,
169  const int bi,
170  const bool gu);
171 
172  // ------------------------------------
173  // Used by legacy linear method algos
174 
175  std::vector<RealType> optdir, optparam;
176 
177  ///Number of iterations maximum before generating new configurations.
179 
181  //-------------------------------------
182 
183  /// Choice of eigenvalue solver
184  std::string eigensolver_;
185 
186  ///Number of iterations maximum before generating new configurations.
189  /// the previous best identity shift
191  /// the previous best overlap shift
193  /// current shift_i, shift_s input values
195  /// accept history, remember the last 2 iterations, value 00, 01, 10, 11
196  std::bitset<2> accept_history;
197  /// Shift_s adjustment base
199 
200  // ------------------------------------
201  // Parameters in this struct are used by one or more of the adaptive LM, descent, or hybrid optimizers
202 
203  //String inputs are listed outside the struct
204  ///name of the current optimization method, updated by processOptXML before run
205  std::string MinMethod;
206 
207  //LMY related input
208  struct LMYOptions
209  {
210  /// number of shifts we will try
211  int num_shifts = 3;
212  /// the maximum relative change in the cost function for the adaptive three-shift scheme
214  ///max amount a parameter may change relative to current wave function weight
216  /// the tolerance to cost function increases when choosing the best shift in the adaptive shift method
218  /// the shift_i value that the adaptive shift method should aim for
220  ///whether we are targeting an excited state
221  bool targetExcited = false;
222  ///whether we are doing block algorithm
223  bool block_lm = false;
224  ///number of blocks used in block algorithm
225  int nblocks = 1;
226  ///number of old updates kept
227  int nolds = 1;
228  ///number of directions kept
229  int nkept = 1;
230  ///number of samples to do in correlated sampling part
231  int nsamp_comp = 0;
232  ///the shift to use when targeting an excited state
234  ///whether to do the first part of block lm
235  bool block_first = true;
236  ///whether to do the second part of block lm
237  bool block_second = false;
238  ///whether to do the third part of block lm
239  bool block_third = false;
240  ///whether to filter parameters for the lm
241  bool filter_param = false;
242  ///whether to filter parameters for the lm
243  bool filter_info = false;
244  ///threshold for filtering parameters for the lm
245  double ratio_threshold = 0.0;
246  ///whether to store samples for the lm
247  bool store_samples = false;
248  ///type of the previous optimization method, updated by processOptXML before run
250  ///type of the current optimization method, updated by processOptXML before run
252  ///whether to use hybrid method
253  bool doHybrid = false;
254  };
255 
256  /// LMY engine related options
258 
259  // ------------------------------------
260 
261  // Test parameter gradients
263 
264  // Output Hamiltonian and overlap matrices
266 
267  // Output Hamiltonian and overlap matrices in HDF format
269 
270  // Flag to open the files on first pass and print header line
272 
275 
276  // Freeze variational parameters. Do not update them during each step.
278 
286 
287  ///xml node to be dumped
288  xmlNodePtr wfNode;
289 
291 
292  ///target cost function to optimize
293  std::unique_ptr<QMCCostFunctionBase> optTarget;
294 
295  ///vmc engine
296  std::unique_ptr<VMCBatched> vmcEngine;
297 
300 
301  /// This is retained in order to construct and reconstruct the vmcEngine.
302  const std::optional<EstimatorManagerInput> global_emi_;
303 };
304 } // namespace qmcplusplus
305 #endif
OptimizerType previous_optimizer_type
type of the previous optimization method, updated by processOptXML before run
helper functions for EinsplineSetBuilder
Definition: Configuration.h:43
Input representation for VMC driver class runtime parameters.
std::vector< double > prepare_shifts(const double central_shift) const
returns a vector of three shift values centered around the provided shift.
OptimizerType current_optimizer_type
type of the current optimization method, updated by processOptXML before run
QMCDriverNew Base class for Unified Drivers.
Definition: QMCDriverNew.h:75
std::ostream & app_log()
Definition: OutputManager.h:65
Declaration of QMCDriverNew.
LinearMethod related functions.
class ProjectData
Definition: ProjectData.h:36
A set of light weight walkers that are carried between driver sections and restart.
void print_cost_summary_header()
prints a header for the summary of each shift&#39;s result
Timer accumulates time and call counts.
Definition: NewTimer.h:135
std::unique_ptr< QMCCostFunctionBase > optTarget
target cost function to optimize
Wrapping information on parallelism.
Definition: Communicate.h:68
RealType shift_i_input
current shift_i, shift_s input values
bool one_shift_run()
Performs one iteration of the linear method using an adaptive scheme that tries three different shift...
RealType max_param_change
max amount a parameter may change relative to current wave function weight
class to handle a set of parameters
Definition: ParameterSet.h:27
bool is_best_cost(const int ii, const std::vector< RealType > &cv, const std::vector< double > &sh, const RealType ic) const
Returns whether the proposed new cost is the best compared to the others.
QMCFixedSampleLinearOptimizeBatched(const ProjectData &project_data, QMCDriverInput &&qmcdriver_input, const std::optional< EstimatorManagerInput > &global_emi, VMCDriverInput &&vmcdriver_input, WalkerConfigurations &wc, MCPopulation &&population, SampleStack &samples, Communicate *comm)
Constructor.
QTBase::ValueType ValueType
Definition: Configuration.h:60
RealType omega_shift
the shift to use when targeting an excited state
void print_cost_summary(const double si, const double ss, const RealType mc, const RealType cv, const int ind, const int bi, const bool gu)
prints a summary of the computed cost for the given shift
RealType cost_increase_tol
the tolerance to cost function increases when choosing the best shift in the adaptive shift method ...
int nsamp_comp
number of samples to do in correlated sampling part
NRCOptimizationFunctionWrapper< QMCFixedSampleLinearOptimizeBatched > objFuncWrapper_
RealType target_shift_i
the shift_i value that the adaptive shift method should aim for
std::string MinMethod
name of the current optimization method, updated by processOptXML before run
void solveShiftsWithoutLMYEngine(const std::vector< double > &shifts_i, const std::vector< double > &shiffts_s, std::vector< std::vector< RealType >> &parameterDirections)
For each set of shifts, solves the linear method eigenproblem by building and diagonalizing the matri...
bool processOptXML(xmlNodePtr cur, const std::string &vmcMove, bool reportH5, bool useGPU)
process xml node value (parameters for both VMC and OPT) for the actual optimization ...
RealType max_relative_cost_change
the maximum relative change in the cost function for the adaptive three-shift scheme ...
LatticeGaussianProduct::ValueType ValueType
QMCTraits::RealType RealType
Definition: QMCDriverNew.h:78
void finish()
common operation to finish optimization, used by the derived classes
int nstabilizers
Number of iterations maximum before generating new configurations.
std::bitset< 2 > accept_history
accept history, remember the last 2 iterations, value 00, 01, 10, 11
int Max_iterations
Number of iterations maximum before generating new configurations.
double ratio_threshold
threshold for filtering parameters for the lm
Output the Hamiltonian and overlap matrices from linear method.
const std::optional< EstimatorManagerInput > global_emi_
This is retained in order to construct and reconstruct the vmcEngine.
void process(xmlNodePtr cur) override
preprocess xml node
Input representation for Driver base class runtime parameters.