DOCI-Exact  1.0
doci.cpp
Go to the documentation of this file.
1 /* Copyright (C) 2014-2015 Ward Poelmans
2 
3 This file is part of DOCI-Exact.
4 
5 Doci-Exact is free software: you can redistribute it and/or modify
6 it under the terms of the GNU General Public License as published by
7 the Free Software Foundation, either version 3 of the License, or
8 (at your option) any later version.
9 
10 Hubbard-GPU is distributed in the hope that it will be useful,
11 but WITHOUT ANY WARRANTY; without even the implied warranty of
12 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 GNU General Public License for more details.
14 
15 You should have received a copy of the GNU General Public License
16 along with DOCI-Exact. If not, see <http://www.gnu.org/licenses/>.
17 */
18 
19 #include <iostream>
20 #include <sstream>
21 #include <iomanip>
22 #include <chrono>
23 #include <getopt.h>
24 
25 #include "Permutation.h"
26 #include "Molecule.h"
27 #include "DOCIHamtilonian.h"
28 #include "DM2.h"
29 
30 // comment out if you don't need it
31 #include "SymMolecule.h"
32 #include "SimulatedAnnealing.h"
33 #include "Hamiltonian.h"
34 #include "OrbitalTransform.h"
35 #include "UnitaryMatrix.h"
36 #include "LocalMinimizer.h"
37 
50 int main(int argc, char **argv)
51 {
52  using namespace doci;
53  using std::cout;
54  using std::endl;
55 
56  cout.precision(10);
57 
58  std::string integralsfile = "mo-integrals.h5";
59  std::string h5name = "rdm.h5";
60  std::string unitary;
61  bool simanneal = false;
62  bool jacobirots = false;
63  bool random = false;
64 
65  struct option long_options[] =
66  {
67  {"integrals", required_argument, 0, 'i'},
68  {"output", required_argument, 0, 'o'},
69  {"unitary", required_argument, 0, 'u'},
70  {"simulated-annealing", no_argument, 0, 's'},
71  {"jacobi-rotations", no_argument, 0, 'j'},
72  {"random", no_argument, 0, 'r'},
73  {"help", no_argument, 0, 'h'},
74  {0, 0, 0, 0}
75  };
76 
77  int i,j;
78 
79  while( (j = getopt_long (argc, argv, "hi:o:su:jr", long_options, &i)) != -1)
80  switch(j)
81  {
82  case 'h':
83  case '?':
84  cout << "Usage: " << argv[0] << " [OPTIONS]\n"
85  "\n"
86  " -i, --integrals=integrals-file Set the input integrals file\n"
87  " -o, --output=h5-file Set the output filename for the RDM\n"
88  " -s, --simulated-annealing Use stimulated Annealing to find lowest energy\n"
89  " -j, --jacobi-rotations Use Jacobi Rotations to find lowest energy\n"
90  " -u, --unitary Use this unitary to calc energy\n"
91  " -r, --random Use a random unitary as start point\n"
92  " -h, --help Display this help\n"
93  "\n";
94  return 0;
95  break;
96  case 'i':
97  integralsfile = optarg;
98  break;
99  case 'o':
100  h5name = optarg;
101  break;
102  case 's':
103  simanneal = true;
104  break;
105  case 'j':
106  jacobirots = true;
107  break;
108  case 'r':
109  random = true;
110  break;
111  case 'u':
112  unitary = optarg;
113  break;
114  }
115 
116  if(simanneal && jacobirots)
117  {
118  cout << "Cannot do both simulated annealing and jacobi rotations, you have to choose!" << endl;
119 
120  return 2;
121  }
122 
123 
124  if(getenv("SAVE_H5_PATH"))
125  {
126  h5name = getenv("SAVE_H5_PATH");
127  h5name += "/rdm.h5";
128  } else
129  // make sure we have a save path, even if it's not specify already
130  // This will not overwrite an already set SAVE_H5_PATH
131  setenv("SAVE_H5_PATH", "./", 0);
132 
133  cout << "Reading: " << integralsfile << endl;
134  Sym_Molecule mol(integralsfile);
135  auto& ham_ints = mol.getHamObject();
136 
137  if(random)
138  {
139  if(!unitary.empty())
140  cout << "Overriding the input unitary " << unitary << " with a random unitary!" << endl;
141 
142  const simanneal::OptIndex opt(ham_ints);
143 
145  X.fill_random();
147  cout << "U=exp(X), X=" << endl;
148  X.print_unitary();
149 
150  simanneal::OrbitalTransform orbtrans(ham_ints);
151  orbtrans.update_unitary(X, false);
152 
153  std::string filename = getenv("SAVE_H5_PATH");
154  filename += "/random-start-unitary.h5";
155 
156  orbtrans.get_unitary().saveU(filename);
157 
158  cout << "Saving random start point to " << filename << endl;
159  unitary = filename;
160  }
161 
162  if(!simanneal && !jacobirots)
163  {
164  if(!unitary.empty())
165  {
166  cout << "Reading unitary " << unitary << endl;
167 
168  simanneal::OrbitalTransform orbtrans(ham_ints);
169 
170  orbtrans.get_unitary().loadU(unitary);
171  orbtrans.fillHamCI(ham_ints);
172  }
173 
174  DOCIHamiltonian ham(mol);
175 
176  auto start = std::chrono::high_resolution_clock::now();
177  ham.Build();
178  auto end = std::chrono::high_resolution_clock::now();
179 
180  cout << "Building took: " << std::fixed << std::chrono::duration_cast<std::chrono::duration<double,std::ratio<1>>>(end-start).count() << " s" << endl;
181 
182  start = std::chrono::high_resolution_clock::now();
183  auto eig2 = ham.Diagonalize();
184 // auto eig3 = ham.CalcEnergy(5);
185  end = std::chrono::high_resolution_clock::now();
186 
187  cout << "Diagonalization took: " << std::chrono::duration_cast<std::chrono::duration<double,std::ratio<1>>>(end-start).count() << " s" << endl;
188 
189 // cout << "Energy levels:" << endl;
190 // for(auto i=0;i<eig3.size();i++)
191 // cout << i << "\t" << eig3[i] + mol.get_nucl_rep() << endl;
192 
193  cout << "E = " << eig2.first + mol.get_nucl_rep() << endl;
194 
195  DM2 rdm(mol);
196  auto perm = ham.getPermutation();
197  start = std::chrono::high_resolution_clock::now();
198  rdm.Build(perm, eig2.second);
199  end = std::chrono::high_resolution_clock::now();
200 
201  cout << "Building 2DM took: " << std::chrono::duration_cast<std::chrono::duration<double,std::ratio<1>>>(end-start).count() << " s" << endl;
202 
203  DM2 rdm_ham(mol);
204 
205  rdm_ham.BuildHamiltonian(mol);
206 
207  cout << "DM2 Energy = " << rdm.Dot(rdm_ham) + mol.get_nucl_rep() << endl;
208  cout << "DM2 Trace = " << rdm.Trace() << endl;
209 
210  cout << "Writing 2DM to " << h5name << endl;
211 
212  rdm.WriteToFile(h5name);
213 
214  return 0;
215  }
216 
217  if(simanneal)
218  {
219  SimulatedAnnealing opt(mol);
220 
221  if(!unitary.empty())
222  {
223  cout << "Reading unitary " << unitary << endl;
224  opt.getOrbitaltf().get_unitary().loadU(unitary);
225  }
226 
227  opt.Set_start_temp(0.1);
228  opt.Set_delta_temp(0.99);
229  opt.Set_max_angle(1.3);
230  opt.Set_delta_angle(0.999);
231 
232  auto start = std::chrono::high_resolution_clock::now();
233 
234  opt.optimize();
235 
236  auto end = std::chrono::high_resolution_clock::now();
237 
238  cout << "The optimal energy is " << opt.get_energy() << std::endl;
239 
240  cout << "Optimization took: " << std::fixed << std::chrono::duration_cast<std::chrono::duration<double,std::ratio<1>>>(end-start).count() << " s" << endl;
241 
242  opt.calc_new_energy();
243 
244  auto &ham = opt.getHam();
245  auto &opt_mol = ham.getMolecule();
246 
247  auto eigs = ham.Diagonalize();
248  cout << "E = " << eigs.first + mol.get_nucl_rep() << endl;
249 
250  DM2 rdm(opt_mol);
251  auto perm = ham.getPermutation();
252  start = std::chrono::high_resolution_clock::now();
253  rdm.Build(perm, eigs.second);
254  end = std::chrono::high_resolution_clock::now();
255 
256  cout << "Building 2DM took: " << std::chrono::duration_cast<std::chrono::duration<double,std::ratio<1>>>(end-start).count() << " s" << endl;
257 
258  DM2 rdm_ham(opt_mol);
259  rdm_ham.BuildHamiltonian(opt_mol);
260 
261  cout << "DM2 Energy = " << rdm.Dot(rdm_ham) + mol.get_nucl_rep() << endl;
262  cout << "DM2 Trace = " << rdm.Trace() << endl;
263 
264  rdm.WriteToFile(h5name);
265  }
266 
267  if(jacobirots)
268  {
269  LocalMinimizer opt(mol);
270 
271  if(!unitary.empty())
272  {
273  cout << "Reading unitary " << unitary << endl;
274  opt.getOrbitalTf().get_unitary().loadU(unitary);
275  }
276 
277  opt.Minimize();
278 
279  cout << "The optimal energy is " << opt.get_energy() << std::endl;
280 
281  DM2 rdm_ham(opt.getHam());
282  rdm_ham.BuildHamiltonian(opt.getHam());
283  const DM2 &rdm = opt.get_DM2();
284 
285  cout << "DM2 Energy = " << rdm.Dot(rdm_ham) + mol.get_nucl_rep() << endl;
286  cout << "DM2 Trace = " << rdm.Trace() << endl;
287 
288 
289  cout << "Writing 2DM to " << h5name << endl;
290 
291  rdm.WriteToFile(h5name);
292  }
293 
294 /* PSI_C1_Molecule mol(integralsfile);
295  *
296  * cout << "RHF energy = " << mol.HF_Energy() + mol.get_nucl_rep() << endl;
297  *
298  * DOCIHamiltonian ham(mol);
299  *
300  * auto start = std::chrono::high_resolution_clock::now();
301  * if(getenv("READ_SPARSE_H5_FILE"))
302  * {
303  * std::stringstream h5_name;
304  * h5_name << getenv("READ_SPARSE_H5_FILE");
305  * ham.ReadFromFile(h5_name.str());
306  * }
307  * else
308  * ham.Build();
309  * auto end = std::chrono::high_resolution_clock::now();
310  *
311  * if(getenv("SAVE_SPARSE_H5_FILE"))
312  * {
313  * std::stringstream h5_name;
314  * h5_name << getenv("SAVE_SPARSE_H5_FILE");
315  * ham.SaveToFile(h5_name.str());
316  * }
317  *
318  * cout << "Building took: " << std::fixed << std::chrono::duration_cast<std::chrono::duration<double,std::ratio<1>>>(end-start).count() << " s" << endl;
319  *
320  * // auto eig = ham.DiagonalizeFull();
321  * //
322  * // for(unsigned int i=0;i<eig.first.size();i++)
323  * // cout << i << "\t" << eig.first[i] + mol.get_nucl_rep() << endl;
324  *
325  * start = std::chrono::high_resolution_clock::now();
326  * auto eig2 = ham.Diagonalize();
327  * end = std::chrono::high_resolution_clock::now();
328  *
329  * cout << "Diagonalization took: " << std::chrono::duration_cast<std::chrono::duration<double,std::ratio<1>>>(end-start).count() << " s" << endl;
330  *
331  * cout << "E = " << eig2.first + mol.get_nucl_rep() << endl;
332  *
333  * DM2 rdm(mol);
334  * auto perm = ham.getPermutation();
335  * start = std::chrono::high_resolution_clock::now();
336  * rdm.Build(perm, eig2.second);
337  * end = std::chrono::high_resolution_clock::now();
338  *
339  * cout << "Building 2DM took: " << std::chrono::duration_cast<std::chrono::duration<double,std::ratio<1>>>(end-start).count() << " s" << endl;
340  *
341  * DM2 rdm_ham(mol);
342  *
343  * rdm_ham.BuildHamiltonian(mol);
344  *
345  * cout << "DM2 Energy = " << rdm.Dot(rdm_ham) + mol.get_nucl_rep() << endl;
346  * cout << "DM2 Trace = " << rdm.Trace() << endl;
347  *
348  * rdm.WriteToFile(h5name);
349  */
350 
351  return 0;
352 }
353 
354 /* vim: set ts=8 sw=4 tw=0 expandtab :*/
simanneal::OrbitalTransform & getOrbitalTf() const
void Minimize(bool dist_choice=false)
double Trace() const
Definition: DM2.cpp:611
void update_unitary(double *step)
std::pair< double, std::vector< double > > Diagonalize() const
doci::DOCIHamiltonian & getHam() const
int main(int argc, char **argv)
Definition: doci.cpp:50
void fillHamCI(CheMPS2::Hamiltonian &HamCI)
void WriteToFile(const std::string) const
Definition: DM2.cpp:213
Definition: DM2.h:24
Molecule const & getMolecule() const
simanneal::OrbitalTransform & getOrbitaltf() const
const doci::DM2 & get_DM2() const
void Build(Permutation &, std::vector< double > &)
Definition: DM2.cpp:373
void BuildHamiltonian(const Molecule &)
Definition: DM2.cpp:537
void saveU(std::string savename) const
Save the unitary to disk.
doci::Sym_Molecule & getHam() const
void loadU(std::string loadname)
Load the unitary from disk.
UnitaryMatrix & get_unitary()
double get_energy() const
double get_nucl_rep() const
Definition: SymMolecule.cpp:48
CheMPS2::Hamiltonian & getHamObject() const
Definition: SymMolecule.cpp:67
double Dot(const DM2 &) const
Definition: DM2.cpp:593
Permutation const & getPermutation() const