HΦ  3.1.0
CalcSpectrum.c
Go to the documentation of this file.
1 /* HPhi - Quantum Lattice Model Simulator */
2 /* Copyright (C) 2015 Takahiro Misawa, Kazuyoshi Yoshimi, Mitsuaki Kawamura, Youhei Yamaji, Synge Todo, Naoki Kawashima */
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 #include "mltply.h"
17 #include "CalcSpectrum.h"
18 #include "CalcSpectrumByLanczos.h"
19 #include "CalcSpectrumByBiCG.h"
20 #include "CalcSpectrumByTPQ.h"
21 #include "CalcSpectrumByFullDiag.h"
22 #include "CalcTime.h"
23 #include "SingleEx.h"
24 #include "PairEx.h"
25 #include "wrapperMPI.h"
26 #include "FileIO.h"
27 #include "mfmemory.h"
28 #include "readdef.h"
29 #include "sz.h"
30 #include "check.h"
31 #include "diagonalcalc.h"
42 int OutputSpectrum(
50  struct EDMainCalStruct *X,
51  int Nomega,
52  double complex *dcSpectrum,
53  double complex *dcomega)
54 {
55  FILE *fp;
56  char sdt[D_FileNameMax];
57  int i;
58 
59  //output spectrum
61  if(childfopenMPI(sdt, "w", &fp)!=0){
62  return FALSE;
63  }
64 
65  for (i = 0; i < Nomega; i++) {
66  fprintf(fp, "%.10lf %.10lf %.10lf %.10lf \n",
67  creal(dcomega[i]-X->Bind.Def.dcOmegaOrg), cimag(dcomega[i]-X->Bind.Def.dcOmegaOrg),
68  creal(dcSpectrum[i]), cimag(dcSpectrum[i]));
69  }/*for (i = 0; i < Nomega; i++)*/
70 
71  fclose(fp);
72  return TRUE;
73 }/*int OutputSpectrum*/
74 
91  struct EDMainCalStruct *X
92  ) {
93  char sdt[D_FileNameMax];
94  char *defname;
95  unsigned long int i;
96  unsigned long int i_max = 0;
97  int i_stp;
98  int iFlagListModified = FALSE;
99  FILE *fp;
100  double dnorm;
101 
102  //ToDo: Nomega should be given as a parameter
103  int Nomega;
104  double complex OmegaMax, OmegaMin;
105  double complex *dcSpectrum;
106  double complex *dcomega;
107  size_t byte_size;
108 
109  //set omega
110  if (SetOmega(&(X->Bind.Def)) != TRUE) {
111  fprintf(stderr, "Error: Fail to set Omega.\n");
112  exitMPI(-1);
113  } else {
114  if (X->Bind.Def.iFlgSpecOmegaOrg == FALSE) {
115  X->Bind.Def.dcOmegaOrg = I*(X->Bind.Def.dcOmegaMax - X->Bind.Def.dcOmegaMin) / (double) X->Bind.Def.iNOmega;
116  }
117  }
118  /*
119  Set & malloc omega grid
120  */
121  Nomega = X->Bind.Def.iNOmega;
122  c_malloc1(dcSpectrum, Nomega);
123  c_malloc1(dcomega, Nomega);
124  OmegaMax = X->Bind.Def.dcOmegaMax + X->Bind.Def.dcOmegaOrg;
125  OmegaMin = X->Bind.Def.dcOmegaMin + X->Bind.Def.dcOmegaOrg;
126  for (i = 0; i < Nomega; i++) {
127  dcomega[i] = (OmegaMax - OmegaMin) / Nomega * i + OmegaMin;
128  }
129  fprintf(stdoutMPI, "\nFrequency range:\n");
130  fprintf(stdoutMPI, " Omega Max. : %15.5e %15.5e\n", creal(OmegaMax), cimag(OmegaMax));
131  fprintf(stdoutMPI, " Omega Min. : %15.5e %15.5e\n", creal(OmegaMin), cimag(OmegaMin));
132  fprintf(stdoutMPI, " Num. of Omega : %d\n", Nomega);
133 
135  fprintf(stderr, "Error: Any excitation operators are not defined.\n");
136  exitMPI(-1);
137  }
138  //Make New Lists
139  if (MakeExcitedList(&(X->Bind), &iFlagListModified) == FALSE) {
140  return FALSE;
141  }
142  X->Bind.Def.iFlagListModified=iFlagListModified;
143 
144  //Set Memory
145  c_malloc1(v1Org, X->Bind.Check.idim_maxOrg+1);
146  for(i=0; i<X->Bind.Check.idim_maxOrg+1; i++){
147  v1Org[i]=0;
148  }
149 
150  //Make excited state
151  StartTimer(6100);
152  if (X->Bind.Def.iFlgCalcSpec == RECALC_NOT ||
153  X->Bind.Def.iFlgCalcSpec == RECALC_OUTPUT_TMComponents_VEC ||
154  (X->Bind.Def.iFlgCalcSpec == RECALC_INOUT_TMComponents_VEC && X->Bind.Def.iCalcType == CG)) {
155  //input eigen vector
156  StartTimer(6101);
157  fprintf(stdoutMPI, " Start: An Eigenvector is inputted in CalcSpectrum.\n");
159  GetFileNameByKW(KWSpectrumVec, &defname);
160  strcat(defname, "_rank_%d.dat");
161 // sprintf(sdt, cFileNameInputEigen, X->Bind.Def.CDataFileHead, X->Bind.Def.k_exct - 1, myrank);
162  sprintf(sdt, defname, myrank);
163  childfopenALL(sdt, "rb", &fp);
164 
165  if (fp == NULL) {
166  fprintf(stderr, "Error: A file of Inputvector does not exist.\n");
167  return -1;
168  }
169 
170  byte_size = fread(&i_stp, sizeof(i_stp), 1, fp);
171  X->Bind.Large.itr = i_stp; //For TPQ
172  byte_size = fread(&i_max, sizeof(i_max), 1, fp);
173  if (i_max != X->Bind.Check.idim_maxOrg) {
174  fprintf(stderr, "Error: myrank=%d, i_max=%ld\n", myrank, i_max);
175  fprintf(stderr, "Error: A file of Inputvector is incorrect.\n");
176  return -1;
177  }
178  byte_size = fread(v1Org, sizeof(complex double), i_max + 1, fp);
179  fclose(fp);
180  StopTimer(6101);
181  if (byte_size == 0) printf("byte_size: %d \n", (int)byte_size);
182 
183  for (i = 0; i <= X->Bind.Check.idim_max; i++) {
184  v0[i] = 0;
185  }
186  fprintf(stdoutMPI, " End: An Inputcector is inputted in CalcSpectrum.\n\n");
189  fprintf(stdoutMPI, " Start: Calculating an excited Eigenvector.\n");
190 
191  //mltply Operator
192  StartTimer(6102);
193  GetExcitedState(&(X->Bind), v0, v1Org);
194  StopTimer(6102);
195 
196  //calculate norm
197  dnorm = NormMPI_dc(X->Bind.Check.idim_max, v0);
198  if (fabs(dnorm) < pow(10.0, -15)) {
199  fprintf(stderr, "Warning: Norm of an excitation vector becomes 0.\n");
200  fprintf(stdoutMPI, " End: Calculating an excited Eigenvector.\n\n");
202  fprintf(stdoutMPI, " End: Calculating a spectrum.\n\n");
204  for (i = 0; i < Nomega; i++) {
205  dcSpectrum[i]=0;
206  }
207  OutputSpectrum(X, Nomega, dcSpectrum, dcomega);
208  return TRUE;
209  }
210  //normalize vector
211 #pragma omp parallel for default(none) private(i) shared(v1, v0) firstprivate(i_max, dnorm, X)
212  for (i = 1; i <= X->Bind.Check.idim_max; i++) {
213  v1[i] = v0[i] / dnorm;
214  }
215 
216  fprintf(stdoutMPI, " End: Calculating an excited Eigenvector.\n\n");
218  }
219  StopTimer(6100);
220  //Reset list_1, list_2_1, list_2_2
221  if (iFlagListModified == TRUE) {
222  free(v1Org);
223  free(list_1_org);
224  free(list_2_1_org);
225  free(list_2_2_org);
226  }
227  //calculate Diagonal term
228  diagonalcalc(&(X->Bind));
229 
230 
231  int iret=TRUE;
232  fprintf(stdoutMPI, " Start: Calculating a spectrum.\n\n");
234  StartTimer(6200);
235  switch (X->Bind.Def.iCalcType) {
236  case Lanczos:
237 
238  iret = CalcSpectrumByLanczos(X, v1, dnorm, Nomega, dcSpectrum, dcomega);
239 
240  if (iret != TRUE) {
241  //Error Message will be added.
242  return FALSE;
243  }
244 
245  break;//Lanczos Spectrum
246 
247  case CG:
248 
249  iret = CalcSpectrumByBiCG(X, v0, v1, vg, Nomega, dcSpectrum, dcomega);
250 
251  if (iret != TRUE) {
252  //Error Message will be added.
253  return FALSE;
254  }
255 
256  break;//Lanczos Spectrum
257 
258  case TPQCalc:
259  fprintf(stderr, " Error: TPQ is not supported for calculating spectrum mode.\n");
260  return FALSE;//TPQ is not supprted.
261 #ifdef _CALCSPEC_TPQ
262  iret = CalcSpectrumByTPQ(X, v1, dnorm, Nomega, dcSpectrum, dcomega);
263  if (iret != TRUE) {
264  //Error Message will be added.
265  return FALSE;
266  }
267 #endif
268 
269  case FullDiag:
270  iret = CalcSpectrumByFullDiag(X, Nomega, dcSpectrum, dcomega);
271  break;
272 
273  // case CalcSpecByShiftedKlyrov will be added
274  default:
275  break;
276  }
277  StopTimer(6200);
278 
279  if (iret != TRUE) {
280  fprintf(stderr, " Error: The selected calculation type is not supported for calculating spectrum mode.\n");
281  return FALSE;
282  }
283  else {
284  fprintf(stdoutMPI, " End: Calculating a spectrum.\n\n");
286  iret = OutputSpectrum(X, Nomega, dcSpectrum, dcomega);
287  return TRUE;
288  }
289 
290  c_free1(dcSpectrum, Nomega);
291  c_free1(dcomega, Nomega);
292 
293 }/*int CalcSpectrum*/
294 
302 int GetExcitedState
303 (
304  struct BindStruct *X,
305  double complex *tmp_v0,
306  double complex *tmp_v1
307 )
308 {
309  if(X->Def.NSingleExcitationOperator > 0 && X->Def.NPairExcitationOperator > 0){
310  fprintf(stderr, "Error: Both single and pair excitation operators exist.\n");
311  return FALSE;
312  }
313 
314 
315  if(X->Def.NSingleExcitationOperator > 0){
316  if(GetSingleExcitedState(X,tmp_v0, tmp_v1)!=TRUE){
317  return FALSE;
318  }
319  }
320  else if(X->Def.NPairExcitationOperator >0){
321  if(GetPairExcitedState(X,tmp_v0, tmp_v1)!=TRUE){
322  return FALSE;
323  }
324  }
325 
326  return TRUE;
327 }
328 
336 int SetOmega
337 (
338  struct DefineList *X
339 ){
340  FILE *fp;
341  char sdt[D_FileNameMax],ctmp[256];
342  int istp=4;
343  double E1, E2, E3, E4, Emax;
344  long unsigned int iline_countMax=2;
345  long unsigned int iline_count=2;
346 
347 
348  if(X->iFlgSpecOmegaMax == TRUE && X->iFlgSpecOmegaMin == TRUE){
349  return TRUE;
350  }
351  else{
352  if (X->iCalcType == Lanczos || X->iCalcType == FullDiag) {
353  sprintf(sdt, cFileNameLanczosStep, X->CDataFileHead);
354  childfopenMPI(sdt, "r", &fp);
355  if (fp == NULL) {
356  fprintf(stdoutMPI, "Error: xx_Lanczos_Step.dat does not exist.\n");
357  return FALSE;
358  }
359  fgetsMPI(ctmp, 256, fp); //1st line is skipped
360  fgetsMPI(ctmp, 256, fp); //2nd line is skipped
361  while (fgetsMPI(ctmp, 256, fp) != NULL) {
362  iline_count++;
363  }
364  iline_countMax = iline_count;
365  iline_count = 2;
366  rewind(fp);
367  fgetsMPI(ctmp, 256, fp); //1st line is skipped
368  fgetsMPI(ctmp, 256, fp); //2nd line is skipped
369 
370  while (fgetsMPI(ctmp, 256, fp) != NULL) {
371  sscanf(ctmp, "stp=%d %lf %lf %lf %lf %lf\n",
372  &istp,
373  &E1,
374  &E2,
375  &E3,
376  &E4,
377  &Emax);
378  iline_count++;
379  if (iline_count == iline_countMax) break;
380  }
381  fclose(fp);
382  if (istp < 4) {
383  fprintf(stdoutMPI, "Error: Lanczos step must be greater than 4 for using spectrum calculation.\n");
384  return FALSE;
385  }
386  }/*if (X->iCalcType == Lanczos || X->iCalcType == FullDiag)*/
387  else
388  {
389  sprintf(sdt, cFileNameEnergy_Lanczos, X->CDataFileHead);
390  childfopenMPI(sdt, "r", &fp);
391  if (fp == NULL) {
392  fprintf(stdoutMPI, "Error: xx_energy.dat does not exist.\n");
393  return FALSE;
394  }/*if (fp == NULL)*/
395  fgetsMPI(ctmp, 256, fp); //1st line is skipped
396  fgetsMPI(ctmp, 256, fp); //1st line is skipped
397  sscanf(ctmp, " Energy %lf \n", &E1);
398  Emax = LargeValue;
399  }
400  //Read Lanczos_Step
401  if(X->iFlgSpecOmegaMax == FALSE){
402  X->dcOmegaMax= Emax*(double)X->Nsite;
403  }
404  if(X->iFlgSpecOmegaMin == FALSE){
405  X->dcOmegaMin= E1;
406  }
407  }/*Omegamax and omegamin is not specified in modpara*/
408 
409  return TRUE;
410 }
411 
422  struct BindStruct *X,
423  int *iFlgListModifed
424 ) {
425  long int j;
426  *iFlgListModifed = FALSE;
427  //To Get Original space
428  if (check(X) == MPIFALSE) {
429  FinalizeMPI();
430  return -1;
431  }
432 
433  X->Check.idim_maxOrg = X->Check.idim_max;
434  X->Check.idim_maxMPIOrg = X->Check.idim_maxMPI;
435 
436  if (X->Def.NSingleExcitationOperator > 0) {
437  switch (X->Def.iCalcModel) {
438  case HubbardGC:
439  break;
440  case HubbardNConserved:
441  case KondoGC:
442  case Hubbard:
443  case Kondo:
444  *iFlgListModifed = TRUE;
445  break;
446  case Spin:
447  case SpinGC:
448  return FALSE;
449  }
450  } else if (X->Def.NPairExcitationOperator > 0) {
451  switch (X->Def.iCalcModel) {
452  case HubbardGC:
453  case SpinGC:
454  case HubbardNConserved:
455  break;
456  case KondoGC:
457  case Hubbard:
458  case Kondo:
459  case Spin:
460  if (X->Def.PairExcitationOperator[0][1] != X->Def.PairExcitationOperator[0][3]) {
461  *iFlgListModifed = TRUE;
462  }
463  break;
464  }
465  } else {
466  return FALSE;
467  }
468 
469  if (*iFlgListModifed == TRUE) {
470  if(GetlistSize(X)==TRUE) {
471  lui_malloc1(list_1_org, X->Check.idim_max + 1);
472 #ifdef MPI
473  lui_malloc1(list_1buf_org, X->Check.idim_maxMPI + 1);
474 #endif // MPI
475  lui_malloc1(list_2_1_org, X->Large.SizeOflist_2_1);
476  lui_malloc1(list_2_2_org, X->Large.SizeOflist_2_2);
477  if(list_1_org==NULL
478  || list_2_1_org==NULL
479  || list_2_2_org==NULL
480  )
481  {
482  return -1;
483  }
484  for(j =0; j<X->Large.SizeOflist_2_1; j++){
485  list_2_1_org[j]=0;
486  }
487  for(j =0; j<X->Large.SizeOflist_2_2; j++){
488  list_2_2_org[j]=0;
489  }
490 
491  }
492 
493  if (sz(X, list_1_org, list_2_1_org, list_2_2_org) != 0) {
494  return FALSE;
495  }
496 
497  if (X->Def.NSingleExcitationOperator > 0) {
498  switch (X->Def.iCalcModel) {
499  case HubbardGC:
500  break;
501  case HubbardNConserved:
502  if (X->Def.SingleExcitationOperator[0][2] == 1) { //cis
503  X->Def.Ne = X->Def.NeMPI + 1;
504  }
505  else{
506  X->Def.Ne = X->Def.NeMPI - 1;
507  }
508  break;
509  case KondoGC:
510  case Hubbard:
511  case Kondo:
512  if (X->Def.SingleExcitationOperator[0][2] == 1) { //cis
513  X->Def.Ne = X->Def.NeMPI + 1;
514  if (X->Def.SingleExcitationOperator[0][1] == 0) {//up
515  X->Def.Nup = X->Def.NupOrg + 1;
516  X->Def.Ndown=X->Def.NdownOrg;
517  } else {//down
518  X->Def.Nup=X->Def.NupOrg;
519  X->Def.Ndown = X->Def.NdownOrg + 1;
520  }
521  } else {//ajt
522  X->Def.Ne = X->Def.NeMPI - 1;
523  if (X->Def.SingleExcitationOperator[0][1] == 0) {//up
524  X->Def.Nup = X->Def.NupOrg - 1;
525  X->Def.Ndown=X->Def.NdownOrg;
526 
527  } else {//down
528  X->Def.Nup=X->Def.NupOrg;
529  X->Def.Ndown = X->Def.NdownOrg - 1;
530  }
531  }
532  break;
533  case Spin:
534  case SpinGC:
535  return FALSE;
536  }
537  } else if (X->Def.NPairExcitationOperator > 0) {
538  X->Def.Ne=X->Def.NeMPI;
539  switch (X->Def.iCalcModel) {
540  case HubbardGC:
541  case SpinGC:
542  case HubbardNConserved:
543  break;
544  case KondoGC:
545  case Hubbard:
546  case Kondo:
547  if (X->Def.PairExcitationOperator[0][1] != X->Def.PairExcitationOperator[0][3]) {
548  if (X->Def.PairExcitationOperator[0][1] == 0) {//up
549  X->Def.Nup = X->Def.NupOrg + 1;
550  X->Def.Ndown = X->Def.NdownOrg - 1;
551  } else {//down
552  X->Def.Nup = X->Def.NupOrg - 1;
553  X->Def.Ndown = X->Def.NdownOrg + 1;
554  }
555  }
556  break;
557  case Spin:
558  if (X->Def.PairExcitationOperator[0][1] != X->Def.PairExcitationOperator[0][3]) {
559  if (X->Def.iFlgGeneralSpin == FALSE) {
560  if (X->Def.PairExcitationOperator[0][1] == 0) {//down
561  X->Def.Nup = X->Def.NupOrg - 1;
562  X->Def.Ndown = X->Def.NdownOrg + 1;
563  } else {//up
564  X->Def.Nup = X->Def.NupOrg + 1;
565  X->Def.Ndown = X->Def.NdownOrg - 1;
566  }
567  }
568  else{//for general spin
569  X->Def.Total2Sz = X->Def.Total2SzMPI+2*(X->Def.PairExcitationOperator[0][1]-X->Def.PairExcitationOperator[0][3]);
570  }
571  }
572  break;
573  }
574  } else {
575  return FALSE;
576  }
577  //Update Infomation
578  X->Def.Nsite=X->Def.NsiteMPI;
579 
580  if (check(X) == MPIFALSE) {
581  FinalizeMPI();
582  return FALSE;
583  }
584  }
585 
586  //set memory
587  if (setmem_large(X) != 0) {
589  exitMPI(-1);
590  }
591 
592  if (sz(X, list_1, list_2_1, list_2_2) != 0) {
593  return FALSE;
594  }
595 
596  if(X->Def.iCalcModel==HubbardNConserved){
597  X->Def.iCalcModel=Hubbard;
598  }
599 
600 #ifdef _DEBUG
601  if (*iFlgListModifed == TRUE) {
602  for(j=1; j<=X->Check.idim_maxOrg; j++){
603  fprintf(stdout, "Debug1: myrank=%d, list_1_org[ %ld] = %ld\n", myrank, j, list_1_org[j]+myrank*X->Def.OrgTpow[2*X->Def.NsiteMPI-1]);
604  }
605 
606  for(j=1; j<=X->Check.idim_max; j++){
607  fprintf(stdout, "Debug2: myrank=%d, list_1[ %ld] = %ld\n", myrank, j, list_1[j]+myrank* 64);
608  }
609  }
610 #endif
611 
612  return TRUE;
613 }
614 
unsigned int NSingleExcitationOperator
Number of single excitaion operator for spectrum.
Definition: struct.h:182
int childfopenALL(const char *_cPathChild, const char *_cmode, FILE **_fp)
All processes open file in output/ directory.
Definition: FileIO.c:50
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
struct DefineList Def
Definision of system (Hamiltonian) etc.
Definition: struct.h:409
int MakeExcitedList(struct BindStruct *X, int *iFlgListModifed)
Make the lists for the excited state; list_1, list_2_1 and list_2_2 (for canonical ensemble)...
Definition: CalcSpectrum.c:421
void StartTimer(int n)
function for initializing elapse time [start]
Definition: time.c:71
double complex * v1Org
Definition: global.h:40
long unsigned int * list_2_1_org
Definition: global.h:55
int itr
Iteration number.
Definition: struct.h:314
unsigned long int idim_max
The dimension of the Hilbert space of this process.
Definition: struct.h:302
int CalcSpectrumByFullDiag(struct EDMainCalStruct *X, int Nomega, double complex *dcSpectrum, double complex *dcomega)
Compute the Green function with the Lehmann representation and FD .
int sz(struct BindStruct *X, long unsigned int *list_1_, long unsigned int *list_2_1_, long unsigned int *list_2_2_)
generating Hilbert space
Definition: sz.c:63
void StopTimer(int n)
function for calculating elapse time [elapse time=StartTimer-StopTimer]
Definition: time.c:83
const char * c_InputEigenVectorEnd
Definition: LogMessage.c:61
int iFlagListModified
When the Hilbert space of excited state differs from the original one.
Definition: struct.h:216
double complex * v1
Definition: global.h:35
double complex dcOmegaMax
Upper limit of the frequency for the spectrum.
Definition: struct.h:208
const char * c_CalcSpectrumEnd
Definition: LogMessage.c:65
#define D_FileNameMax
Definition: global.h:23
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.h:411
int childfopenMPI(const char *_cPathChild, const char *_cmode, FILE **_fp)
Only the root process open file in output/ directory.
Definition: FileIO.c:27
const char * c_CalcExcitedStateEnd
Definition: LogMessage.c:63
const char * cFileNameLanczosStep
Definition: global.c:39
int check(struct BindStruct *X)
A program to check size of dimension for Hilbert-space.
Definition: check.c:51
int setmem_large(struct BindStruct *X)
Set size of memories for Hamiltonian (Ham, L_vec), 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.c:224
#define TRUE
Definition: global.h:26
long unsigned int * list_1buf_org
Definition: global.h:54
int GetPairExcitedState(struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
Calculating the pair excited state by the pair operator; , where indicates a creation (anti-creat...
Definition: PairEx.c:47
Bind.
Definition: struct.h:418
const char * c_CalcSpectrumStart
Definition: LogMessage.c:64
int GetExcitedState(struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
Parent function to calculate the excited state.
Definition: CalcSpectrum.c:303
double NormMPI_dc(unsigned long int idim, double complex *_v1)
Compute norm of process-distributed vector .
Definition: wrapperMPI.c:289
int iErrCodeMem
Error Message in HPhiMain.c.
Definition: ErrorMessage.c:20
long unsigned int * list_2_2_org
Definition: global.h:56
const char * cFileNameCalcDynamicalGreen
Definition: global.c:51
int CalcSpectrumByBiCG(struct EDMainCalStruct *X, double complex *vrhs, double complex *v2, double complex *v4, int Nomega, double complex *dcSpectrum, double complex *dcomega)
A main function to calculate spectrum by BiCG method In this function, the library is used...
const char * c_InputEigenVectorStart
Definition: LogMessage.c:60
int diagonalcalc(struct BindStruct *X)
Calculate diagonal components and obtain the list, list_diagonal.
Definition: diagonalcalc.c:85
Bind.
Definition: struct.h:408
int CalcSpectrumByLanczos(struct EDMainCalStruct *X, double complex *tmp_v1, double dnorm, int Nomega, double complex *dcSpectrum, double complex *dcomega)
A main function to calculate spectrum by continued fraction expansions.
unsigned long int idim_maxOrg
The local Hilbert-space dimention of original state for the spectrum.
Definition: struct.h:304
int CalcSpectrumByTPQ(struct EDMainCalStruct *X, double complex *tmp_v1, double dnorm, int Nomega, double complex *dcSpectrum, double complex *dcomega)
A main function to calculate spectrum by TPQ (Note: This method is trial)
#define MPIFALSE
Definition: global.h:24
const char * c_CalcExcitedStateStart
Definition: LogMessage.c:62
long unsigned int * list_1_org
Definition: global.h:53
long unsigned int * list_2_1
Definition: global.h:49
int GetlistSize(struct BindStruct *X)
Set size of lists for the canonical ensemble.
Definition: xsetmem.c:379
void FinalizeMPI()
MPI Finitialization wrapper.
Definition: wrapperMPI.c:74
int CalcSpectrum(struct EDMainCalStruct *X)
A main function to calculate spectrum.
Definition: CalcSpectrum.c:90
char * cErrLargeMem
Definition: ErrorMessage.c:35
int SetOmega(struct DefineList *X)
Set target frequencies.
Definition: CalcSpectrum.c:337
int iNOmega
Number of frequencies for spectrum.
Definition: struct.h:211
#define FALSE
Definition: global.h:25
long unsigned int * list_1
Definition: global.h:47
int iFlgSpecOmegaOrg
Whether DefineList::dcOmegaOrg is input or not.
Definition: struct.h:214
unsigned int NPairExcitationOperator
Number of pair excitaion operator for spectrum.
Definition: struct.h:188
const char * cFileNameTimeKeep
Definition: global.c:23
double complex * vg
Definition: global.h:41
long unsigned int * list_2_2
Definition: global.h:50
double complex dcOmegaOrg
Origin limit of the frequency for the spectrum.
Definition: struct.h:210
struct EDMainCalStruct X
Definition: struct.h:431
double complex dcOmegaMin
Lower limit of the frequency for the spectrum.
Definition: struct.h:209
int GetSingleExcitedState(struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
Calculation of single excited state Target System: Hubbard, Kondo.
Definition: SingleEx.c:30
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
int iFlgCalcSpec
Input parameter CalcSpec in teh CalcMod file.
Definition: struct.h:215
int OutputSpectrum(struct EDMainCalStruct *X, int Nomega, double complex *dcSpectrum, double complex *dcomega)
Output spectrum.
Definition: CalcSpectrum.c:49
int GetFileNameByKW(int iKWidx, char **FileName)
function of getting file name labeled by the keyword
Definition: readdef.c:2711
char * fgetsMPI(char *InputString, int maxcount, FILE *fp)
MPI file I/O (get a line, fgets) wrapper. Only the root node (myrank = 0) reads and broadcast string...
Definition: wrapperMPI.c:122
double complex * v0
Definition: global.h:34
const char * cFileNameEnergy_Lanczos
Definition: global.c:40
struct CheckList Check
Size of the Hilbert space.
Definition: struct.h:410
char * CDataFileHead
Read from Calcmod in readdef.h. Header of output file such as Green&#39;s function.
Definition: struct.h:42
double LargeValue
Definition: global.h:64
struct BindStruct Bind
Binded struct.
Definition: struct.h:419
Definision of system (Hamiltonian) etc.
Definition: struct.h:41
int TimeKeeper(struct BindStruct *X, const char *cFileName, const char *cTimeKeeper_Message, const char *cWriteType)
Functions for writing a time log.
Definition: log.c:42
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:164