HPhi++  3.1.0
xsetmem.cpp
Go to the documentation of this file.
1 /* HPhi - Quantum Lattice Model Simulator */
2 /* Copyright (C) 2015 The University of Tokyo */
3 
4 /* This program is free software: you can redistribute it and/or modify */
5 /* it under the terms of the GNU General Public License as published by */
6 /* the Free Software Foundation, either version 3 of the License, or */
7 /* (at your option) any later version. */
8 
9 /* This program is distributed in the hope that it will be useful, */
10 /* but WITHOUT ANY WARRANTY; without even the implied warranty of */
11 /* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the */
12 /* GNU General Public License for more details. */
13 
14 /* You should have received a copy of the GNU General Public License */
15 /* along with this program. If not, see <http://www.gnu.org/licenses/>. */
16 
17 #include "Common.hpp"
18 #include "common/setmemory.hpp"
19 #include "xsetmem.hpp"
20 #include "wrapperMPI.hpp"
21 #include <iostream>
33 void setmem_HEAD
39 (
40  struct BindStruct *X
41  )
42 {
43  X->Def.CDataFileHead = (char*)malloc(D_FileNameMax*sizeof(char));
44  X->Def.CParaFileHead = (char*)malloc(D_FileNameMax*sizeof(char));
45 }
46 
52 void setmem_def
53 (
54  struct BindStruct *X,
55  struct BoostList *xBoost
56  ) {
57  X->Def.Tpow = li_1d_allocate(2 * X->Def.Nsite + 2);
58  X->Def.OrgTpow = li_1d_allocate(2 * X->Def.Nsite + 2);
59  X->Def.SiteToBit = li_1d_allocate(X->Def.Nsite + 1);
60  X->Def.LocSpn = i_1d_allocate(X->Def.Nsite);
61  X->Phys.spin_real_cor = d_1d_allocate(X->Def.Nsite * X->Def.Nsite);
62  X->Phys.charge_real_cor = d_1d_allocate(X->Def.Nsite * X->Def.Nsite);
63  X->Phys.loc_spin_z = d_1d_allocate(X->Def.Nsite * X->Def.Nsite);
64  X->Def.EDChemi = i_1d_allocate(X->Def.NInterAll + X->Def.NTransfer);
65  X->Def.EDSpinChemi = i_1d_allocate(X->Def.NInterAll + X->Def.NTransfer);
66  X->Def.EDParaChemi = d_1d_allocate(X->Def.NInterAll + X->Def.NTransfer);
67  X->Def.GeneralTransfer = i_2d_allocate(X->Def.NTransfer, 4);
68  X->Def.ParaGeneralTransfer = cd_1d_allocate(X->Def.NTransfer);
69 
70  if (X->Def.iCalcType == TimeEvolution) {
71  X->Def.EDGeneralTransfer = i_2d_allocate(X->Def.NTransfer + X->Def.NTETransferMax, 4);
72  X->Def.EDParaGeneralTransfer = cd_1d_allocate(X->Def.NTransfer + X->Def.NTETransferMax);
73  } else {
74  X->Def.EDGeneralTransfer = i_2d_allocate(X->Def.NTransfer, 4);
75  X->Def.EDParaGeneralTransfer = cd_1d_allocate(X->Def.NTransfer);
76  }
77 
78  X->Def.CoulombIntra = i_2d_allocate(X->Def.NCoulombIntra, 1);
79  X->Def.ParaCoulombIntra = d_1d_allocate(X->Def.NCoulombIntra);
80  X->Def.CoulombInter = i_2d_allocate(X->Def.NCoulombInter + X->Def.NIsingCoupling, 2);
81  X->Def.ParaCoulombInter = d_1d_allocate(X->Def.NCoulombInter + X->Def.NIsingCoupling);
82  X->Def.HundCoupling = i_2d_allocate(X->Def.NHundCoupling + X->Def.NIsingCoupling, 2);
83  X->Def.ParaHundCoupling = d_1d_allocate(X->Def.NHundCoupling + X->Def.NIsingCoupling);
84  X->Def.PairHopping = i_2d_allocate(X->Def.NPairHopping, 2);
85  X->Def.ParaPairHopping = d_1d_allocate(X->Def.NPairHopping);
86  X->Def.ExchangeCoupling = i_2d_allocate(X->Def.NExchangeCoupling, 2);
87  X->Def.ParaExchangeCoupling = d_1d_allocate(X->Def.NExchangeCoupling);
88  X->Def.PairLiftCoupling = i_2d_allocate(X->Def.NPairLiftCoupling, 2);
89  X->Def.ParaPairLiftCoupling = d_1d_allocate(X->Def.NPairLiftCoupling);
90 
91  X->Def.InterAll = i_2d_allocate(X->Def.NInterAll, 8);
92  X->Def.ParaInterAll = cd_1d_allocate(X->Def.NInterAll);
93 
94  X->Def.CisAjt = i_2d_allocate(X->Def.NCisAjt, 4);
95  X->Def.CisAjtCkuAlvDC = i_2d_allocate(X->Def.NCisAjtCkuAlvDC, 8);
96 
98  X->Def.SingleExcitationOperator = (int***)malloc(sizeof(int**)*X->Def.NNSingleExcitationOperator);
99  X->Def.ParaSingleExcitationOperator = (std::complex<double>**)malloc(
100  sizeof(std::complex<double>*)*X->Def.NNSingleExcitationOperator);
102  X->Def.PairExcitationOperator = (int***)malloc(sizeof(int**)*X->Def.NNPairExcitationOperator);
103  X->Def.ParaPairExcitationOperator = (std::complex<double>**)malloc(
104  sizeof(std::complex<double>*)*X->Def.NNPairExcitationOperator);
105 
106  X->Def.ParaLaser = d_1d_allocate(X->Def.NLaser);
107 
108  xBoost->list_6spin_star = i_2d_allocate(xBoost->R0 * xBoost->num_pivot, 7);
109  xBoost->list_6spin_pair = i_3d_allocate(xBoost->R0 * xBoost->num_pivot, 7, 15);
110  xBoost->arrayJ = cd_3d_allocate(xBoost->NumarrayJ, 3, 3);
111 
112  int NInterAllSet;
113  NInterAllSet = (X->Def.iCalcType == TimeEvolution) ? X->Def.NInterAll + X->Def.NTEInterAllMax : X->Def.NInterAll;
114  X->Def.InterAll_OffDiagonal = i_2d_allocate(NInterAllSet, 8);
115  X->Def.ParaInterAll_OffDiagonal = cd_1d_allocate(NInterAllSet);
116  X->Def.InterAll_Diagonal = i_2d_allocate(NInterAllSet, 4);
117  X->Def.ParaInterAll_Diagonal = d_1d_allocate(NInterAllSet);
118 
119  if (X->Def.iCalcType == TimeEvolution) {
120  X->Def.TETime = d_1d_allocate(X->Def.NTETimeSteps);
121  //Time-dependent Transfer
122  X->Def.NTETransfer = i_1d_allocate(X->Def.NTETimeSteps);
123  X->Def.NTETransferDiagonal = i_1d_allocate(X->Def.NTETimeSteps);
124  X->Def.TETransfer = i_3d_allocate(X->Def.NTETimeSteps, X->Def.NTETransferMax, 4);
125  X->Def.TETransferDiagonal = i_3d_allocate(X->Def.NTETimeSteps, X->Def.NTETransferMax, 2);
126  X->Def.ParaTETransfer = cd_2d_allocate(X->Def.NTETimeSteps, X->Def.NTETransferMax);
127  X->Def.ParaTETransferDiagonal = d_2d_allocate(X->Def.NTETimeSteps, X->Def.NTETransferMax);
128  //Time-dependent InterAll
129  X->Def.NTEInterAll = i_1d_allocate(X->Def.NTETimeSteps);
130  X->Def.NTEInterAllDiagonal = i_1d_allocate(X->Def.NTETimeSteps);
131  X->Def.TEInterAll = i_3d_allocate(X->Def.NTETimeSteps, X->Def.NTEInterAllMax, 8);
132  X->Def.TEInterAllDiagonal = i_3d_allocate(X->Def.NTETimeSteps, X->Def.NTEInterAllMax, 4);
133  X->Def.ParaTEInterAll = cd_2d_allocate(X->Def.NTETimeSteps, X->Def.NTEInterAllMax);
134  X->Def.ParaTEInterAllDiagonal = d_2d_allocate(X->Def.NTETimeSteps, X->Def.NTEInterAllMax);
135  X->Def.NTEInterAllOffDiagonal = i_1d_allocate(X->Def.NTETimeSteps);
136  X->Def.TEInterAllOffDiagonal = i_3d_allocate(X->Def.NTETimeSteps, X->Def.NTEInterAllMax, 8);
137  X->Def.ParaTEInterAllOffDiagonal = cd_2d_allocate(X->Def.NTETimeSteps, X->Def.NTEInterAllMax);
138  //Time-dependent Chemi generated by InterAll diagonal components
139  X->Def.NTEChemi = i_1d_allocate(X->Def.NTETimeSteps);
140  X->Def.TEChemi = i_2d_allocate(X->Def.NTETimeSteps, X->Def.NTEInterAllMax);
141  X->Def.SpinTEChemi = i_2d_allocate(X->Def.NTETimeSteps, X->Def.NTEInterAllMax);
142  X->Def.ParaTEChemi = d_2d_allocate(X->Def.NTETimeSteps, X->Def.NTEInterAllMax);
143  }
144 }
145 
146 
153 int setmem_large
154 (
155  struct BindStruct *X
156 ) {
157  int nstate;
158 
159  if (GetlistSize(X) == TRUE) {
160  list_1 = li_1d_allocate(X->Check.idim_max + 1);
161  list_2_1 = li_1d_allocate(X->Large.SizeOflist_2_1);
162  list_2_2 = li_1d_allocate(X->Large.SizeOflist_2_2);
163  if (list_1 == NULL
164  || list_2_1 == NULL
165  || list_2_2 == NULL
166  ) {
167  return -1;
168  }
169  }
170 
171  list_Diagonal = d_1d_allocate(X->Check.idim_max + 1);
172 
173  if (X->Def.iCalcType == FullDiag) {
174  nstate = X->Check.idim_max;
175  }
176  else if (X->Def.iCalcType == CG) {
177  nstate = X->Def.k_exct;
178  }
179  else if (X->Def.iCalcType == TPQCalc) {
180  nstate = NumAve;
181  }
182  else {
183  nstate = 1;
184  }
185  v0 = cd_2d_allocate(X->Check.idim_max + 1, nstate);
186  v1 = cd_2d_allocate(X->Check.idim_max + 1, nstate);
187 #ifdef __MPI
188  long int MAXidim_max;
189  MAXidim_max = MaxMPI_li(X->Check.idim_max);
190  if (GetlistSize(X) == TRUE) list_1buf = li_1d_allocate(MAXidim_max + 1);
191  v1buf = cd_2d_allocate(MAXidim_max + 1, nstate);
192 #else
193  if (X->Def.iCalcType == CG) v1buf = cd_2d_allocate(X->Check.idim_max + 1, nstate);
194 #endif // MPI
195 
196  X->Phys.num_down = d_1d_allocate(nstate);
197  X->Phys.num_up = d_1d_allocate(nstate);
198  X->Phys.num = d_1d_allocate(nstate);
199  X->Phys.num2 = d_1d_allocate(nstate);
200  X->Phys.energy = d_1d_allocate(nstate);
201  X->Phys.var = d_1d_allocate(nstate);
202  X->Phys.doublon = d_1d_allocate(nstate);
203  X->Phys.doublon2 = d_1d_allocate(nstate);
204  X->Phys.Sz = d_1d_allocate(nstate);
205  X->Phys.Sz2 = d_1d_allocate(nstate);
206  X->Phys.s2 = d_1d_allocate(nstate);
207 
208  fprintf(stdoutMPI, "%s", "\n###### LARGE ALLOCATE FINISH ! ######\n\n");
209  return 0;
210 }
221  struct BindStruct *X
222 ) {
223  switch (X->Def.iCalcModel) {
224  case Spin:
225  case Hubbard:
226  case HubbardNConserved:
227  case Kondo:
228  case KondoGC:
229  if (X->Def.iFlgGeneralSpin == FALSE) {
230  if (X->Def.iCalcModel == Spin && X->Def.Nsite % 2 == 1) {
231  X->Large.SizeOflist_2_1 = X->Check.sdim * 2 + 2;
232  }
233  else {
234  X->Large.SizeOflist_2_1 = X->Check.sdim + 2;
235  }
236  X->Large.SizeOflist_2_2 = X->Check.sdim + 2;
237  X->Large.SizeOflistjb = X->Check.sdim + 2;
238  }
239  else {//for spin-canonical general spin
240  X->Large.SizeOflist_2_1 = X->Check.sdim + 2;
241  X->Large.SizeOflist_2_2 =
242  X->Def.Tpow[X->Def.Nsite - 1] * X->Def.SiteToBit[X->Def.Nsite - 1] / X->Check.sdim + 2;
243  X->Large.SizeOflistjb =
244  X->Def.Tpow[X->Def.Nsite - 1] * X->Def.SiteToBit[X->Def.Nsite - 1] / X->Check.sdim + 2;
245  }
246  break;
247  default:
248  X->Large.SizeOflistjb = 1;
249  return FALSE;
250  }
251  return TRUE;
252 }
double * doublon
Expectation value of the Doublon.
Definition: struct.hpp:357
int * NTETransferDiagonal
Definition: struct.hpp:260
struct DefineList Def
Definision of system (Hamiltonian) etc.
Definition: struct.hpp:395
std::complex< double > * EDParaGeneralTransfer
Value of general transfer integrals by a def file. malloc in setmem_def(). Data Format [DefineList::N...
Definition: struct.hpp:116
double * num_down
Expectation value of the number of down-spin electtrons.
Definition: struct.hpp:364
std::complex< double > * ParaInterAll
[DefineList::NInterAll] Coupling constant of inter-all term. malloc in setmem_def().
Definition: struct.hpp:166
double * var
Expectation value of the Energy variance.
Definition: struct.hpp:367
int *** TEInterAll
Definition: struct.hpp:280
FILE * stdoutMPI
File pointer to the standard output defined in InitializeMPI()
Definition: global.cpp:75
int NCoulombInter
Number of off-site Coulomb interaction.
Definition: struct.hpp:127
long int * list_2_1
Definition: global.cpp:27
std::complex< double > ** ParaTETransfer
Definition: struct.hpp:266
int * EDSpinChemi
[DefineList::Nsite]
Definition: struct.hpp:99
long int * list_2_2
Definition: global.cpp:28
int ** list_6spin_star
Definition: struct.hpp:388
int ** SpinTEChemi
Definition: struct.hpp:297
std::complex< double > ** ParaTEInterAll
Definition: struct.hpp:288
long int * OrgTpow
[2 * DefineList::NsiteMPI] malloc in setmem_def().
Definition: struct.hpp:92
int ** ExchangeCoupling
[DefineList::NExchangeCoupling][2] Index of exchange term. malloc in setmem_def().
Definition: struct.hpp:146
std::complex< double > ** v1buf
Definition: global.cpp:22
int NumAve
Definition: global.cpp:43
double * ParaInterAll_Diagonal
[DefineList::NInterAll_Diagonal] Coupling constant of diagonal inter-all term. malloc in setmem_def()...
Definition: struct.hpp:168
int NPairLiftCoupling
Number of pair-lift term.
Definition: struct.hpp:153
std::complex< double > ** ParaSingleExcitationOperator
[DefineList::NSingleExcitationOperator] Coefficient of single excitaion operator for spectrum...
Definition: struct.hpp:184
std::complex< double > ** v0
Definition: global.cpp:20
int * NSingleExcitationOperator
Number of single excitaion operator for spectrum.
Definition: struct.hpp:183
int *** TETransferDiagonal
Definition: struct.hpp:264
double * s2
Expectation value of the square of the total S.
Definition: struct.hpp:365
int Nsite
Number of sites in the INTRA process region.
Definition: struct.hpp:56
int NumarrayJ
Definition: struct.hpp:385
double * TETime
Definition: struct.hpp:249
int ** TEChemi
Definition: struct.hpp:295
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.hpp:397
double * num2
Expectation value of the quare of the number of electrons.
Definition: struct.hpp:360
int * LocSpn
[DefineList::NLocSpn] Flag (and size) of the local spin. malloc in setmem_def().
Definition: struct.hpp:82
int ** CisAjtCkuAlvDC
[DefineList::NCisAjtCkuAlvDC][4] Indices of two-body correlation function. malloc in setmem_def()...
Definition: struct.hpp:177
int ** EDGeneralTransfer
Index of transfer integrals for calculation. malloc in setmem_def(). Data Format [DefineList::NTransf...
Definition: struct.hpp:110
void setmem_HEAD(struct BindStruct *X)
Set size of memories headers of output files.
Definition: xsetmem.cpp:39
struct PhysList Phys
Physical quantities.
Definition: struct.hpp:398
double * Sz2
Expectation value of the Square of total Sz.
Definition: struct.hpp:362
long int SizeOflist_2_1
Size of list_2_1.
Definition: struct.hpp:319
long int MaxMPI_li(long int idim)
MPI wrapper function to obtain maximum unsigned long integer across processes.
Definition: wrapperMPI.cpp:171
int * NPairExcitationOperator
Number of pair excitaion operator for spectrum.
Definition: struct.hpp:190
int *** list_6spin_pair
Definition: struct.hpp:389
std::complex< double > ** v1
Definition: global.cpp:21
int * NTEChemi
Definition: struct.hpp:296
double * ParaLaser
Definition: struct.hpp:253
double ** ParaTEChemi
Definition: struct.hpp:298
double * Sz
Expectation value of the Total Sz.
Definition: struct.hpp:361
int NHundCoupling
Number of Hund coupling.
Definition: struct.hpp:133
int NIsingCoupling
Number of Ising term.
Definition: struct.hpp:151
std::complex< double > ** ParaPairExcitationOperator
[DefineList::NPairExcitationOperator] Coefficient of pair excitaion operator for spectrum. malloc in setmem_def().
Definition: struct.hpp:191
int NInterAll
Total Number of Interacted quartet.
Definition: struct.hpp:163
int ** InterAll_OffDiagonal
[DefineList::NinterAll_OffDiagonal][8] Interacted quartet
Definition: struct.hpp:161
int ** InterAll_Diagonal
[DefineList::NinterAll_Diagonal][4] Interacted quartet
Definition: struct.hpp:162
double * loc_spin_z
Malloc, but Not used ???
Definition: struct.hpp:372
int ** GeneralTransfer
Index of transfer integrals obtained by a def file. malloc in setmem_def(). Data Format [DefineList::...
Definition: struct.hpp:106
int * NTEInterAllOffDiagonal
Definition: struct.hpp:275
int ** PairHopping
[DefineList::NPairHopping][2] Index of pair-hopping. malloc in setmem_def().
Definition: struct.hpp:140
long int * list_1buf
Definition: global.cpp:26
double ** ParaTETransferDiagonal
Definition: struct.hpp:268
int NNSingleExcitationOperator
Number of single excitaion operator for spectrum.
Definition: struct.hpp:182
int ** PairLiftCoupling
[DefineList::NPairHopping][2] Index of pair-lift term. malloc in setmem_def().
Definition: struct.hpp:154
int ** HundCoupling
[DefineList::NHundCoupling][2] Index of Hund coupling. malloc in setmem_def().
Definition: struct.hpp:134
int * NTEInterAll
Definition: struct.hpp:273
int NCisAjt
Number of indices of two-body correlation function.
Definition: struct.hpp:175
double * num
Expectation value of the Number of electrons.
Definition: struct.hpp:359
Bind.
Definition: struct.hpp:394
double * ParaExchangeCoupling
[DefineList::NExchangeCoupling] Coupling constant of exchange term. malloc in setmem_def().
Definition: struct.hpp:148
double * ParaCoulombInter
[DefineList::NCoulombInter]Coupling constant of off-site Coulomb interaction. malloc in setmem_def()...
Definition: struct.hpp:130
int NLaser
Definition: struct.hpp:252
double * ParaPairLiftCoupling
[DefineList::NPairHopping] Coupling constant of pair-lift term. malloc in setmem_def().
Definition: struct.hpp:156
int * NTETransfer
Definition: struct.hpp:258
int ** CoulombInter
Definition: struct.hpp:128
void setmem_def(struct BindStruct *X, struct BoostList *xBoost)
Set size of memories for Def and Phys in BindStruct.
Definition: xsetmem.cpp:53
double * EDParaChemi
[DefineList::Nsite] On-site potential parameter. malloc in setmem_def().
Definition: struct.hpp:100
double * energy
Expectation value of the total energy.
Definition: struct.hpp:356
int NExchangeCoupling
Number of exchange term.
Definition: struct.hpp:145
int * NTEInterAllDiagonal
Definition: struct.hpp:278
int *** TETransfer
Definition: struct.hpp:262
std::complex< double > * ParaInterAll_OffDiagonal
[DefineList::NInterAll_OffDiagonal] Coupling constant of off-diagonal inter-all term. malloc in setmem_def().
Definition: struct.hpp:170
char * CParaFileHead
Read from Calcmod in readdef.h. It is not used. Just for the compatibility to mVMC.
Definition: struct.hpp:44
int NTransfer
Number of transfer integrals obtained by a def file.
Definition: struct.hpp:104
int iFlgGeneralSpin
Flag for the general (Sz/=1/2) spin.
Definition: struct.hpp:86
int NTETransferMax
Definition: struct.hpp:257
double * list_Diagonal
Definition: global.cpp:24
long int * SiteToBit
[DefineList::NsiteMPI] Similar to DefineList::Tpow. For general spin.
Definition: struct.hpp:94
double * doublon2
Expectation value of the Square of doublon.
Definition: struct.hpp:358
long int * Tpow
[2 * DefineList::NsiteMPI] malloc in setmem_def().
Definition: struct.hpp:90
int ** CisAjt
[DefineList::NCisAjt][4] Indices of one-body correlation function. malloc in setmem_def().
Definition: struct.hpp:174
double ** ParaTEInterAllDiagonal
Definition: struct.hpp:293
double * ParaPairHopping
[DefineList::NPairHopping] Coupling constant of pair-hopping term. malloc in setmem_def().
Definition: struct.hpp:142
int ** CoulombIntra
Definition: struct.hpp:122
long int num_pivot
Definition: struct.hpp:383
double * charge_real_cor
Malloc, but Not used ???
Definition: struct.hpp:371
For Boost.
Definition: struct.hpp:379
double * ParaHundCoupling
[DefineList::NHundCoupling] Hund coupling constant. malloc in setmem_def().
Definition: struct.hpp:136
std::complex< double > * ParaGeneralTransfer
Value of general transfer integrals by a def file. malloc in setmem_def(). Data Format [DefineList::N...
Definition: struct.hpp:113
int NPairHopping
Number of pair-hopping term.
Definition: struct.hpp:139
int NCoulombIntra
Definition: struct.hpp:121
int NTETimeSteps
Definition: struct.hpp:248
double * spin_real_cor
Malloc, but Not used ???
Definition: struct.hpp:370
int setmem_large(struct BindStruct *X)
Set size of memories for vectors(vg, v0, v1, v2, vec, alpha, beta), lists (list_1, list_2_1, list_2_2, list_Diagonal) and Phys(BindStruct.PhysList) struct in the case of Full Diag mode.
Definition: xsetmem.cpp:154
double * num_up
Expectation value of the number of up-spin electtrons.
Definition: struct.hpp:363
int iCalcModel
Switch for model. 0:Hubbard, 1:Spin, 2:Kondo, 3:HubbardGC, 4:SpinGC, 5:KondoGC, 6:HubbardNConserved.
Definition: struct.hpp:200
long int R0
Definition: struct.hpp:381
int *** TEInterAllOffDiagonal
Definition: struct.hpp:283
double * ParaCoulombIntra
Definition: struct.hpp:124
int GetlistSize(struct BindStruct *X)
Set size of lists for the canonical ensemble.
Definition: xsetmem.cpp:220
int *** TEInterAllDiagonal
Definition: struct.hpp:286
int k_exct
Read from Calcmod in readdef.h.
Definition: struct.hpp:47
int ** InterAll
[DefineList::NinterAll][8] Interacted quartet
Definition: struct.hpp:160
long int SizeOflistjb
Used for computing Sz.
Definition: struct.hpp:321
int * EDChemi
[DefineList::Nsite] Chemical potential. malloc in setmem_def().
Definition: struct.hpp:98
std::complex< double > *** arrayJ
Definition: struct.hpp:386
struct CheckList Check
Size of the Hilbert space.
Definition: struct.hpp:396
char * CDataFileHead
Read from Calcmod in readdef.h. Header of output file such as Green&#39;s function.
Definition: struct.hpp:42
long int sdim
Dimension for Ogata-Lin ???
Definition: struct.hpp:309
int NNPairExcitationOperator
Number of pair excitaion operator for spectrum.
Definition: struct.hpp:189
int NCisAjtCkuAlvDC
Number of indices of two-body correlation function.
Definition: struct.hpp:178
int *** SingleExcitationOperator
[DefineList::NSingleExcitationOperator][3] Indices of single excitaion operator for spectrum...
Definition: struct.hpp:180
long int * list_1
Definition: global.cpp:25
long int idim_max
The dimension of the Hilbert space of this process.
Definition: struct.hpp:305
long int SizeOflist_2_2
Size of list_2_2.
Definition: struct.hpp:320
int *** PairExcitationOperator
[DefineList::NPairExcitationOperator][5] Indices of pair excitaion operator for spectrum. malloc in setmem_def().
Definition: struct.hpp:187
int NTEInterAllMax
Definition: struct.hpp:272
std::complex< double > ** ParaTEInterAllOffDiagonal
Definition: struct.hpp:290
int iCalcType
Switch for calculation type. 0:Lanczos, 1:TPQCalc, 2:FullDiag.
Definition: struct.hpp:194