ALPS MPS Codes
Reference documentation.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
super_models_none.hpp
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  * Michele Dolfi <dolfim@phys.ethz.ch>
8  * 2012 by Jan Gukelberger <gukelberger@phys.ethz.ch>
9  *
10  * This software is part of the ALPS Applications, published under the ALPS
11  * Application License; you can use, redistribute it and/or modify it under
12  * the terms of the license, either version 1 or (at your option) any later
13  * version.
14  *
15  * You should have received a copy of the ALPS Application License along with
16  * the ALPS Applications; see the file LICENSE.txt. If not, the license is also
17  * available from http://alps.comp-phys.org/.
18  *
19  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
20  * IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
21  * FITNESS FOR A PARTICULAR PURPOSE, TITLE AND NON-INFRINGEMENT. IN NO EVENT
22  * SHALL THE COPYRIGHT HOLDERS OR ANYONE DISTRIBUTING THE SOFTWARE BE LIABLE
23  * FOR ANY DAMAGES OR OTHER LIABILITY, WHETHER IN CONTRACT, TORT OR OTHERWISE,
24  * ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
25  * DEALINGS IN THE SOFTWARE.
26  *
27  *****************************************************************************/
28 
29 #ifndef SUPER_MODELS_NONE_HPP
30 #define SUPER_MODELS_NONE_HPP
31 
32 #include <sstream>
33 
34 #include "dmrg/models/model.h"
36 
37 
38 /*
39  * Time evolution for Bose-Hubbard density operators.
40  *
41  * The density operator rho is represented as a superstate whose time evolution
42  * is generated by a superoperator according to the Lindblad master equation
43  * d rho/dt = (-i ad H + lind L) rho .
44  * In this framework the time evolution of the density operator is computed
45  * just like that of a pure state
46  * rho(t) = exp{ -i (ad H + i lind L) t } rho(0) .
47  *
48  * Operators are represented in the occupation number basis
49  * O_{ij} = <i|O|j>
50  * with basis states {|i>}, i=0,...,N-1 for each site.
51  * The corresponding superstate is obtained by the mapping
52  * O_{ij} -> O_n, n=i*N+j ,
53  * i.e. the basis for operator states is given by the matrix elements
54  * |i><j|
55  *
56  * Superoperators represent the map
57  * rho -> d rho/dt .
58  * We mostly need the adjoint Hamiltonian superoperator
59  * ad H: rho -> [H,rho]
60  * and the dissipative superoperator
61  * lind L: rho -> 2 L rho L^\dag - L^\dag L rho - rho L^\dag L .
62  * The superoperator representation of any linear map is conveniently computed
63  * by applying the map to all N^2 operator basis states, filling the
64  * corresponding columns of the superoperator.
65  */
66 
67 
68 /// Take a Hamilton operator H and construct its adjoint Hamiltonian
69 /// superoperator, corresponding to the map
70 /// ad H: rho -> [H,rho]
71 template<class Matrix>
72 Matrix adjoint_hamiltonian(const Matrix& h)
73 {
74  const size_t N = num_rows(h);
75  const size_t N2 = N*N;
76  assert( num_cols(h) == N );
77 
78  Matrix adH(N2,N2);
79  for( size_t i = 0; i < N; ++i )
80  {
81  for( size_t j = 0; j < N; ++j )
82  {
83  const size_t n = i + j*N;
84  Matrix phi (N,N); phi(i,j) = 1.; // n'th operator basis state phi
85  Matrix hphi(N,N); gemm(h,phi,hphi); // H.phi
86  Matrix phih(N,N); gemm(phi,h,phih); // phi.H
87  Matrix comm = hphi - phih; // [H,phi]
88 
89  // operator state adH.phi -> n'th column of adH
90  for( size_t ii = 0; ii < N; ++ii )
91  for( size_t jj = 0; jj < N; ++jj )
92  adH(ii+jj*N, n) = comm(ii,jj);
93  }
94  }
95 
96  return adH;
97 }
98 
99 /// Take a Lindblad operator L and construct its contribution to the Liouville
100 /// superoperator, i.e. the map
101 /// lind L: rho -> 2 L rho L^\dag - L^\dag L rho - rho L^\dag L
102 template<class Matrix>
103 Matrix super_lindblad(const Matrix& l)
104 {
105  const size_t N = num_rows(l);
106  const size_t N2 = N*N;
107  assert( num_cols(l) == N );
108 
109  Matrix ldag = transpose(conj(l)); // L^\dag
110  Matrix ldagl(N,N); gemm(ldag,l,ldagl); // L^\dag L
111  Matrix lindL(N2,N2);
112 
113  for( size_t i = 0; i < N; ++i )
114  {
115  for( size_t j = 0; j < N; ++j )
116  {
117  const size_t n = i + j*N;
118  Matrix phi (N,N); phi(i,j) = 1.; // n'th operator basis state phi
119 
120  // phi -> lindL.phi = 2*l.phi.ldag - ldag.l.phi - phi.ldag.l
121  Matrix lphi (N,N); gemm(l,phi,lphi);
122  Matrix lphildag(N,N); gemm(lphi,ldag,lphildag);
123  Matrix ldaglphi(N,N); gemm(ldag,lphi,ldaglphi);
124  Matrix phildagl(N,N); gemm(phi,ldagl,phildagl);
125  Matrix ll = 2*lphildag - ldaglphi - phildagl;
126 
127  // operator state lindL.phi -> n'th column of lindL
128  for( size_t ii = 0; ii < N; ++ii )
129  for( size_t jj = 0; jj < N; ++jj )
130  lindL(ii+jj*N, n) = ll(ii,jj);
131  }
132  }
133 
134  return lindL;
135 }
136 
137 /// Take a Lindblad operator L and construct the map corresponding to left
138 /// multiplication with L
139 /// rho -> L rho
140 template<class Matrix>
141 Matrix super_left(const Matrix& l)
142 {
143  const size_t N = num_rows(l);
144  const size_t N2 = N*N;
145  assert( num_cols(l) == N );
146 
147  Matrix superL(N2,N2);
148 
149  for( size_t i = 0; i < N; ++i )
150  {
151  for( size_t j = 0; j < N; ++j )
152  {
153  const size_t n = i + j*N;
154  Matrix phi(N,N); phi(i,j) = 1.; // n'th operator basis state phi
155 
156  // phi -> L.phi
157  Matrix lphi(N,N); gemm(l,phi,lphi);
158  for( size_t ii = 0; ii < N; ++ii )
159  for( size_t jj = 0; jj < N; ++jj )
160  superL(ii+jj*N, n) = lphi(ii,jj);
161  }
162  }
163 
164  return superL;
165 }
166 
167 /// Take a Lindblad operator L and construct the map corresponding to right
168 /// multiplication with L
169 /// rho -> rho L
170 template<class Matrix>
171 Matrix super_right(const Matrix& l)
172 {
173  const size_t N = num_rows(l);
174  const size_t N2 = N*N;
175  assert( num_cols(l) == N );
176 
177  Matrix superL(N2,N2);
178 
179  for( size_t i = 0; i < N; ++i )
180  {
181  for( size_t j = 0; j < N; ++j )
182  {
183  const size_t n = i + j*N;
184  Matrix phi (N,N); phi(i,j) = 1.; // n'th operator basis state phi
185 
186  // phi -> phi.L
187  Matrix phil(N,N); gemm(phi,l,phil);
188  for( size_t ii = 0; ii < N; ++ii )
189  for( size_t jj = 0; jj < N; ++jj )
190  superL(ii+jj*N, n) = phil(ii,jj);
191  }
192  }
193 
194  return superL;
195 }
196 
197 
198 /// Fuse indices n[i] into one p = \sum_i n[i] d^i
199 template<class T, class A>
200 T fuse(const A& ind, T d)
201 {
202  T fused = 0;
203  T stride = 1;
204  for( T i = 0; i < ind.size(); ++i )
205  {
206  fused += ind[i] * stride;
207  stride *= d;
208  }
209  return fused;
210 }
211 
212 ///
213 template<class T, class A>
214 void unfuse(T fused, T d, A& ind)
215 {
216  for( T i = 0; i < ind.size(); ++i )
217  {
218  ind[i] = fused % d;
219  fused /= d;
220  }
221  assert( fused == 0 );
222 }
223 
224 
225 
226 ///
227 template<class Matrix>
228 Matrix reshape_bond2site(const Matrix& a)
229 {
230  typedef typename Matrix::size_type size_type;
231  size_type d4 = num_rows(a);
232  size_type d = sqrt(sqrt(double(d4)));
233  assert( d4 == num_cols(a) );
234  assert( d4 == d*d*d*d );
235 
236  Matrix b(d4,d4);
237  boost::array<size_type,4> ii, jj, kk, ll;
238  for( size_type i = 0; i < d4; ++i )
239  {
240  unfuse(i,d,ii);
241  for( size_type j = 0; j < d4; ++j )
242  {
243  unfuse(j,d,jj);
244  kk[0] = ii[0]; kk[1] = ii[2]; kk[2] = jj[0]; kk[3] = jj[2];
245  ll[0] = ii[1]; ll[1] = ii[3]; ll[2] = jj[1]; ll[3] = jj[3];
246  b(fuse(kk,d),fuse(ll,d)) = a(i,j);
247  }
248  }
249  return b;
250 }
251 
252 ///
253 template<class Op,class Matrix>
254 std::vector< std::pair<Op,Op> > decompose_bond_super(const Matrix& bondop, const Index<TrivialGroup>& phys)
255 {
256  Matrix rbond = reshape_bond2site(bondop);
257 
258  Matrix U, V;
259  typename alps::numeric::associated_real_diagonal_matrix<Matrix>::type S, Ssqrt;
260  svd(rbond, U, V, S);
261  Ssqrt = sqrt(S);
262 
264  Op left, right;
265  {
266  Matrix tmp;
267  gemm(U, Ssqrt, tmp);
268  left.insert_block(tmp, C, C);
269  }
270  {
271  Matrix tmp;
272  gemm(Ssqrt, V, tmp);
273  right.insert_block(tmp, C, C);
274  }
275 
276  std::vector<Op> leftops = reshape_right_to_list(phys, left);
277  std::vector<Op> rightops = reshape_left_to_list (phys, right);
278  assert(leftops.size() == rightops.size());
279 
280  // discard terms with no weight
281  std::vector< std::pair<Op,Op> > terms;
282  for( unsigned i = 0; i < num_rows(S) && std::abs(S(i,i)) > 1e-10; ++i )
283  {
284  leftops[i].transpose_inplace();
285  rightops[i].transpose_inplace();
286  terms.push_back(std::pair<Op,Op>( leftops[i], rightops[i] ));
287  }
288  return terms;
289 }
290 
291 
292 /* ****************** BOSE-HUBBARD */
293 template<class Matrix>
294 class SuperBoseHubbardNone : public model_impl<Matrix, TrivialGroup>
295 {
297 
298  typedef typename base::table_type table_type;
299  typedef typename base::table_ptr table_ptr;
300  typedef typename base::tag_type tag_type;
301 
302  typedef typename base::term_descriptor term_descriptor;
303  typedef typename base::terms_type terms_type;
304  typedef typename base::op_t op_t;
305  typedef typename base::measurements_type measurements_type;
306 
307  typedef typename base::size_t size_t;
308  typedef typename Matrix::value_type value_type;
309 
310 public:
311  // Dissipation needs complex types, that's why we forward to do_init with a tag
313  : model(model_)
314  , lattice(lat)
315  , tag_handler(new table_type())
316  {
317  do_init(lat,model,typename Matrix::value_type());
318  }
319 
320  void do_init(const Lattice& lat, BaseParameters & model_, double)
321  {
322  throw std::runtime_error("need complex value type");
323  }
324 
325  void do_init(const Lattice& lat, BaseParameters & model_, std::complex<double>)
326  {
327  // retrieve model parameters
328  int Nmax = model["Nmax"];
329  double t = model["t"];
330  double U = model["U"];
331  double mu = model["mu"];
332  double omega = model["omega"];
333  double V = model["V"];
334  double Lambda = model["Lambda"];
335  double Delta = model["Delta"];
336  double Gamma1a = model["Gamma1a"];
337  double Gamma1b = model["Gamma1b"];
338  double Gamma2 = model["Gamma2"];
339  double nbar = model["nbar"];
340 
342  const size_t N = Nmax+1;
343  const size_t N2 = N*N;
344  const std::complex<double> I(0,1);
345 
346  phys_psi.insert(std::make_pair(C, N));
347 
348  phys.insert(std::make_pair(C, N2));
349  ident.insert_block(Matrix::identity_matrix(N2), C, C);
350 
351  // construct basic on-site operators
352  mcount.resize(N,N); minteraction.resize(N,N); mcreate.resize(N,N); mdestroy.resize(N,N);
353  for( int n = 1; n <= Nmax; ++n )
354  {
355  mcount(n,n) = n;
356  if ((n*n-n) != 0)
357  minteraction(n,n) = n*n-n;
358 
359  mcreate(n,n-1) = std::sqrt(value_type(n)); // input n-1, output n
360  mdestroy(n-1,n) = std::sqrt(value_type(n)); // input n, output n-1
361  }
362  Matrix mcreate2 (N,N); gemm(mcreate ,mcreate ,mcreate2 );
363  Matrix mdestroy2(N,N); gemm(mdestroy,mdestroy,mdestroy2);
364 
365  // construct on-site superoperators
366  Matrix screate = adjoint_hamiltonian(mcreate);
367  Matrix sdestroy = adjoint_hamiltonian(mdestroy);
368  Matrix sdrive = adjoint_hamiltonian(mcreate+mdestroy);
369  Matrix spump = adjoint_hamiltonian(mcreate2+mdestroy2);
370  Matrix scount = adjoint_hamiltonian(mcount);
371  Matrix sinteraction = adjoint_hamiltonian(minteraction);
372 
373  Matrix ldestroy = super_lindblad(mdestroy );
374  Matrix lcreate = super_lindblad(mcreate );
375  Matrix ldestroy2 = super_lindblad(mdestroy2);
376 
377  // cast superoperators to op_t
378  create .insert_block(transpose(screate ), C,C);
379  destroy .insert_block(transpose(sdestroy ), C,C);
380  drive .insert_block(transpose(sdrive ), C,C);
381  pump .insert_block(transpose(spump ), C,C);
382  count .insert_block(transpose(scount ), C,C);
383  interaction.insert_block(transpose(sinteraction), C,C);
384 
385  lindDestroy .insert_block(transpose(ldestroy ), C,C);
386  lindCreate .insert_block(transpose(lcreate ), C,C);
387  lindDestroy2.insert_block(transpose(ldestroy2 ), C,C);
388 
389  leftDestroy.insert_block(transpose(super_left(mdestroy)), C,C);
390  rightCreate.insert_block(transpose(super_right(mcreate)), C,C);
391 
392  std::vector< std::pair<op_t,op_t> > hopops = decompose_bond_super<op_t>(adjoint_hamiltonian(kron(mcreate, mdestroy)),phys);
393  std::vector< std::pair<op_t,op_t> > Vops = decompose_bond_super<op_t>(adjoint_hamiltonian(kron(mcount, mcount)),phys);
394 
395 
396 #define REGISTER(op, kind) op ## _tag = tag_handler->register_op(op, kind);
399  REGISTER(destroy, tag_detail::bosonic)
403  REGISTER(interaction, tag_detail::bosonic)
404  REGISTER(lindDestroy, tag_detail::bosonic)
405  REGISTER(lindCreate, tag_detail::bosonic)
406  REGISTER(lindDestroy2,tag_detail::bosonic)
407  REGISTER(leftDestroy, tag_detail::bosonic)
408  REGISTER(rightCreate, tag_detail::bosonic)
409 
410  std::vector< std::pair<tag_type,tag_type> > hopops_tag(hopops.size());
411  for (size_t i=0; i<hopops.size(); ++i)
412  hopops_tag[i] = std::make_pair( tag_handler->register_op(hopops[i].first, tag_detail::bosonic),
413  tag_handler->register_op(hopops[i].second, tag_detail::bosonic) );
414  std::vector< std::pair<tag_type,tag_type> > Vops_tag(Vops.size());
415  for (size_t i=0; i<Vops.size(); ++i)
416  Vops_tag[i] = std::make_pair( tag_handler->register_op(Vops[i].first, tag_detail::bosonic),
417  tag_handler->register_op(Vops[i].second, tag_detail::bosonic) );
418 
419 
420 #undef REGISTER
421 
422  // insert superoperators for each site
423  for( int p=0; p < lat.size(); ++p )
424  {
425  // interaction H_U = U/2 n_i (n_i - 1)
426  if( U != 0 )
427  {
428  term_descriptor term;
429  term.coeff = 0.5*U;
430  term.push_back( boost::make_tuple(p, interaction_tag) );
431  this->terms_.push_back(term);
432  }
433 
434  // site-dependent chemical potential H_mu = -mu n_i
435  // harmonic trap H_omega = (omega x_i)^2/2 n_i
436  double x = p - 0.5*lat.size();
437  double mup = -mu + 0.5*omega*omega*x*x;
438  if( mup != 0 )
439  {
440  term_descriptor term;
441  term.coeff = mup;
442  term.push_back( boost::make_tuple(p, count_tag) );
443  this->terms_.push_back(term);
444  }
445 
446  // drive H_Lambda = Lambda (b_i^\dag + b_i)
447  if( Lambda != 0 )
448  {
449  term_descriptor term;
450  term.coeff = Lambda;
451  term.push_back( boost::make_tuple(p, drive_tag) );
452  this->terms_.push_back(term);
453  }
454 
455  // pump H_Delta = Delta (b_i^\dag^2 + b_i^2)
456  if( Delta != 0 )
457  {
458  term_descriptor term;
459  term.coeff = Delta;
460  term.push_back( boost::make_tuple(p, pump_tag) );
461  this->terms_.push_back(term);
462  }
463 
464  // one-boson dissipation L_{1a} = Gamma_{1a} (1+\bar{n}) lind b_i
465  // = Gamma_{1a} (1+\bar{n}) (2 b_i rho b_i^\dag - b_i^\dag b_i rho - rho b_i^\dag b_i)
466  if( Gamma1a != 0 )
467  {
468  term_descriptor term;
469  term.coeff = I*Gamma1a*(1+nbar);
470  term.push_back( boost::make_tuple(p, lindDestroy_tag) );
471  this->terms_.push_back(term);
472  }
473 
474  // one-boson dissipation at finite temperature (thermal population \bar{n})
475  // L_{1a} = Gamma_{1a} \bar{n} lind b_i^\dag
476  // = Gamma_{1a} \bar{n} (2 b_i^\dag rho b_i - b_i b_i^\dag rho - rho b_i b_i^\dag)
477  if( Gamma1a*nbar != 0 )
478  {
479  term_descriptor term;
480  term.coeff = I*Gamma1a*nbar;
481  term.push_back( boost::make_tuple(p, lindCreate_tag) );
482  this->terms_.push_back(term);
483  }
484 
485  // two-boson dissipation L_2 = Gamma_2/2 lind b_i^2
486  // = Gamma_2/2 (2 b_i^2 rho b_i^\dag^2 - b_i^\dag^2 b_i^2 rho - rho b_i^\dag^2 b_i^2)
487  if( Gamma2 != 0 )
488  {
489  term_descriptor term;
490  term.coeff = I*0.5*Gamma2;
491  term.push_back( boost::make_tuple(p, lindDestroy2_tag) );
492  this->terms_.push_back(term);
493  }
494 
495  // bond terms
496  std::vector<int> neighs = lat.forward(p);
497  for( int n = 0; n < neighs.size(); ++n )
498  {
499  // hopping H_J = -J (b_i^\dag b_{i+1} + b_i b_{i+1}^\dag)
500  for( unsigned i = 0; i < hopops_tag.size(); ++i )
501  {
502  {
503  term_descriptor term;
504  term.coeff = -t;
505  term.push_back( boost::make_tuple(p, hopops_tag[i].first) );
506  term.push_back( boost::make_tuple(neighs[n], hopops_tag[i].second) );
507  this->terms_.push_back(term);
508  }
509  {
510  term_descriptor term;
511  term.coeff = -t;
512  term.push_back( boost::make_tuple(p, hopops_tag[i].second) );
513  term.push_back( boost::make_tuple(neighs[n], hopops_tag[i].first) );
514  this->terms_.push_back(term);
515  }
516  }
517 
518  // nearest-neighbor interaction H_V = V n_i n_{i+1}
519  if( V != 0 )
520  {
521  for( unsigned i = 0; i < Vops_tag.size(); ++i )
522  {
523  {
524  term_descriptor term;
525  term.coeff = V;
526  term.push_back( boost::make_tuple(p, Vops_tag[i].first) );
527  term.push_back( boost::make_tuple(neighs[n],Vops_tag[i].second) );
528  this->terms_.push_back(term);
529  }
530  {
531  term_descriptor term;
532  term.coeff = V;
533  term.push_back( boost::make_tuple(p, Vops_tag[i].second) );
534  term.push_back( boost::make_tuple(neighs[n],Vops_tag[i].first) );
535  this->terms_.push_back(term);
536  }
537  }
538  }
539 
540  // one-boson dissipation L_{1b} = -Gamma_{1b}/2 (
541  // 2 b_{i+1} rho b_i^\dag - b_i^\dag b_{i+1} rho - rho b_i^\dag b_{i+1}
542  // + 2 b_i rho b_{i+1}^\dag - b_i b_{i+1}^\dag rho - rho b_i b_{i+1}^\dag )
543  // = Gamma_{1b}/2 (
544  // (ad b^\dag)_i (b rho)_{i+1} + (b rho)_i (ad b^\dag)_{i+1}
545  // - (ad b)_i (rho b^\dag)_{i+1} - (rho b^\dag)_i (ad b)_{i+1} )
546  if( Gamma1b != 0 )
547  {
548  term_descriptor term;
549  term.coeff = I*Gamma1b/2.;
550  term.push_back( boost::make_tuple(p, create_tag) );
551  term.push_back( boost::make_tuple(neighs[n], leftDestroy_tag) );
552  this->terms_.push_back(term);
553  }
554  if( Gamma1b != 0 )
555  {
556  term_descriptor term;
557  term.coeff = I*Gamma1b/2.;
558  term.push_back( boost::make_tuple(p, leftDestroy_tag) );
559  term.push_back( boost::make_tuple(neighs[n], create_tag) );
560  this->terms_.push_back(term);
561  }
562  if( Gamma1b != 0 )
563  {
564  term_descriptor term;
565  term.coeff = -I*Gamma1b/2.;
566  term.push_back( boost::make_tuple(p, destroy_tag) );
567  term.push_back( boost::make_tuple(neighs[n], rightCreate_tag) );
568  this->terms_.push_back(term);
569  }
570  if( Gamma1b != 0 )
571  {
572  term_descriptor term;
573  term.coeff = -I*Gamma1b/2.;
574  term.push_back( boost::make_tuple(p, rightCreate_tag) );
575  term.push_back( boost::make_tuple(neighs[n], destroy_tag) );
576  this->terms_.push_back(term);
577  }
578  }
579  }
580 
581  }
582 
583  void update(BaseParameters const& p)
584  {
585  // TODO: update this->terms_ with the new parameters
586  throw std::runtime_error("update() not yet implemented for this model.");
587  return;
588  }
589 
590  Index<TrivialGroup> const& phys_dim(size_t type) const
591  {
592  return phys;
593  }
594  tag_type identity_matrix_tag(size_t type) const
595  {
596  return ident_tag;
597  }
598  tag_type filling_matrix_tag(size_t type) const
599  {
600  return identity_matrix_tag(type);
601  }
603  {
605  }
606 
607  measurements_type measurements() const
608  {
610 
611  op_t ident_psi = identity_matrix<Matrix>(phys_psi);
612  op_t count_psi, create_psi, destroy_psi;
613  count_psi.insert_block(mcount, C, C);
614  create_psi.insert_block(transpose(mcreate), C, C);
615  destroy_psi.insert_block(transpose(mdestroy), C, C);
616 
617  typedef std::vector<block_matrix<Matrix, TrivialGroup> > op_vec;
618  typedef std::vector<std::pair<op_vec, bool> > bond_element;
619 
620  measurements_type meas;
621 
622  if (model["MEASURE[Density]"]) {
623  meas.push_back( new measurements::average<Matrix, TrivialGroup>("Density", lattice,
624  op_vec(1,ident_psi), op_vec(1,ident_psi),
625  op_vec(1,count_psi)) );
626  meas[meas.size()-1].set_super_meas(phys_psi);
627  }
628 
629  if (model["MEASURE[Local density]"]) {
630  meas.push_back( new measurements::local<Matrix, TrivialGroup>("Local density", lattice,
631  op_vec(1,ident_psi), op_vec(1,ident_psi),
632  op_vec(1,count_psi)) );
633  meas[meas.size()-1].set_super_meas(phys_psi);
634  }
635 
636  if (model["MEASURE[Local density^2]"]) {
637  op_t count2_psi;
638  gemm(count_psi, count_psi, count2_psi);
639  meas.push_back( new measurements::local<Matrix, TrivialGroup>("Local density^2", lattice,
640  op_vec(1,ident_psi), op_vec(1,ident_psi),
641  op_vec(1,count2_psi)) );
642  meas[meas.size()-1].set_super_meas(phys_psi);
643  }
644 
645  if (model["MEASURE[Onebody density matrix]"]) {
646  bond_element ops;
647  ops.push_back( std::make_pair(op_vec(1,create_psi), false) );
648  ops.push_back( std::make_pair(op_vec(1,destroy_psi), false) );
649  meas.push_back( new measurements::correlations<Matrix, TrivialGroup>("Onebody density matrix", lattice,
650  op_vec(1,ident_psi), op_vec(1,ident_psi),
651  ops, true, false) );
652  meas[meas.size()-1].set_super_meas(phys_psi);
653  }
654 
655  if (model["MEASURE[Density correlation]"]) {
656  bond_element ops;
657  ops.push_back( std::make_pair(op_vec(1,count_psi), false) );
658  ops.push_back( std::make_pair(op_vec(1,count_psi), false) );
659  meas.push_back( new measurements::correlations<Matrix, TrivialGroup>("Density correlation", lattice,
660  op_vec(1,ident_psi), op_vec(1,ident_psi),
661  ops, true, false) );
662  meas[meas.size()-1].set_super_meas(phys_psi);
663  }
664 
665  return meas;
666  }
667 
668  tag_type get_operator_tag(std::string const & name, size_t type) const
669  {
670  throw std::runtime_error("Operator not defined for this model.");
671  return 0;
672  }
673 
674  table_ptr operators_table() const
675  {
676  return tag_handler;
677  }
678 
679 private:
680  BaseParameters & model;
681  Lattice lattice;
682 
683  op_t ident;
684  Matrix mcount, minteraction, mcreate, mdestroy;
685  op_t create, destroy, drive, pump, count, interaction;
686  op_t lindDestroy, lindCreate, lindDestroy2, leftDestroy, rightCreate;
687  Index<TrivialGroup> phys, phys_psi;
688 
689  boost::shared_ptr<TagHandler<Matrix, TrivialGroup> > tag_handler;
690  tag_type ident_tag;
691  tag_type create_tag, destroy_tag, drive_tag, pump_tag, count_tag, interaction_tag;
692  tag_type lindDestroy_tag, lindCreate_tag, lindDestroy2_tag, leftDestroy_tag, rightCreate_tag;
693 };
694 
695 
696 
697 #endif
Matrix adjoint_hamiltonian(const Matrix &h)
boost::ptr_vector< measurement< Matrix, TrivialGroup > > measurements_type
Definition: model.h:58
void update(BaseParameters const &p)
tag_type get_operator_tag(std::string const &name, size_t type) const
tag_type identity_matrix_tag(size_t type) const
definition of Model base class
block_matrix< Matrix, SymmGroup > identity_matrix(Index< SymmGroup > const &size)
size_type insert_block(Matrix const &, charge, charge)
boost::shared_ptr< table_type > table_ptr
Definition: model.h:52
std::vector< block_matrix< Matrix, SymmGroup > > reshape_right_to_list(Index< SymmGroup > const &phys, block_matrix< Matrix, SymmGroup > const &A)
Definition: reshapes.h:498
table_type::tag_type tag_type
Definition: model.h:53
base type for all models
Definition: model.h:47
Matrix super_left(const Matrix &l)
std::vector< term_descriptor > terms_type
Definition: model.h:56
std::size_t num_rows(maquis::dmrg::one_matrix< T > const &m)
Definition: one_matrix.hpp:99
block_matrix< Matrix, SymmGroup > sqrt(block_matrix< Matrix, SymmGroup > m)
Matrix super_lindblad(const Matrix &l)
std::vector< std::pair< Op, Op > > decompose_bond_super(const Matrix &bondop, const Index< TrivialGroup > &phys)
std::vector< block_matrix< Matrix, SymmGroup > > reshape_left_to_list(Index< SymmGroup > const &phys, block_matrix< Matrix, SymmGroup > const &A)
Definition: reshapes.h:562
void do_init(const Lattice &lat, BaseParameters &model_, double)
::term_descriptor< typename Matrix::value_type > term_descriptor
Definition: model.h:55
void do_init(const Lattice &lat, BaseParameters &model_, std::complex< double >)
void svd(block_matrix< Matrix, SymmGroup > const &M, block_matrix< Matrix, SymmGroup > &U, block_matrix< Matrix, SymmGroup > &V, block_matrix< DiagMatrix, SymmGroup > &S)
std::size_t num_cols(maquis::dmrg::one_matrix< T > const &m)
Definition: one_matrix.hpp:102
table_ptr operators_table() const
measurements_type measurements() const
Matrix reshape_bond2site(const Matrix &a)
static const charge IdentityCharge
Definition: none.h:44
#define REGISTER(op, kind)
Matrix super_right(const Matrix &l)
std::vector< pos_t > forward(pos_t site) const
Definition: lattice.h:100
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)
block_matrix< typename maquis::traits::transpose_view< Matrix >::type, SymmGroup > transpose(block_matrix< Matrix, SymmGroup > const &m)
TrivialGroup::charge total_quantum_numbers(BaseParameters &parms) const
pimpl resolved Lattice
Definition: lattice.h:84
tag_type filling_matrix_tag(size_t type) const
void unfuse(T fused, T d, A &ind)
pos_t size() const
Definition: lattice.h:115
SuperBoseHubbardNone(const Lattice &lat, BaseParameters &model_)
Index< TrivialGroup > const & phys_dim(size_t type) const
T fuse(const A &ind, T d)
Fuse indices n[i] into one p = n[i] d^i.