HPhi++  3.1.0
check.cpp File Reference

File for giving a function of calculating size of Hilbert space. More...

#include "bitcalc.hpp"
#include "sz.hpp"
#include "FileIO.hpp"
#include "common/setmemory.hpp"
#include "check.hpp"
#include "wrapperMPI.hpp"
#include "CheckMPI.hpp"

Go to the source code of this file.

Functions

int check (struct BindStruct *X)
 A program to check size of dimension for Hilbert-space. More...
 

Detailed Description

File for giving a function of calculating size of Hilbert space.

Version
0.1, 0.2
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)

Definition in file check.cpp.

Function Documentation

◆ check()

int check ( struct BindStruct X)

A program to check size of dimension for Hilbert-space.

Parameters
[in,out]XCommon data set used in HPhi.
Return values
TRUEnormal termination
FALSEabnormal termination
MPIFALSECheckMPI abnormal termination
Version
0.2

add function of calculating Hilbert space for canonical ensemble.

Version
0.1
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)

Definition at line 51 of file check.cpp.

References Binomial(), BindStruct::Check, CheckMPI(), CheckMPI_Summary(), childfopenMPI(), BindStruct::Def, GetLocal2Sz(), GetSplitBitForGeneralSpin(), DefineList::iCalcModel, DefineList::iCalcType, CheckList::idim_max, DefineList::iFlgCalcSpec, DefineList::iFlgGeneralSpin, DefineList::iFlgScaLAPACK, DefineList::k_exct, CheckList::max_mem, MaxMPI_d(), MaxMPI_li(), DefineList::Ndown, DefineList::Ne, DefineList::NLocSpn, DefineList::Nsite, DefineList::NsiteMPI, NumAve, DefineList::Nup, CheckList::sdim, DefineList::SiteToBit, stdoutMPI, DefineList::Total2Sz, DefineList::Total2SzMPI, and DefineList::Tpow.

Referenced by main(), and MakeExcitedList().

