HΦ  3.2.0
check.c File Reference

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

#include "bitcalc.h"
#include "sz.h"
#include "FileIO.h"
#include "common/setmemory.h"
#include "check.h"
#include "wrapperMPI.h"
#include "CheckMPI.h"
+ Include dependency graph for check.c:

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.c.

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.c.

References Binomial(), cErrNoModel, cFileNameCheckMemory, cFileNameCheckSdim, BindStruct::Check, CheckMPI(), CheckMPI_Summary(), childfopenMPI(), BindStruct::Def, FALSE, 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(), mfint, MPIFALSE, DefineList::Ndown, DefineList::Ne, DefineList::NLocSpn, DefineList::Nsite, DefineList::NsiteMPI, DefineList::Nup, CheckList::sdim, DefineList::SiteToBit, stdoutMPI, DefineList::Total2Sz, DefineList::Total2SzMPI, DefineList::Tpow, and TRUE.

Referenced by main(), and MakeExcitedList().

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