ALPS MPS Codes
Reference documentation.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
ietl_jacobi_davidson.h
Go to the documentation of this file.
1 /*****************************************************************************
2  *
3  * ALPS MPS DMRG Project
4  *
5  * Copyright (C) 2013 Institute for Theoretical Physics, ETH Zurich
6  * 2011-2011 by Bela Bauer <bauerb@phys.ethz.ch>
7  *
8  * This software is part of the ALPS Applications, published under the ALPS
9  * Application License; you can use, redistribute it and/or modify it under
10  * the terms of the license, either version 1 or (at your option) any later
11  * version.
12  *
13  * You should have received a copy of the ALPS Application License along with
14  * the ALPS Applications; see the file LICENSE.txt. If not, the license is also
15  * available from http://alps.comp-phys.org/.
16  *
17  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
18  * IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
19  * FITNESS FOR A PARTICULAR PURPOSE, TITLE AND NON-INFRINGEMENT. IN NO EVENT
20  * SHALL THE COPYRIGHT HOLDERS OR ANYONE DISTRIBUTING THE SOFTWARE BE LIABLE
21  * FOR ANY DAMAGES OR OTHER LIABILITY, WHETHER IN CONTRACT, TORT OR OTHERWISE,
22  * ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
23  * DEALINGS IN THE SOFTWARE.
24  *
25  *****************************************************************************/
26 
27 #ifndef IETL_JD_SOLVER_H
28 #define IETL_JD_SOLVER_H
29 
31 
32 #include "ietl_lanczos_solver.h"
33 
34 #include "ietl/jacobi.h"
35 #include "ietl/jd.h"
36 
37 template<class Matrix, class SymmGroup>
38 std::pair<double, MPSTensor<Matrix, SymmGroup> >
40  MPSTensor<Matrix, SymmGroup> const & initial,
41  BaseParameters & params,
42  std::vector<MPSTensor<Matrix, SymmGroup> > ortho_vecs = std::vector<MPSTensor<Matrix, SymmGroup> >())
43 {
44  if (initial.num_elements() <= ortho_vecs.size())
45  ortho_vecs.resize(initial.num_elements()-1);
46  // Gram-Schmidt the ortho_vecs
47  for (int n = 1; n < ortho_vecs.size(); ++n)
48  for (int n0 = 0; n0 < n; ++n0)
49  ortho_vecs[n] -= ietl::dot(ortho_vecs[n0], ortho_vecs[n])/ietl::dot(ortho_vecs[n0],ortho_vecs[n0])*ortho_vecs[n0];
50 
51  typedef MPSTensor<Matrix, SymmGroup> Vector;
52  SingleSiteVS<Matrix, SymmGroup> vs(initial, ortho_vecs);
53 
54  ietl::jcd_gmres_solver<SiteProblem<Matrix, SymmGroup>, SingleSiteVS<Matrix, SymmGroup> >
55  jcd_gmres(sp, vs, params["ietl_jcd_gmres"]);
56 
57  ietl::jacobi_davidson<SiteProblem<Matrix, SymmGroup>, SingleSiteVS<Matrix, SymmGroup> >
58  jd(sp, vs, ietl::Smallest);
59 
60  double tol = params["ietl_jcd_tol"];
61  ietl::basic_iteration<double> iter(params["ietl_jcd_maxiter"], tol, tol);
62 
63 // maquis::cout << "Ortho vecs " << ortho_vecs.size() << std::endl;
64  for (int n = 0; n < ortho_vecs.size(); ++n) {
65 // maquis::cout << "Ortho norm " << n << ": " << ietl::two_norm(ortho_vecs[n]) << std::endl;
66  maquis::cout << "Input <MPS|O[" << n << "]> : " << ietl::dot(initial, ortho_vecs[n]) << std::endl;
67  }
68 
69  std::pair<double, Vector> r0 = jd.calculate_eigenvalue(initial, jcd_gmres, iter);
70 
71  for (int n = 0; n < ortho_vecs.size(); ++n)
72  maquis::cout << "Output <MPS|O[" << n << "]> : " << ietl::dot(r0.second, ortho_vecs[n]) << std::endl;
73 
74  maquis::cout << "JCD used " << iter.iterations() << " iterations." << std::endl;
75 
76  return r0;
77 }
78 
79 /*template<class Matrix, class SymmGroup>
80 std::pair<double, MPSTensor<Matrix, SymmGroup> >
81 solve_ietl_new_jd(SiteProblem<Matrix, SymmGroup> & sp,
82  MPSTensor<Matrix, SymmGroup> const & initial,
83  BaseParameters & params)
84 {
85  typedef MPSTensor<Matrix, SymmGroup> Vector;
86  typedef SiteProblem<Matrix, SymmGroup> Operator;
87  typedef SingleSiteVS<Matrix, SymmGroup> Vecspace;
88  Vecspace vs(initial);
89 
90  int n_evals = 1; // number of eigenpairs to be calculated
91  int max_iter = params["ietl_jcd_maxiter"]; // maximal number of iterations
92 
93  int m_max = 40;
94  int m_min = 20;
95 
96  // tolerance
97  double rel_tol = params["ietl_jcd_tol"];
98  double abs_tol = rel_tol;
99 
100  // maximal iterations for the correction equation solver
101  unsigned max_cor_iter = params["ietl_jcd_gmres"];
102 
103  ietl::jd_iteration<double> iter(max_iter, m_min, m_max, rel_tol, abs_tol);
104  ietl::jd<Operator, Vecspace> jd(sp, vs);
105 
106  // the correction equation solver must be an function object
107  ietl::ietl_gmres gmres(max_cor_iter);
108 
109  jd.eigensystem(iter, initial, n_evals, gmres);
110 
111  std::vector<double> evals = jd.eigenvalues();
112 
113  return std::make_pair(evals[0], jd.eigenvector(0));
114 }*/
115 
116 #endif
call IETL Lanczos solver
MPSTensor< Matrix, SymmGroup >::scalar_type dot(MPSTensor< Matrix, SymmGroup > const &x, MPSTensor< Matrix, SymmGroup > const &y)
std::pair< double, MPSTensor< Matrix, SymmGroup > > solve_ietl_jcd(SiteProblem< Matrix, SymmGroup > &sp, MPSTensor< Matrix, SymmGroup > const &initial, BaseParameters &params, std::vector< MPSTensor< Matrix, SymmGroup > > ortho_vecs=std::vector< MPSTensor< Matrix, SymmGroup > >())
std::size_t num_elements() const
Definition: mpstensor.hpp:639