51  {
52 
53  FILE *fp;
54  long int i,tmp_sdim;
55  int NLocSpn, NCond;
56  int Nup, Ndown;
57  long int u_tmp;
58  long int tmp;
59  long int Ns,comb_1,comb_2,comb_3,comb_sum, comb_up, comb_down;
60  int u_loc;
61  long int **comb;
62  long int idimmax=0;
63  long int idim=0;
64  long int isite=0;
65  int tmp_sz=0;
66  int iMinup=0;
67  if(X->Def.iCalcModel ==Spin ||X->Def.iCalcModel ==SpinGC )
68  {
69  X->Def.Ne=X->Def.Nup;
70  }
71 
72  int iAllup=X->Def.Ne;
73 
74  if(X->Def.iFlgScaLAPACK == 0) {
75  /*
76  Set Site number per MPI process
77  */
78  if (CheckMPI(X) != TRUE) {
79  return MPIFALSE;
80  }
81  }
82  else{
83  X->Def.NsiteMPI = X->Def.Nsite;
84  X->Def.Total2SzMPI = X->Def.Total2Sz;
85  }
86 
87  Ns = X->Def.Nsite;
88 
89  comb = li_2d_allocate(Ns+1,Ns+1);
90 
91  //idim_max
92  switch(X->Def.iCalcModel){
93  case HubbardGC:
94  //comb_sum = 2^(2*Ns)=4^Ns
95  comb_sum = 1;
96  for(i=0;i<2*X->Def.Nsite;i++){
97  comb_sum= 2*comb_sum;
98  }
99  break;
100  case SpinGC:
101  //comb_sum = 2^(Ns)
102  comb_sum = 1;
103  if(X->Def.iFlgGeneralSpin ==FALSE){
104  for(i=0;i<X->Def.Nsite;i++){
105  comb_sum= 2*comb_sum;
106  }
107  }
108  else{
109  for(i=0; i<X->Def.Nsite;i++){
110  comb_sum=comb_sum*X->Def.SiteToBit[i];
111  }
112  }
113  break;
114 
115  case Hubbard:
116  comb_up= Binomial(Ns, X->Def.Nup, comb, Ns);
117  comb_down= Binomial(Ns, X->Def.Ndown, comb, Ns);
118  comb_sum=comb_up*comb_down;
119  break;
120 
121  case HubbardNConserved:
122  comb_sum=0;
123  if(X->Def.Ne > X->Def.Nsite){
124  iMinup = X->Def.Ne-X->Def.Nsite;
125  iAllup = X->Def.Nsite;
126  }
127 
128  for(i=iMinup; i<= iAllup; i++){
129  comb_up= Binomial(Ns, i, comb, Ns);
130  comb_down= Binomial(Ns, X->Def.Ne-i, comb, Ns);
131  comb_sum +=comb_up*comb_down;
132  }
133  break;
134 
135  case Kondo:
136  Nup = X->Def.Nup;
137  Ndown = X->Def.Ndown;
138  NCond = X->Def.Nsite-X->Def.NLocSpn;
139  NLocSpn = X->Def.NLocSpn;
140  comb_sum = 0;
141  for(u_loc=0;u_loc<=X->Def.Nup;u_loc++){
142  comb_1 = Binomial(NLocSpn,u_loc,comb,Ns);
143  comb_2 = Binomial(NCond,Nup-u_loc,comb,Ns);
144  comb_3 = Binomial(NCond,Ndown+u_loc-NLocSpn,comb,Ns);
145  comb_sum += comb_1*comb_2*comb_3;
146  }
147  break;
148  case KondoGC:
149  comb_sum = 1;
150  NCond = X->Def.Nsite-X->Def.NLocSpn;
151  NLocSpn = X->Def.NLocSpn;
152  //4^Nc*2^Ns
153  for(u_loc=0;u_loc <(2*NCond+NLocSpn); u_loc++){
154  comb_sum= 2*comb_sum;
155  }
156  break;
157  case Spin:
158 
159  if(X->Def.iFlgGeneralSpin ==FALSE){
160  if(X->Def.Nup+X->Def.Ndown != X->Def.Nsite){
161  fprintf(stderr, " 2Sz is incorrect.\n");
162  return FALSE;
163  }
164  //comb_sum= Binomial(Ns, X->Def.Ne, comb, Ns);
165  comb_sum= Binomial(Ns, X->Def.Nup, comb, Ns);
166  }
167  else{
168  idimmax = 1;
169  X->Def.Tpow[0]=idimmax;
170  for(isite=0; isite<X->Def.Nsite;isite++){
171  idimmax=idimmax*X->Def.SiteToBit[isite];
172  X->Def.Tpow[isite+1]=idimmax;
173  }
174  comb_sum=0;
175 #pragma omp parallel for default(none) reduction(+:comb_sum) private(tmp_sz, isite) firstprivate(idimmax, X)
176  for(idim=0; idim<idimmax; idim++){
177  tmp_sz=0;
178  for(isite=0; isite<X->Def.Nsite;isite++){
179  tmp_sz += GetLocal2Sz(isite+1,idim, X->Def.SiteToBit, X->Def.Tpow );
180  }
181  if(tmp_sz == X->Def.Total2Sz){
182  comb_sum +=1;
183  }
184  }
185 
186  }
187 
188  break;
189  default:
190  fprintf(stderr, "Error: CalcModel %d is incorrect.\n", X->Def.iCalcModel);
191  free_li_2d_allocate(comb);
192  return FALSE;
193  }
194 
195  //fprintf(stdoutMPI, "Debug: comb_sum= %ld \n",comb_sum);
196 
197  X->Check.idim_max = comb_sum;
198  switch(X->Def.iCalcType) {
199  case CG:
200  switch (X->Def.iCalcModel) {
201  case Hubbard:
202  case HubbardNConserved:
203  case Kondo:
204  case KondoGC:
205  case Spin:
206  X->Check.max_mem = (7 * X->Def.k_exct + 1.5) * X->Check.idim_max * 16.0 / (pow(10, 9));
207  break;
208  case HubbardGC:
209  case SpinGC:
210  X->Check.max_mem = (7 * X->Def.k_exct + 1.0) * X->Check.idim_max * 16.0 / (pow(10, 9));
211  break;
212  }
213  break;
214  case TPQCalc:
215  switch (X->Def.iCalcModel) {
216  case Hubbard:
217  case HubbardNConserved:
218  case Kondo:
219  case KondoGC:
220  case Spin:
221  if (X->Def.iFlgCalcSpec != CALCSPEC_NOT) {
222  X->Check.max_mem = NumAve * 3 * X->Check.idim_max * 16.0 / (pow(10, 9));
223  } else {
224  X->Check.max_mem = 4.5 * X->Check.idim_max * 16.0 / (pow(10, 9));
225  }
226  break;
227  case HubbardGC:
228  case SpinGC:
229  if (X->Def.iFlgCalcSpec != CALCSPEC_NOT) {
230  X->Check.max_mem = NumAve * 3 * X->Check.idim_max * 16.0 / (pow(10, 9));
231  } else {
232  X->Check.max_mem = 3.5 * X->Check.idim_max * 16.0 / (pow(10, 9));
233  }
234  break;
235  }
236  break;
237  case FullDiag:
238  X->Check.max_mem = X->Check.idim_max * 8.0 * X->Check.idim_max * 8.0 / (pow(10, 9));
239  break;
240  case TimeEvolution:
241  X->Check.max_mem = (4 + 2 + 1) * X->Check.idim_max * 16.0 / (pow(10, 9));
242  break;
243  default:
244  return FALSE;
245  //break;
246  }
247 
248  //fprintf(stdoutMPI, " MAX DIMENSION idim_max=%ld \n",X->Check.idim_max);
249  //fprintf(stdoutMPI, " APPROXIMATE REQUIRED MEMORY max_mem=%lf GB \n",X->Check.max_mem);
250  long int li_dim_max=MaxMPI_li(X->Check.idim_max);
251  fprintf(stdoutMPI, " MAX DIMENSION idim_max=%ld \n",li_dim_max);
252  double dmax_mem=MaxMPI_d(X->Check.max_mem);
253  fprintf(stdoutMPI, " APPROXIMATE REQUIRED MEMORY max_mem=%lf GB \n",dmax_mem);
254  if(childfopenMPI("CHECK_Memory.dat","w", &fp)!=0){
255  free_li_2d_allocate(comb);
256  return FALSE;
257  }
258  fprintf(fp," MAX DIMENSION idim_max=%ld \n", li_dim_max);
259  fprintf(fp," APPROXIMATE REQUIRED MEMORY max_mem=%lf GB \n", dmax_mem);
260 
261 
262  /*
263  fprintf(fp," MAX DIMENSION idim_max=%ld \n",X->Check.idim_max);
264  fprintf(fp," APPROXIMATE REQUIRED MEMORY max_mem=%lf GB \n",X->Check.max_mem);
265  */
266  fclose(fp);
267 
268  //sdim
269  tmp=1;
270  tmp_sdim=1;
271 
272  switch(X->Def.iCalcModel){
273  case HubbardGC:
274  case KondoGC:
275  case HubbardNConserved:
276  case Hubbard:
277  case Kondo:
278  while(tmp <= X->Def.Nsite){
279  tmp_sdim=tmp_sdim*2;
280  tmp+=1;
281  }
282  break;
283  case Spin:
284  case SpinGC:
285  if(X->Def.iFlgGeneralSpin==FALSE){
286  while(tmp <= X->Def.Nsite/2){
287  tmp_sdim=tmp_sdim*2;
288  tmp+=1;
289  }
290  }
291  else{
292  GetSplitBitForGeneralSpin(X->Def.Nsite, &tmp_sdim, X->Def.SiteToBit);
293  }
294  break;
295  default:
296  fprintf(stdoutMPI, "Error: CalcModel %d is incorrect.\n", X->Def.iCalcModel);
297  free_li_2d_allocate(comb);
298  return FALSE;
299  }
300  X->Check.sdim=tmp_sdim;
301 
302  if(childfopenMPI("CHECK_Sdim.dat","w", &fp)!=0){
303  free_li_2d_allocate(comb);
304  return FALSE;
305  }
306 
307  switch(X->Def.iCalcModel){
308  case HubbardGC:
309  case KondoGC:
310  case HubbardNConserved:
311  case Hubbard:
312  case Kondo:
313  //fprintf(stdoutMPI, "sdim=%ld =2^%d\n",X->Check.sdim,X->Def.Nsite);
314  fprintf(fp,"sdim=%ld =2^%d\n",X->Check.sdim,X->Def.Nsite);
315  break;
316  case Spin:
317  case SpinGC:
318  if(X->Def.iFlgGeneralSpin==FALSE){
319  //fprintf(stdoutMPI, "sdim=%ld =2^%d\n",X->Check.sdim,X->Def.Nsite/2);
320  fprintf(fp,"sdim=%ld =2^%d\n",X->Check.sdim,X->Def.Nsite/2);
321  }
322  break;
323  default:
324  break;
325  }
326 
327  free_li_2d_allocate(comb);
328 
329  u_tmp=1;
330  X->Def.Tpow[0]=u_tmp;
331  switch(X->Def.iCalcModel){
332  case HubbardGC:
333  case KondoGC:
334  for(i=1;i<=2*X->Def.Nsite;i++){
335  u_tmp=u_tmp*2;
336  X->Def.Tpow[i]=u_tmp;
337  fprintf(fp,"%ld %ld \n",i,u_tmp);
338  }
339  break;
340  case HubbardNConserved:
341  case Hubbard:
342  case Kondo:
343  for(i=1;i<=2*X->Def.Nsite-1;i++){
344  u_tmp=u_tmp*2;
345  X->Def.Tpow[i]=u_tmp;
346  fprintf(fp,"%ld %ld \n",i,u_tmp);
347  }
348  break;
349  case SpinGC:
350  if(X->Def.iFlgGeneralSpin==FALSE){
351  for(i=1;i<=X->Def.Nsite;i++){
352  u_tmp=u_tmp*2;
353  X->Def.Tpow[i]=u_tmp;
354  fprintf(fp,"%ld %ld \n",i,u_tmp);
355  }
356  }
357  else{
358  X->Def.Tpow[0]=u_tmp;
359  fprintf(fp,"%d %ld \n", 0, u_tmp);
360  for(i=1;i<X->Def.Nsite;i++){
361  u_tmp=u_tmp*X->Def.SiteToBit[i-1];
362  X->Def.Tpow[i]=u_tmp;
363  fprintf(fp,"%ld %ld \n",i,u_tmp);
364  }
365  }
366  break;
367  case Spin:
368  if(X->Def.iFlgGeneralSpin==FALSE){
369  for(i=1;i<=X->Def.Nsite-1;i++){
370  u_tmp=u_tmp*2;
371  X->Def.Tpow[i]=u_tmp;
372  fprintf(fp,"%ld %ld \n",i,u_tmp);
373  }
374  }
375  else{
376  for(i=0;i<X->Def.Nsite;i++){
377  fprintf(fp,"%ld %ld \n",i,X->Def.Tpow[i]);
378  }
379  }
380  break;
381  default:
382  fprintf(stdoutMPI, "Error: CalcModel %d is incorrect.\n", X->Def.iCalcModel);
383  free_li_2d_allocate(comb);
384  return FALSE;
385  }
386  fclose(fp);
387  /*
388  Print MPI-site information and Modify Tpow
389  in the inter process region.
390  */
391  CheckMPI_Summary(X);
392 
393  return TRUE;
394 }
long int Binomial(int n, int k, long int **comb, int Nsite)
Definition: sz.cpp:178
int Nup
Number of spin-up electrons in this process.
Definition: struct.hpp:58
struct DefineList Def
Definision of system (Hamiltonian) etc.
Definition: struct.hpp:395
FILE * stdoutMPI
File pointer to the standard output defined in InitializeMPI()
Definition: global.cpp:75
int NumAve
Definition: global.cpp:43
int Total2Sz
Total in this process.
Definition: struct.hpp:69
int Nsite
Number of sites in the INTRA process region.
Definition: struct.hpp:56
double MaxMPI_d(double dvalue)
MPI wrapper function to obtain maximum Double across processes.
Definition: wrapperMPI.cpp:188
long int MaxMPI_li(long int idim)
MPI wrapper function to obtain maximum unsigned long integer across processes.
Definition: wrapperMPI.cpp:171
int iFlgScaLAPACK
ScaLAPACK mode ( only for FullDiag )
Definition: struct.hpp:237
int CheckMPI(struct BindStruct *X)
Define the number of sites in each PE (DefineList.Nsite). Reduce the number of electrons (DefineList...
Definition: CheckMPI.cpp:27
int NsiteMPI
Total number of sites, differ from DefineList::Nsite.
Definition: struct.hpp:57
int Ne
Number of electrons in this process.
Definition: struct.hpp:71
int GetSplitBitForGeneralSpin(const int Nsite, long int *ihfbit, const long int *SiteToBit)
function of getting right, left and half bits corresponding to a original space.
Definition: bitcalc.cpp:124
int GetLocal2Sz(const int isite, const long int org_bit, const long int *SiteToBit, const long int *Tpow)
get 2sz at a site for general spin
Definition: bitcalc.cpp:448
int iFlgGeneralSpin
Flag for the general (Sz/=1/2) spin.
Definition: struct.hpp:86
long int * SiteToBit
[DefineList::NsiteMPI] Similar to DefineList::Tpow. For general spin.
Definition: struct.hpp:94
long int * Tpow
[2 * DefineList::NsiteMPI] malloc in setmem_def().
Definition: struct.hpp:90
int iCalcModel
Switch for model. 0:Hubbard, 1:Spin, 2:Kondo, 3:HubbardGC, 4:SpinGC, 5:KondoGC, 6:HubbardNConserved.
Definition: struct.hpp:200
int iFlgCalcSpec
Input parameter CalcSpec in teh CalcMod file.
Definition: struct.hpp:218
double max_mem
Estimated memory size.
Definition: struct.hpp:310
int k_exct
Read from Calcmod in readdef.h.
Definition: struct.hpp:47
int Ndown
Number of spin-down electrons in this process.
Definition: struct.hpp:59
int NLocSpn
Number of local spins.
Definition: struct.hpp:84
struct CheckList Check
Size of the Hilbert space.
Definition: struct.hpp:396
long int sdim
Dimension for Ogata-Lin ???
Definition: struct.hpp:309
void CheckMPI_Summary(struct BindStruct *X)
Print infomation of MPI parallelization Modify Definelist::Tpow in the inter process region...
Definition: CheckMPI.cpp:323
long int idim_max
The dimension of the Hilbert space of this process.
Definition: struct.hpp:305
int Total2SzMPI
Total across processes.
Definition: struct.hpp:70
int iCalcType
Switch for calculation type. 0:Lanczos, 1:TPQCalc, 2:FullDiag.
Definition: struct.hpp:194
int childfopenMPI(const char *_cPathChild, const char *_cmode, FILE **_fp)
Only the root process open file in output/ directory.
Definition: FileIO.cpp:27