ALPS MPS Codes
Reference documentation.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
mps_mpo_ops.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 MPS_MPO_OPS_H
28 #define MPS_MPO_OPS_H
29 
30 #include "dmrg/mp_tensors/mps.h"
31 #include "dmrg/mp_tensors/mpo.h"
32 
35 
36 #include "dmrg/utils/utils.hpp"
37 #include "utils/traits.hpp"
38 
39 template<class Matrix, class SymmGroup>
40 std::vector<Boundary<Matrix, SymmGroup> >
42 {
43  assert(mpo.length() == mps.length());
44  std::size_t L = mps.length();
45 
46  std::vector<Boundary<Matrix, SymmGroup> > left_(L+1);
47  left_[0] = mps.left_boundary();
48 
49  for (int i = 0; i < L; ++i) {
50  left_[i+1] = contraction::overlap_mpo_left_step(mps[i], mps[i], left_[i], mpo[i]);
51  }
52  return left_;
53 }
54 
55 template<class Matrix, class SymmGroup>
56 std::vector<Boundary<Matrix, SymmGroup> >
58 {
59  assert(mpo.length() == mps.length());
60  std::size_t L = mps.length();
61 
62  std::vector<Boundary<Matrix, SymmGroup> > right_(L+1);
63  right_[L] = mps.right_boundary();
64 
65  for (int i = L-1; i >= 0; --i) {
66  right_[i] = contraction::overlap_mpo_right_step(mps[i], mps[i], right_[i+1], mpo[i]);
67  }
68  return right_;
69 }
70 
71 template<class Matrix, class SymmGroup>
72 double expval(MPS<Matrix, SymmGroup> const & mps, MPO<Matrix, SymmGroup> const & mpo, int d)
73 {
74  if (d == 0) {
75  std::vector<Boundary<Matrix, SymmGroup> > left_ = left_mpo_overlaps(mps, mpo);
76  assert( check_real(left_[mps.length()][0].trace()) );
77  return maquis::real(left_[mps.length()][0].trace());
78  } else {
79  std::vector<Boundary<Matrix, SymmGroup> > right_ = right_mpo_overlaps(mps, mpo);
80  assert( check_real(right_[0][0].trace()) );
81  return maquis::real(right_[0][0].trace());
82  }
83 }
84 
85 template<class Matrix, class SymmGroup>
86 double expval(MPS<Matrix, SymmGroup> const & mps, MPO<Matrix, SymmGroup> const & mpo,
87  bool verbose = false)
88 {
89  assert(mpo.length() == mps.length());
90  std::size_t L = mps.length();
91 
93 
94  semi_parallel_for (/*locale::compact(L)*/, std::size_t i = 0; i < L; ++i) {
95  if (verbose)
96  maquis::cout << "expval site " << (size_t)i << std::endl;
97  left = contraction::overlap_mpo_left_step(mps[i], mps[i], left, mpo[i]);
98  }
99 
100  return maquis::real(left[0].trace());
101 }
102 
103 template<class Matrix, class SymmGroup>
104 std::vector<typename MPS<Matrix, SymmGroup>::scalar_type> multi_expval(MPS<Matrix, SymmGroup> const & mps,
105  MPO<Matrix, SymmGroup> const & mpo)
106 {
107  assert(mpo.length() == mps.length());
108  std::size_t L = mps.length();
109 
111 
112  for (int i = 0; i < L; ++i) {
113  left = contraction::overlap_mpo_left_step(mps[i], mps[i], left, mpo[i]);
114  }
115 
116  return left.traces();
117 }
118 
119 template<class Matrix, class SymmGroup>
121 {
122  std::size_t L = mps.length();
123 
125  left.insert_block(Matrix(1, 1, 1), SymmGroup::IdentityCharge, SymmGroup::IdentityCharge);
126 
127  semi_parallel_for (/*removed...*/, std::size_t i = 0; i < L; ++i) {
128  MPSTensor<Matrix, SymmGroup> cpy = mps[i];
129  left = contraction::overlap_left_step(mps[i], cpy, left); // serial
130  }
131 
132  return trace(left);
133 }
134 
135 template<class Matrix, class SymmGroup>
137  MPS<Matrix, SymmGroup> const & mps2)
138 {
139  assert(mps1.length() == mps2.length());
140 
141  std::size_t L = mps1.length();
142 
144  left.insert_block(Matrix(1, 1, 1), SymmGroup::IdentityCharge, SymmGroup::IdentityCharge);
145 
146  semi_parallel_for (/*locale::compact(L)*/, std::size_t i = 0; i < L; ++i) {
147  left = contraction::overlap_left_step(mps1[i], mps2[i], left);
148  }
149 
150  return trace(left);
151 }
152 
153 template<class Matrix, class SymmGroup>
154 std::vector<typename MPS<Matrix, SymmGroup>::scalar_type> multi_overlap(MPS<Matrix, SymmGroup> const & mps1,
155  MPS<Matrix, SymmGroup> const & mps2)
156 {
157  // assuming mps2 to have `correct` shape, i.e. left size=1, right size=1
158  // mps1 more generic, i.e. left size=1, right size arbitrary
159 
160 
161  assert(mps1.length() == mps2.length());
162 
163  std::size_t L = mps1.length();
164 
166  left.insert_block(Matrix(1, 1, 1), SymmGroup::IdentityCharge, SymmGroup::IdentityCharge);
167 
168  for (int i = 0; i < L; ++i) {
169  left = contraction::overlap_left_step(mps1[i], mps2[i], left);
170  }
171 
172  assert(left.right_basis().sum_of_sizes() == 1);
173  std::vector<typename MPS<Matrix, SymmGroup>::scalar_type> vals;
174  vals.reserve(left.left_basis().sum_of_sizes());
175  for (int n=0; n<left.n_blocks(); ++n)
176  for (int i=0; i<left.left_basis()[n].second; ++i)
177  vals.push_back( left[n](i,0) );
178 
179  return vals;
180 }
181 
182 //typedef std::vector< std::vector< std::pair<std::string, double> > > entanglement_spectrum_type;
183 typedef std::vector< std::pair<std::vector<std::string>, std::vector<double> > > entanglement_spectrum_type;
184 template<class Matrix, class SymmGroup>
185 std::vector<double>
187  std::vector<int> * measure_es_where = NULL,
188  entanglement_spectrum_type * spectra = NULL) // to be optimized later
189 {
190  std::size_t L = mps.length();
191  std::vector<double> ret;
192 
193  MPS<Matrix, SymmGroup> const& constmps = mps;
194 
196 
197  if (spectra != NULL)
198  spectra->clear();
199 
200  mps.canonize(0);
201  for (std::size_t p = 1; p < L; ++p)
202  {
205 
206  constmps[p-1].make_left_paired();
207  constmps[p].make_right_paired();
208 
209  gemm(constmps[p-1].data(), constmps[p].data(), t);
210 
211  svd(t, u, v, s);
212 
213  std::vector<double> sv = maquis::dmrg::detail::bond_renyi_entropies(s);
214 
215  if (spectra != NULL && measure_es_where != NULL
216  && std::find(measure_es_where->begin(), measure_es_where->end(), p) != measure_es_where->end()) {
217  std::vector< std::string > labels;
218  std::vector< double > values;
219  for (std::size_t k = 0; k < s.n_blocks(); ++k) {
220  std::ostringstream oss_c;
221  oss_c << s.left_basis()[k].first;
222  std::string c_str = oss_c.str();
223  for (std::size_t l = 0; l < s.left_basis()[k].second; ++l) {
224  labels.push_back( c_str );
225  values.push_back( s[k](l,l) );
226  }
227  }
228  spectra->push_back(std::make_pair(labels, values));
229  }
230 
231  double S = 0;
232  if (n == 1) {
233  for (std::vector<double>::const_iterator it = sv.begin();
234  it != sv.end(); ++it)
235  S += *it * log(*it);
236  ret.push_back(-S);
237  } else {
238  for (std::vector<double>::const_iterator it = sv.begin();
239  it != sv.end(); ++it)
240  S += pow(*it, n);
241  ret.push_back(1/(1-n)*log(S));
242  }
243 
244  mps.move_normalization_l2r(p-1, p, DefaultSolver());
245  }
246 
247  return ret;
248 }
249 
250 template<class Matrix, class SymmGroup>
251 std::vector<double>
253 {
254  return calculate_bond_renyi_entropies(mps, 1, NULL);
255 }
256 
257 template<class Matrix, class SymmGroup>
259 {
260  typedef typename SymmGroup::charge charge;
261  charge I = SymmGroup::IdentityCharge;
262  size_t L = mps.length();
263 
264  Index<SymmGroup> phys_rho = phys_psi * adjoin(phys_psi);
265  ProductBasis<SymmGroup> pb(phys_psi, phys_psi, boost::lambda::bind(static_cast<charge(*)(charge, charge)>(SymmGroup::fuse),
266  boost::lambda::_1, -boost::lambda::_2));
267 
268  Matrix identblock(phys_rho.size_of_block(I), 1, 0.);
269  for (int s=0; s<phys_psi.size(); ++s)
270  for (int ss=0; ss<phys_psi[s].second; ++ss) {
271  identblock(pb(phys_psi[s].first, phys_psi[s].first) + ss*phys_psi[s].second+ss, 0) = 1.;
272  }
274  ident.insert_block(identblock, I, I);
275 
276  Index<SymmGroup> trivial_i;
277  trivial_i.insert(std::make_pair(I, 1));
278  MPSTensor<Matrix, SymmGroup> mident(phys_rho, trivial_i, trivial_i);
279  mident.data() = ident;
280 
281  MPS<Matrix,SymmGroup> mps_ident(L);
282  for (int p=0; p<L; ++p)
283  mps_ident[p] = mident;
284 
285  return overlap(mps, mps_ident);
286 }
287 
288 
289 // Specific to Fermi-Hubbard on a Ladder!!
290 template<class Matrix, class SymmGroup>
291 void fix_density(MPS<Matrix, SymmGroup> & mps, std::vector<block_matrix<Matrix, SymmGroup> > const & dens_ops, std::vector<std::vector<double> > const & dens)
292 {
293  assert( mps.size() == dens[0].size() );
294  assert( dens_ops.size() == dens.size() );
295  size_t L = mps.size();
296 
297  mps.normalize_left();
298  mps.canonize(0);
299  for (int p=0; p<L; ++p)
300  {
301 
302 
303  Index<SymmGroup> phys = mps[p].site_dim();
304  typename SymmGroup::charge empty, up, down, updown;
305  empty[0] = 0; empty[1] = 0;
306  up[0] = 1; up[1] = 0;
307  down[0] = 0; down[1] = 1;
308  updown[0] = 1; updown[1] = 1;
309 
310 
312 
313  for (size_t j=0; j<dens.size(); ++j) {
314 
315  MPSTensor<Matrix, SymmGroup> tmp = contraction::local_op(mps[p], dens_ops[j]);
316  double cur_dens = mps[p].scalar_overlap(tmp);
317  maquis::cout << "Density[" << j << "] (before) = " << cur_dens << std::endl;
318  }
319 
320  double a = trace(rho(down, down)) * trace(rho(updown, updown));
321  double b = trace(rho(up, up)) * trace(rho(down, down)) + dens[0][p] * trace(rho(updown, updown)) - dens[1][p] * trace(rho(updown, updown));
322  double c = - dens[1][p] * trace(rho(up, up));
323  double k2 = ( -b + sqrt(b*b - 4*a*c) ) / (2*a);
324 
325  double k1 = dens[0][p] / ( trace(rho(up, up)) + k2*trace(rho(updown,updown)) );
326 
327  double t0 = 0.;
328  t0 += k1*trace( rho(up, up) );
329  t0 += k2*trace( rho(down, down) );
330  t0 += k1*k2*trace( rho(updown, updown) );
331  double k0 = (1.-t0) / trace(rho(empty, empty));
332 
333  maquis::cout << "k0 = " << k0 << std::endl;
334  maquis::cout << "k1 = " << k1 << std::endl;
335  maquis::cout << "k2 = " << k2 << std::endl;
336  assert( k0 > 0 ); // not always the case!!!
337 
338  block_matrix<Matrix, SymmGroup> rescale = identity_matrix<Matrix>(phys);
339  rescale(empty, empty) *= std::sqrt(k0);
340  rescale(up, up) *= std::sqrt(k1);
341  rescale(down, down) *= std::sqrt(k2);
342  rescale(updown, updown) *= std::sqrt(k1*k2);
343 
344  mps[p] = contraction::local_op(mps[p], rescale);
345 
346  {
347  for (size_t j=0; j<dens.size(); ++j) {
348  MPSTensor<Matrix, SymmGroup> tmp = contraction::local_op(mps[p], dens_ops[j]);
349  double meas_dens = mps[p].scalar_overlap(tmp) / mps[p].scalar_norm();
350  maquis::cout << "Density[" << j << "] (after) = " << meas_dens << ", should be " << dens[j][p] << std::endl;
351  }
352  }
353 
354  block_matrix<Matrix, SymmGroup> t_norm = mps[p].normalize_left(DefaultSolver());
355  if (p < L-1)
356  mps[p+1].multiply_from_left(t_norm);
357 
358  }
359 
360 }
361 
362 
363 #endif
void move_normalization_l2r(size_t p1, size_t p2, DecompMethod method=DefaultSolver())
Definition: mps.hpp:167
double expval(MPS< Matrix, SymmGroup > const &mps, MPO< Matrix, SymmGroup > const &mpo, int d)
Definition: mps_mpo_ops.h:72
void normalize_left()
Definition: mps.hpp:124
std::vector< typename MPS< Matrix, SymmGroup >::scalar_type > multi_overlap(MPS< Matrix, SymmGroup > const &mps1, MPS< Matrix, SymmGroup > const &mps2)
Definition: mps_mpo_ops.h:154
size_type n_blocks() const
declaration of the MPS class (vector of MPSTensor)
static MPSTensor< Matrix, SymmGroup > local_op(MPSTensor< Matrix, SymmGroup > const &mps, block_matrix< Matrix, SymmGroup > const &op)
Definition: contractions.h:774
std::vector< double > calculate_bond_entropies(MPS< Matrix, SymmGroup > &mps)
Definition: mps_mpo_ops.h:252
std::vector< double > bond_renyi_entropies(const block_matrix< alps::numeric::diagonal_matrix< T >, SymmGroup > &set)
std::vector< Boundary< Matrix, SymmGroup > > left_mpo_overlaps(MPS< Matrix, SymmGroup > const &mps, MPO< Matrix, SymmGroup > const &mpo)
Definition: mps_mpo_ops.h:41
std::vector< typename MPS< Matrix, SymmGroup >::scalar_type > multi_expval(MPS< Matrix, SymmGroup > const &mps, MPO< Matrix, SymmGroup > const &mpo)
Definition: mps_mpo_ops.h:104
size_type insert_block(Matrix const &, charge, charge)
std::size_t size_of_block(charge c) const
Boundary< Matrix, SymmGroup > left_boundary() const
Definition: mps.hpp:253
static Boundary< OtherMatrix, SymmGroup > overlap_mpo_right_step(MPSTensor< Matrix, SymmGroup > const &bra_tensor, MPSTensor< Matrix, SymmGroup > const &ket_tensor, Boundary< OtherMatrix, SymmGroup > const &right, MPOTensor< Matrix, SymmGroup > const &mpo)
Definition: contractions.h:386
data_type data
static block_matrix< Matrix, SymmGroup > density_matrix(MPSTensor< Matrix, SymmGroup > const &bra_tensor, MPSTensor< Matrix, SymmGroup > const &ket_tensor)
Definition: contractions.h:851
Definition: mpo.h:36
std::vector< std::pair< std::vector< std::string >, std::vector< double > > > entanglement_spectrum_type
Definition: mps_mpo_ops.h:183
Boundary< Matrix, SymmGroup > right_boundary() const
Definition: mps.hpp:266
bool check_real(T x)
Definition: utils.hpp:43
#define semi_parallel_for(constraint,...)
block_matrix< Matrix, SymmGroup > sqrt(block_matrix< Matrix, SymmGroup > m)
std::vector< Boundary< Matrix, SymmGroup > > right_mpo_overlaps(MPS< Matrix, SymmGroup > const &mps, MPO< Matrix, SymmGroup > const &mpo)
Definition: mps_mpo_ops.h:57
static block_matrix< OtherMatrix, SymmGroup > overlap_left_step(MPSTensor< Matrix, SymmGroup > const &bra_tensor, MPSTensor< Matrix, SymmGroup > const &ket_tensor, block_matrix< OtherMatrix, SymmGroup > const &left, block_matrix< OtherMatrix, SymmGroup > *localop=NULL)
Definition: contractions.h:41
block_matrix< Matrix, SymmGroup >::scalar_type trace(block_matrix< Matrix, SymmGroup > const &m)
block_matrix< Matrix, SymmGroup > adjoin(block_matrix< Matrix, SymmGroup > const &m)
MPS< Matrix, SymmGroup >::scalar_type overlap(MPS< Matrix, SymmGroup > const &mps1, MPS< Matrix, SymmGroup > const &mps2)
Definition: mps_mpo_ops.h:136
size_t size() const
Definition: mps.h:57
definition of MPO class (vector of MPOTensor)
Definition: mps.h:40
void svd(block_matrix< Matrix, SymmGroup > const &M, block_matrix< Matrix, SymmGroup > &U, block_matrix< Matrix, SymmGroup > &V, block_matrix< DiagMatrix, SymmGroup > &S)
MPSTensor< Matrix, SymmGroup >::scalar_type scalar_type
Definition: mps.h:51
MPS< Matrix, SymmGroup >::scalar_type norm(MPS< Matrix, SymmGroup > const &mps)
Definition: mps_mpo_ops.h:120
Index< SymmGroup > const & right_basis() const
static Boundary< OtherMatrix, SymmGroup > overlap_mpo_left_step(MPSTensor< Matrix, SymmGroup > const &bra_tensor, MPSTensor< Matrix, SymmGroup > const &ket_tensor, Boundary< OtherMatrix, SymmGroup > const &left, MPOTensor< Matrix, SymmGroup > const &mpo)
Definition: contractions.h:340
Logger< storage::archive > log
Definition: utils.cpp:40
size_t length() const
Definition: mps.h:58
void fix_density(MPS< Matrix, SymmGroup > &mps, std::vector< block_matrix< Matrix, SymmGroup > > const &dens_ops, std::vector< std::vector< double > > const &dens)
Definition: mps_mpo_ops.h:291
std::vector< scalar_type > traces() const
Definition: boundary.h:84
std::size_t insert(std::pair< charge, std::size_t > const &x)
void gemm(block_matrix< Matrix1, SymmGroup > const &A, block_matrix< Matrix2, SymmGroup > const &B, block_matrix< Matrix3, SymmGroup > &C)
void canonize(size_t center, DecompMethod method=DefaultSolver())
Definition: mps.hpp:147
Index< SymmGroup > const & site_dim(size_t i) const
Definition: mps.h:59
std::vector< double > calculate_bond_renyi_entropies(MPS< Matrix, SymmGroup > &mps, double n, std::vector< int > *measure_es_where=NULL, entanglement_spectrum_type *spectra=NULL)
Definition: mps_mpo_ops.h:186
Index< SymmGroup > const & left_basis() const
block_matrix< Matrix, SymmGroup > & data()
Definition: mpstensor.hpp:423
std::size_t size() const
functions to contract tensor network states
alps::numeric::real_type< T >::type real(T f)
Definition: bindings.hpp:38
MPS< Matrix, SymmGroup >::scalar_type dm_trace(MPS< Matrix, SymmGroup > const &mps, Index< SymmGroup > const &phys_psi)
Definition: mps_mpo_ops.h:258
T fuse(const A &ind, T d)
Fuse indices n[i] into one p = n[i] d^i.