v2DM-DOCI  1.0
doci.cpp
Go to the documentation of this file.
1 /*
2  * @BEGIN LICENSE
3  *
4  * Copyright (C) 2014-2015 Ward Poelmans
5  *
6  * This file is part of v2DM-DOCI.
7  *
8  * v2DM-DOCI is free software: you can redistribute it and/or modify
9  * it under the terms of the GNU General Public License as published by
10  * the Free Software Foundation, either version 3 of the License, or
11  * (at your option) any later version.
12  *
13  * Foobar is distributed in the hope that it will be useful,
14  * but WITHOUT ANY WARRANTY; without even the implied warranty of
15  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  * GNU General Public License for more details.
17  *
18  * You should have received a copy of the GNU General Public License
19  * along with Foobar. If not, see <http://www.gnu.org/licenses/>.
20  *
21  * @END LICENSE
22  */
23 
24 #include <iostream>
25 #include <chrono>
26 #include <getopt.h>
27 #include <mpi.h>
28 #include <signal.h>
29 
30 #include "include.h"
31 #include "BoundaryPoint.h"
32 #include "PotentialReducation.h"
33 // from CheMPS2
34 #include "Hamiltonian.h"
35 
36 #include "SimulatedAnnealing.h"
37 
38 sig_atomic_t stopping = 0;
39 sig_atomic_t stopping_min = 0;
40 
41 int main(int argc,char **argv)
42 {
43  using std::cout;
44  using std::endl;
45  using namespace doci2DM;
47 
48  cout.precision(10);
49 
50  std::string integralsfile = "mo-integrals.h5";
51  bool bp = false;
52  bool pr = false;
53  std::string startinput;
54 
55  struct option long_options[] =
56  {
57  {"integrals", required_argument, 0, 'i'},
58  {"start", required_argument, 0, 's'},
59  {"boundary-point", no_argument, 0, 'b'},
60  {"potential-reduction", no_argument, 0, 'p'},
61  {"help", no_argument, 0, 'h'},
62  {0, 0, 0, 0}
63  };
64 
65  int i,j;
66 
67  while( (j = getopt_long (argc, argv, "hi:bps:", long_options, &i)) != -1)
68  switch(j)
69  {
70  case 'h':
71  case '?':
72  cout << "Usage: " << argv[0] << " [OPTIONS]\n"
73  "\n"
74  " -i, --integrals=integrals-file Set the input integrals file\n"
75  " -b, --boundary-point Use the boundary point method as solver (default)\n"
76  " -p, --potential-reduction Use the potential reduction method as solver\n"
77  " -s, --start Use this a start point for the Simulated Annealing\n"
78  " -h, --help Display this help\n"
79  "\n";
80  return 0;
81  break;
82  case 'i':
83  integralsfile = optarg;
84  break;
85  case 's':
86  startinput = optarg;
87  break;
88  case 'b':
89  bp = true;
90  break;
91  case 'p':
92  pr = true;
93  break;
94  }
95 
96  if(!bp && !pr)
97  {
98  std::cout << "Tell me what to do..." << std::endl;
99  return 0;
100  }
101 
102  cout << "Reading: " << integralsfile << endl;
103 
104  // make sure we have a save path, even if it's not specify already
105  // This will not overwrite an already set SAVE_H5_PATH
106  setenv("SAVE_H5_PATH", "./", 0);
107 
108  cout << "Using save path: " << getenv("SAVE_H5_PATH") << endl;
109 
110  MPI_Init(&argc,&argv);
111 
112  SimulatedAnnealing opt(CheMPS2::Hamiltonian::CreateFromH5(integralsfile));
113 
114  if(bp)
115  opt.UseBoundaryPoint();
116  else if(pr)
117  opt.UsePotentialReduction();
118 
119  auto& ham = opt.getHam();
120 
121  const auto L = ham.getL(); //dim sp hilbert space
122  const auto N = ham.getNe(); //nr of particles
123 
124  cout << "Starting with L=" << L << " N=" << N << endl;
125  auto start = std::chrono::high_resolution_clock::now();
126 
127  opt.Set_start_temp(0.1);
128  opt.Set_delta_temp(0.99);
129  opt.Set_max_angle(1.3);
130  opt.Set_delta_angle(0.999);
131 
132 // opt.get_Optimal_Unitary().loadU("optimale-uni.h5");
133 // opt.calc_new_energy();
134 
135  if(bp)
136  {
137  try {
138  opt.getMethod_BP().set_use_prev_result(true);
139  } catch (const std::bad_cast err)
140  {
141  std::cerr << "Shit just hit the fan!\t" << err.what() << std::endl;
142  }
143  }
144 
145  if(! startinput.empty())
146  opt.getMethod().getRDM().ReadFromFile(startinput.c_str());
147  else
148  opt.calc_energy();
149 
150  opt.optimize_mpi();
151 
152  auto end = std::chrono::high_resolution_clock::now();
153 
154  cout << "Total Runtime: " << std::fixed << std::chrono::duration_cast<std::chrono::duration<double,std::ratio<1>>>(end-start).count() << " s" << endl;
155 
156  int rank;
157  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
158  if(rank == 0)
159  {
160  std::stringstream h5_name1;
161  h5_name1 << getenv("SAVE_H5_PATH") << "/rdm.h5";
162  opt.getMethod().getRDM().WriteToFile(h5_name1.str().c_str());
163 
164  std::stringstream h5_name2;
165  h5_name2 << getenv("SAVE_H5_PATH") << "/optimale-uni.h5";
166  opt.get_Optimal_Unitary().saveU(h5_name2.str());
167  }
168 
169  MPI_Finalize();
170 
171  return 0;
172 }
173 
174 /* vim: set ts=3 sw=3 expandtab :*/
static Hamiltonian CreateFromH5(const string filename)
sig_atomic_t stopping_min
Definition: doci.cpp:39
sig_atomic_t stopping
Definition: doci.cpp:38
int main(int argc, char **argv)
Definition: doci.cpp:41