CheMPS2
TensorT.cpp
1 /*
2  CheMPS2: a spin-adapted implementation of DMRG for ab initio quantum chemistry
3  Copyright (C) 2013-2016 Sebastian Wouters
4 
5  This program 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 2 of the License, or
8  (at your option) any later version.
9 
10  This program 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 along
16  with this program; if not, write to the Free Software Foundation, Inc.,
17  51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
18 */
19 
20 #include <stdlib.h> /*rand*/
21 #include <algorithm>
22 #include <math.h>
23 
24 #include "TensorT.h"
25 #include "Lapack.h"
26 
27 using std::min;
28 
29 CheMPS2::TensorT::TensorT( const int site_index, const SyBookkeeper * denBK ) : Tensor(){
30 
31  this->index = site_index; //left boundary = index ; right boundary = index+1
32  this->denBK = denBK;
33 
34  AllocateAllArrays();
35 
36 }
37 
38 void CheMPS2::TensorT::AllocateAllArrays(){
39 
40  nKappa = 0;
41  for ( int NL = denBK->gNmin( index ); NL <= denBK->gNmax( index ); NL++ ){
42  for ( int TwoSL = denBK->gTwoSmin( index, NL ); TwoSL <= denBK->gTwoSmax( index, NL ); TwoSL += 2 ){
43  for ( int IL = 0; IL < denBK->getNumberOfIrreps(); IL++ ){
44  const int dimL = denBK->gCurrentDim( index, NL, TwoSL, IL );
45  if ( dimL > 0 ){
46  for ( int NR = NL; NR <= NL+2; NR++ ){
47  const int TwoJ = (( NR == NL + 1 ) ? 1 : 0 );
48  for ( int TwoSR = TwoSL - TwoJ; TwoSR <= TwoSL + TwoJ; TwoSR += 2 ){
49  if ( TwoSR >= 0 ){
50  int IR = (( NR == NL + 1 ) ? Irreps::directProd( IL, denBK->gIrrep( index ) ) : IL );
51  const int dimR = denBK->gCurrentDim( index + 1, NR, TwoSR, IR );
52  if ( dimR > 0 ){
53  nKappa++;
54  }
55  }
56  }
57  }
58  }
59  }
60  }
61  }
62 
63  sectorNL = new int[ nKappa ];
64  sectorNR = new int[ nKappa ];
65  sectorIL = new int[ nKappa ];
66  sectorIR = new int[ nKappa ];
67  sectorTwoSL = new int[ nKappa ];
68  sectorTwoSR = new int[ nKappa ];
69  kappa2index = new int[ nKappa + 1 ];
70  kappa2index[ 0 ] = 0;
71 
72  nKappa = 0;
73  for ( int NL = denBK->gNmin( index ); NL <= denBK->gNmax( index ); NL++ ){
74  for ( int TwoSL = denBK->gTwoSmin( index, NL ); TwoSL <= denBK->gTwoSmax( index, NL ); TwoSL += 2 ){
75  for ( int IL = 0; IL < denBK->getNumberOfIrreps(); IL++ ){
76  const int dimL = denBK->gCurrentDim( index, NL, TwoSL, IL );
77  if ( dimL > 0 ){
78  for ( int NR = NL; NR <= NL+2; NR++ ){
79  const int TwoJ = (( NR == NL + 1 ) ? 1 : 0 );
80  for ( int TwoSR = TwoSL - TwoJ; TwoSR <= TwoSL + TwoJ; TwoSR += 2 ){
81  if ( TwoSR >= 0 ){
82  int IR = (( NR == NL + 1 ) ? Irreps::directProd( IL, denBK->gIrrep( index ) ) : IL );
83  const int dimR = denBK->gCurrentDim( index + 1, NR, TwoSR, IR );
84  if ( dimR > 0 ){
85  sectorNL[ nKappa ] = NL;
86  sectorNR[ nKappa ] = NR;
87  sectorIL[ nKappa ] = IL;
88  sectorIR[ nKappa ] = IR;
89  sectorTwoSL[ nKappa ] = TwoSL;
90  sectorTwoSR[ nKappa ] = TwoSR;
91  kappa2index[ nKappa + 1 ] = kappa2index[ nKappa ] + dimL * dimR;
92  nKappa++;
93  }
94  }
95  }
96  }
97  }
98  }
99  }
100  }
101 
102  storage = new double[ kappa2index[ nKappa ] ];
103 
104 }
105 
107 
108  DeleteAllArrays();
109 
110 }
111 
112 void CheMPS2::TensorT::DeleteAllArrays(){
113 
114  delete [] sectorNL;
115  delete [] sectorNR;
116  delete [] sectorIL;
117  delete [] sectorIR;
118  delete [] sectorTwoSL;
119  delete [] sectorTwoSR;
120  delete [] kappa2index;
121  delete [] storage;
122 
123 }
124 
126 
127  DeleteAllArrays();
128  AllocateAllArrays();
129 
130 }
131 
132 int CheMPS2::TensorT::gNKappa() const { return nKappa; }
133 
134 double * CheMPS2::TensorT::gStorage() { return storage; }
135 
136 int CheMPS2::TensorT::gKappa( const int N1, const int TwoS1, const int I1, const int N2, const int TwoS2, const int I2 ) const{
137 
138  for ( int cnt = 0; cnt < nKappa; cnt++ ){
139  if (( sectorNL[ cnt ] == N1 ) &&
140  ( sectorNR[ cnt ] == N2 ) &&
141  ( sectorIL[ cnt ] == I1 ) &&
142  ( sectorIR[ cnt ] == I2 ) &&
143  ( sectorTwoSL[ cnt ] == TwoS1 ) &&
144  ( sectorTwoSR[ cnt ] == TwoS2 )){ return cnt; }
145  }
146 
147  return -1;
148 
149 }
150 
151 int CheMPS2::TensorT::gKappa2index( const int kappa ) const{ return kappa2index[ kappa ]; }
152 
153 double * CheMPS2::TensorT::gStorage( const int N1, const int TwoS1, const int I1, const int N2, const int TwoS2, const int I2 ){
154 
155  int kappa = gKappa( N1, TwoS1, I1, N2, TwoS2, I2 );
156  if ( kappa == -1 ){ return NULL; }
157  return storage + kappa2index[ kappa ];
158 
159 }
160 
161 int CheMPS2::TensorT::gIndex() const { return index; }
162 
163 const CheMPS2::SyBookkeeper * CheMPS2::TensorT::gBK() const{ return denBK; }
164 
165 void CheMPS2::TensorT::sBK( const SyBookkeeper * newBK ){ denBK = newBK; }
166 
168 
169  for ( int cnt = 0; cnt < kappa2index[ nKappa ]; cnt++ ){
170  storage[ cnt ] = ((double) rand()) / RAND_MAX;
171  }
172 
173 }
174 
175 void CheMPS2::TensorT::number_operator( const double alpha, const double beta ){
176 
177  #pragma omp parallel for schedule(dynamic)
178  for ( int ikappa = 0; ikappa < nKappa; ikappa++ ){
179  int size = kappa2index[ ikappa + 1 ] - kappa2index[ ikappa ];
180  double * array = storage + kappa2index[ ikappa ];
181  double factor = beta + alpha * ( sectorNR[ ikappa ] - sectorNL[ ikappa ] );
182  int inc1 = 1;
183  dscal_( &size, &factor, array, &inc1 );
184  }
185 
186 }
187 
188 void CheMPS2::TensorT::QR(Tensor * Rstorage){
189 
190  //Left normalization occurs in T-convention: no pre or after multiplication
191  //Work per right symmetry sector
192 
193  //PARALLEL
194  #pragma omp parallel for schedule(dynamic)
195  for (int NR=denBK->gNmin(index+1); NR<=denBK->gNmax(index+1); NR++){
196  for (int TwoSR=denBK->gTwoSmin(index+1,NR); TwoSR<=denBK->gTwoSmax(index+1,NR); TwoSR+=2){
197  for (int IR=0; IR<denBK->getNumberOfIrreps(); IR++){
198  int dimR = denBK->gCurrentDim(index+1,NR,TwoSR,IR);
199  if (dimR>0){
200 
201  //Find out the total left dimension
202  int dimLtotal = 0;
203  for (int ikappa=0; ikappa<nKappa; ikappa++){
204  if ((NR==sectorNR[ikappa])&&(TwoSR==sectorTwoSR[ikappa])&&(IR==sectorIR[ikappa])){
205  dimLtotal += denBK->gCurrentDim(index,sectorNL[ikappa],sectorTwoSL[ikappa],sectorIL[ikappa]);
206  }
207  }
208 
209  if (dimLtotal>0){ //Due to the initial truncation, it is possible that dimLtotal is temporarily smaller than dimR ...
210 
211  double * mem = new double[dimLtotal*dimR];
212  //Copy the relevant parts from storage to mem
213  int dimLtotal2 = 0;
214  for (int ikappa=0; ikappa<nKappa; ikappa++){
215  if ((NR==sectorNR[ikappa])&&(TwoSR==sectorTwoSR[ikappa])&&(IR==sectorIR[ikappa])){
216  int dimL = denBK->gCurrentDim(index,sectorNL[ikappa],sectorTwoSL[ikappa],sectorIL[ikappa]);
217  if (dimL>0){
218  for (int l=0; l<dimL; l++){
219  for (int r=0; r<dimR; r++){
220  mem[dimLtotal2 + l + dimLtotal * r] = storage[kappa2index[ikappa] + l + dimL * r];
221  }
222  }
223  dimLtotal2 += dimL;
224  }
225  }
226  }
227 
228  //QR mem --> m = dimLtotal ; n = dimR
229  int info;
230  int minofdims = min(dimR,dimLtotal);
231  double * tau = new double[minofdims];
232  double * work = new double[dimR];
233  dgeqrf_(&dimLtotal,&dimR,mem,&dimLtotal,tau,work,&dimR,&info);
234 
235  //Copy R to Rstorage
236  double * wheretoput = Rstorage->gStorage(NR,TwoSR,IR,NR,TwoSR,IR); //dimR x dimR
237 
238  for (int irow = 0; irow<minofdims; irow++){
239  for (int icol = 0; icol<irow; icol++){
240  wheretoput[irow + dimR * icol] = 0.0;
241  }
242  for (int icol = irow; icol<dimR; icol++){
243  wheretoput[irow + dimR * icol] = mem[irow + dimLtotal * icol];
244  }
245  }
246  for (int irow = minofdims; irow<dimR; irow++){
247  for (int icol = 0; icol<dimR; icol++){
248  wheretoput[irow + dimR * icol] = 0.0;
249  }
250  }
251 
252  //Construct Q
253  dorgqr_(&dimLtotal,&minofdims,&minofdims,mem,&dimLtotal,tau,work,&dimR,&info);
254  if (dimLtotal < dimR){ //if number of cols larger than number of rows, rest of cols zero.
255  for (int irow=0; irow<dimLtotal; irow++){
256  for (int icol=dimLtotal; icol<dimR; icol++){
257  mem[irow + dimLtotal * icol] = 0.0;
258  }
259  }
260  }
261 
262  //Copy from mem to storage
263  dimLtotal2 = 0;
264  for (int ikappa=0; ikappa<nKappa; ikappa++){
265  if ((NR==sectorNR[ikappa])&&(TwoSR==sectorTwoSR[ikappa])&&(IR==sectorIR[ikappa])){
266  int dimL = denBK->gCurrentDim(index,sectorNL[ikappa],sectorTwoSL[ikappa],sectorIL[ikappa]);
267  if (dimL>0){
268  for (int l=0; l<dimL; l++){
269  for (int r=0; r<dimR; r++){
270  storage[kappa2index[ikappa] + l + dimL * r] = mem[dimLtotal2 + l + dimLtotal * r];
271  }
272  }
273  dimLtotal2 += dimL;
274  }
275  }
276  }
277 
278  //Clear the memory
279  delete [] work;
280  delete [] tau;
281  delete [] mem;
282 
283  }
284  }
285  }
286  }
287  }
288 
289 }
290 
291 void CheMPS2::TensorT::LQ(Tensor * Lstorage){
292 
293  //Right normalization occurs in U-convention: pre-multiplication with sqrt{2jR+1/2jL+1} and after multiplication with sqrt{2jL+1/2jR+1}
294  //Work per left symmetry sector
295 
296  //PARALLEL
297  #pragma omp parallel for schedule(dynamic)
298  for (int NL=denBK->gNmin(index); NL<=denBK->gNmax(index); NL++){
299  for (int TwoSL=denBK->gTwoSmin(index,NL); TwoSL<=denBK->gTwoSmax(index,NL); TwoSL+=2){
300  for (int IL=0; IL<denBK->getNumberOfIrreps(); IL++){
301  int dimL = denBK->gCurrentDim(index,NL,TwoSL,IL);
302  if (dimL>0){
303 
304  //Find out the total right dimension
305  int dimRtotal = 0;
306  for (int ikappa=0; ikappa<nKappa; ikappa++){
307  if ((NL==sectorNL[ikappa])&&(TwoSL==sectorTwoSL[ikappa])&&(IL==sectorIL[ikappa])){
308  dimRtotal += denBK->gCurrentDim(index+1,sectorNR[ikappa],sectorTwoSR[ikappa],sectorIR[ikappa]);
309  }
310  }
311 
312  if (dimRtotal>0){ //Due to the initial truncation, it is possible that dimRtotal is temporarily smaller than dimL ...
313 
314  double * mem = new double[dimRtotal*dimL];
315  //Copy the relevant parts from storage to mem & multiply with factor !!
316  int dimRtotal2 = 0;
317  for (int ikappa=0; ikappa<nKappa; ikappa++){
318  if ((NL==sectorNL[ikappa])&&(TwoSL==sectorTwoSL[ikappa])&&(IL==sectorIL[ikappa])){
319  int dimR = denBK->gCurrentDim(index+1,sectorNR[ikappa],sectorTwoSR[ikappa],sectorIR[ikappa]);
320  if (dimR>0){
321  double factor = sqrt((sectorTwoSR[ikappa]+1.0)/(TwoSL+1.0));
322  for (int l=0; l<dimL; l++){
323  for (int r=0; r<dimR; r++){
324  mem[l + dimL * (dimRtotal2 + r)] = factor * storage[kappa2index[ikappa] + l + dimL * r];
325  }
326  }
327  dimRtotal2 += dimR;
328  }
329  }
330  }
331 
332  //LQ mem --> m = dimL ; n = dimRtotal
333  int info;
334  int minofdims = min(dimL,dimRtotal);
335  double * tau = new double[minofdims];
336  double * work = new double[dimL];
337  dgelqf_(&dimL,&dimRtotal,mem,&dimL,tau,work,&dimL,&info);
338 
339  //Copy L to Lstorage
340  double * wheretoput = Lstorage->gStorage(NL,TwoSL,IL,NL,TwoSL,IL); //dimL x dimL
341 
342  for (int irow = 0; irow<dimL; irow++){
343  for (int icol = 0; icol<min(irow+1,dimRtotal); icol++){ //icol can be max. irow and max. dimRtotal-1
344  wheretoput[irow + dimL * icol] = mem[irow + dimL * icol];
345  }
346  for (int icol = min(irow+1,dimRtotal); icol<dimL; icol++){
347  wheretoput[irow + dimL * icol] = 0.0;
348  }
349  }
350 
351  //Construct Q
352  dorglq_(&minofdims,&dimRtotal,&minofdims,mem,&dimL,tau,work,&dimL,&info);
353  if (dimRtotal < dimL){ //if number of rows larger than number of cols, rest of rows zero.
354  for (int irow=dimRtotal; irow<dimL; irow++){
355  for (int icol=0; icol<dimRtotal; icol++){
356  mem[irow + dimL * icol] = 0.0;
357  }
358  }
359  }
360 
361  //Copy from mem to storage & multiply with factor !!
362  dimRtotal2 = 0;
363  for (int ikappa=0; ikappa<nKappa; ikappa++){
364  if ((NL==sectorNL[ikappa])&&(TwoSL==sectorTwoSL[ikappa])&&(IL==sectorIL[ikappa])){
365  int dimR = denBK->gCurrentDim(index+1,sectorNR[ikappa],sectorTwoSR[ikappa],sectorIR[ikappa]);
366  if (dimR>0){
367  double factor = sqrt((TwoSL+1.0)/(sectorTwoSR[ikappa]+1.0));
368  for (int l=0; l<dimL; l++){
369  for (int r=0; r<dimR; r++){
370  storage[kappa2index[ikappa] + l + dimL * r] = factor * mem[l + dimL * (r + dimRtotal2)];
371  }
372  }
373  dimRtotal2 += dimR;
374  }
375  }
376  }
377 
378  //Clear the memory
379  delete [] work;
380  delete [] tau;
381  delete [] mem;
382 
383  }
384  }
385  }
386  }
387  }
388 
389 }
390 
392 
393  //PARALLEL
394  #pragma omp parallel for schedule(dynamic)
395  for (int ikappa=0; ikappa<nKappa; ikappa++){
396  int dimL = denBK->gCurrentDim(index,sectorNL[ikappa],sectorTwoSL[ikappa],sectorIL[ikappa]);
397  int dimR = denBK->gCurrentDim(index+1,sectorNR[ikappa],sectorTwoSR[ikappa],sectorIR[ikappa]);
398  double * MxBlock = Mx->gStorage(sectorNL[ikappa],sectorTwoSL[ikappa],sectorIL[ikappa],sectorNL[ikappa],sectorTwoSL[ikappa],sectorIL[ikappa]);
399  char notrans = 'N';
400  double one = 1.0;
401  double zero = 0.0;
402  int dim = dimL*dimR;
403  double * mem = new double[dim];
404  dgemm_(&notrans,&notrans,&dimL,&dimR,&dimL,&one,MxBlock,&dimL,storage+kappa2index[ikappa],&dimL,&zero,mem,&dimL);
405  int inc = 1;
406  dcopy_(&dim,mem,&inc,storage+kappa2index[ikappa],&inc);
407  delete [] mem;
408  }
409 
410 }
411 
413 
414  //PARALLEL
415  #pragma omp parallel for schedule(dynamic)
416  for (int ikappa=0; ikappa<nKappa; ikappa++){
417  int dimL = denBK->gCurrentDim(index,sectorNL[ikappa],sectorTwoSL[ikappa],sectorIL[ikappa]);
418  int dimR = denBK->gCurrentDim(index+1,sectorNR[ikappa],sectorTwoSR[ikappa],sectorIR[ikappa]);
419  double * MxBlock = Mx->gStorage(sectorNR[ikappa],sectorTwoSR[ikappa],sectorIR[ikappa],sectorNR[ikappa],sectorTwoSR[ikappa],sectorIR[ikappa]);
420  char notrans = 'N';
421  double one = 1.0;
422  double zero = 0.0;
423  int dim = dimL*dimR;
424  double * mem = new double[dim];
425  dgemm_(&notrans,&notrans,&dimL,&dimR,&dimR,&one,storage+kappa2index[ikappa],&dimL,MxBlock,&dimR,&zero,mem,&dimL);
426  int inc = 1;
427  dcopy_(&dim,mem,&inc,storage+kappa2index[ikappa],&inc);
428  delete [] mem;
429  }
430 
431 }
432 
434 
435  bool isLeftNormal = true;
436 
437  for (int NR=denBK->gNmin(index+1); NR<=denBK->gNmax(index+1); NR++){
438  for (int TwoSR=denBK->gTwoSmin(index+1,NR); TwoSR<=denBK->gTwoSmax(index+1,NR); TwoSR+=2){
439  for (int IR=0; IR<denBK->getNumberOfIrreps(); IR++){
440  int dimR = denBK->gCurrentDim(index+1,NR,TwoSR,IR);
441  if (dimR>0){
442  double * result = new double[dimR*dimR];
443  bool firsttime = true;
444  for (int NL=NR-2; NL<=NR; NL++){
445  for (int TwoSL=TwoSR-((NR==NL+1)?1:0); TwoSL<TwoSR+2; TwoSL+=2){
446  int IL = (NR==NL+1)?(Irreps::directProd(denBK->gIrrep(index),IR)):IR;
447  int dimL = denBK->gCurrentDim(index,NL,TwoSL,IL);
448  if (dimL>0){
449  double * Block = storage + kappa2index[gKappa(NL,TwoSL,IL,NR,TwoSR,IR)];
450  char trans = 'T';
451  char notrans = 'N';
452  double one = 1.0;
453  double beta = (firsttime)?0.0:1.0;
454  dgemm_(&trans,&notrans,&dimR,&dimR,&dimL,&one,Block,&dimL,Block,&dimL,&beta,result,&dimR);
455  firsttime = false;
456  }
457  }
458  }
459  for (int cnt=0; cnt<dimR; cnt++) result[(dimR+1)*cnt] -= 1.0;
460  char norm = 'F'; //Frobenius norm
461  char uplo = 'U'; //Doesn't matter as result is fully filled
462  double TwoNorm = dlansy_(&norm,&uplo,&dimR,result,&dimR,result); //last result is work space, not referenced as norm='F';
463  if (TwoNorm>CheMPS2::TENSORT_orthoComparison) isLeftNormal = false;
464  delete [] result;
465  }
466  }
467  }
468  }
469 
470  return isLeftNormal;
471 
472 }
473 
475 
476  bool isRightNormal = true;
477 
478  for (int NL=denBK->gNmin(index); NL<=denBK->gNmax(index); NL++){
479  for (int TwoSL=denBK->gTwoSmin(index,NL); TwoSL<=denBK->gTwoSmax(index,NL); TwoSL+=2){
480  for (int IL=0; IL<denBK->getNumberOfIrreps(); IL++){
481  int dimL = denBK->gCurrentDim(index,NL,TwoSL,IL);
482  if (dimL>0){
483  double * result = new double[dimL*dimL];
484  bool firsttime = true;
485  for (int NR=NL; NR<=NL+2; NR++){
486  for (int TwoSR=TwoSL-((NR==NL+1)?1:0); TwoSR<TwoSL+2; TwoSR+=2){
487  int IR = (NR==NL+1)?(Irreps::directProd(denBK->gIrrep(index),IL)):IL;
488  int dimR = denBK->gCurrentDim(index+1,NR,TwoSR,IR);
489  if (dimR>0){
490  double * Block = storage + kappa2index[gKappa(NL,TwoSL,IL,NR,TwoSR,IR)];
491  char trans = 'T';
492  char notrans = 'N';
493  double alpha = (TwoSR+1.0)/(TwoSL+1.0);
494  double beta = (firsttime)?0.0:1.0;
495 
496  dgemm_(&notrans,&trans,&dimL,&dimL,&dimR,&alpha,Block,&dimL,Block,&dimL,&beta,result,&dimL);
497  firsttime = false;
498  }
499  }
500  }
501  for (int cnt=0; cnt<dimL; cnt++) result[(dimL+1)*cnt] -= 1.0;
502  char norm = 'F'; //Frobenius norm
503  char uplo = 'U'; //Doesn't matter as result is fully filled
504  double TwoNorm = dlansy_(&norm,&uplo,&dimL,result,&dimL,result); //last result is work space, not referenced as norm='F';
505  if (TwoNorm>CheMPS2::TENSORT_orthoComparison) isRightNormal = false;
506  delete [] result;
507  }
508  }
509  }
510  }
511 
512  return isRightNormal;
513 
514 }
515 
516 
void random()
Fill storage with random numbers 0 < val < 1.
Definition: TensorT.cpp:167
int gNmin(const int boundary) const
Get the min. possible particle number for a certain boundary.
void RightMultiply(Tensor *Mx)
Multiply at the right with a diagonal TensorOperator.
Definition: TensorT.cpp:412
int index
Index of the Tensor object. For TensorT: a site index; for other tensors: a boundary index...
Definition: Tensor.h:75
virtual ~TensorT()
Destructor.
Definition: TensorT.cpp:106
void LQ(Tensor *Lstorage)
Right-normalization.
Definition: TensorT.cpp:291
TensorT(const int site_index, const SyBookkeeper *denBK)
Constructor.
Definition: TensorT.cpp:29
int gIndex() const
Get the location index.
Definition: TensorT.cpp:161
int gCurrentDim(const int boundary, const int N, const int TwoS, const int irrep) const
Get the current virtual dimensions.
int gNKappa() const
Get the number of symmetry blocks.
Definition: TensorT.cpp:132
virtual double * gStorage()=0
Get the pointer to the storage.
int gTwoSmin(const int boundary, const int N) const
Get the minimum possible spin value for a certain boundary and particle number.
int gTwoSmax(const int boundary, const int N) const
Get the maximum possible spin value for a certain boundary and particle number.
void sBK(const SyBookkeeper *newBK)
Set the pointer to the symmetry bookkeeper.
Definition: TensorT.cpp:165
double * gStorage()
Get the pointer to the storage.
Definition: TensorT.cpp:134
static int directProd(const int Irrep1, const int Irrep2)
Get the direct product of the irreps with numbers Irrep1 and Irrep2: a bitwise XOR for psi4&#39;s convent...
Definition: Irreps.h:123
int gKappa(const int N1, const int TwoS1, const int I1, const int N2, const int TwoS2, const int I2) const
Get the index corresponding to a certain tensor block.
Definition: TensorT.cpp:136
void QR(Tensor *Rstorage)
Left-normalization.
Definition: TensorT.cpp:188
void LeftMultiply(Tensor *Mx)
Multiply at the left with a diagonal TensorOperator.
Definition: TensorT.cpp:391
const SyBookkeeper * gBK() const
Get the pointer to the symmetry bookkeeper.
Definition: TensorT.cpp:163
int getNumberOfIrreps() const
Get the total number of irreps.
int * kappa2index
kappa2index[kappa] indicates the start of tensor block kappa in storage. kappa2index[nKappa] gives th...
Definition: Tensor.h:84
int gNmax(const int boundary) const
Get the max. possible particle number for a certain boundary.
void number_operator(const double alpha, const double beta)
Apply alpha * ( number operator ) + beta to the MPS tensor.
Definition: TensorT.cpp:175
bool CheckLeftNormal() const
Check whether the TensorT is left-normal.
Definition: TensorT.cpp:433
int nKappa
Number of Tensor blocks.
Definition: Tensor.h:81
int gIrrep(const int orbital) const
Get an orbital irrep.
bool CheckRightNormal() const
Check whether the TensorT is right-normal.
Definition: TensorT.cpp:474
void Reset()
Reset the TensorT (if virtual dimensions are changed)
Definition: TensorT.cpp:125
int gKappa2index(const int kappa) const
Get the storage jump corresponding to a certain tensor block.
Definition: TensorT.cpp:151
double * storage
The actual variables. Tensor block kappa begins at storage+kappa2index[kappa] and ends at storage+kap...
Definition: Tensor.h:78