HΦ  3.1.0
readdef.c
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  *[ver.2008.11.4]
18  * Read Definition files
19  *-------------------------------------------------------------
20  * Copyright (C) 2007-2009 Daisuke Tahara. All rights reserved.
21  *-------------------------------------------------------------*/
22 
23 
24 /*=================================================================================================*/
25 
35 #include "Common.h"
36 #include "readdef.h"
37 #include <ctype.h>
38 #include "LogMessage.h"
39 #include "wrapperMPI.h"
40 #include "mfmemory.h"
41 
46 static char cKWListOfFileNameList[][D_CharTmpReadDef]={
47  "CalcMod",
48  "ModPara",
49  "LocSpin",
50  "Trans",
51  "CoulombIntra",
52  "CoulombInter",
53  "Hund",
54  "PairHop",
55  "Exchange",
56  "InterAll",
57  "OneBodyG",
58  "TwoBodyG",
59  "PairLift",
60  "Ising",
61  "Boost",
62  "SingleExcitation",
63  "PairExcitation",
64  "SpectrumVec",
65  "Laser",
66  "TEOneBody",
67  "TETwoBody"
68 };
69 
71 
75 static char (*cFileNameListFile)[D_CharTmpReadDef];
76 
78  const int NTETrans,
79  const int idx);
80 
81 
83  int iCalcModel,
84  int Nsite,
85  int iFlgGeneralSpin,
86  int *iLocSpin,
87  int isite1, int isigma1,
88  int isite2, int isigma2,
89  int isite3, int isigma3,
90  int isite4, int isigma4
91 );
92 
94  int *icnt_interall,
95  int **iInterAllInfo,
96  double complex *cInterAllValue,
97  int isite1, int isigma1,
98  int isite2, int isigma2,
99  int isite3, int isigma3,
100  int isite4, int isigma4,
101  double re_value, double im_value
102 );
103 
104 
113  const char *defname
114  ){
115  fprintf(stdoutMPI, cErrReadDefFile, defname);
116  return (-1);
117 }
118 
131  const int icheckValue,
132  const int ilowestValue,
133  const int iHighestValue
134  ){
135 
136  if(icheckValue < ilowestValue || icheckValue > iHighestValue){
137  return(-1);
138  }
139  return 0;
140 }
141 
155  const char* cKW,
156  char cKWList[][D_CharTmpReadDef],
157  int iSizeOfKWidx,
158  int* iKWidx
159  ){
160  *iKWidx=-1;
161  int itmpKWidx;
162  int iret=-1;
163  for(itmpKWidx=0; itmpKWidx<iSizeOfKWidx; itmpKWidx++){
164  if(strcmp(cKW,"")==0){
165  break;
166  }
167  else if(CheckWords(cKW, cKWList[itmpKWidx])==0){
168  iret=0;
169  *iKWidx=itmpKWidx;
170  }
171  }
172  return iret;
173 }
174 
188  char *ctmpLine,
189  char *ctmp,
190  int *itmp
191  )
192 {
193  char *ctmpRead;
194  char *cerror;
195  char csplit[] = " ,.\t\n";
196  if(*ctmpLine=='\n') return 1;
197  ctmpRead = strtok(ctmpLine, csplit);
198  if(strncmp(ctmpRead, "=", 1)==0 || strncmp(ctmpRead, "#", 1)==0 || ctmpRead==NULL){
199  return 1;
200  }
201  strcpy(ctmp, ctmpRead);
202 
203  ctmpRead = strtok( NULL, csplit );
204  *itmp = strtol(ctmpRead, &cerror, 0);
205  //if ctmpRead is not integer type
206  if(*cerror != '\0'){
207  fprintf(stdoutMPI, cErrDefFileFormat, cerror);
208  return(-1);
209  }
210 
211  ctmpRead = strtok( NULL, csplit );
212  if(ctmpRead != NULL){
213  fprintf(stdoutMPI, cErrDefFileFormat, ctmpRead);
214  return(-1);
215  }
216 
217  return 0;
218 }
219 
231  const char *defname,
232  struct DefineList *X
233  )
234 {
235  FILE *fp;
236  int itmp, iret;
237  char ctmpLine[D_CharTmpReadDef+D_CharKWDMAX];
238  char ctmp[D_CharKWDMAX];
239  X->iCalcType=0;
240  X->iFlgFiniteTemperature=0;
241  X->iCalcModel=0;
242  X->iOutputMode=0;
243  X->iCalcEigenVec=0;
244  X->iInitialVecType=0;
245  X->iOutputEigenVec=0;
246  X->iInputEigenVec=0;
247  X->iOutputHam=0;
248  X->iInputHam=0;
249  X->iFlgCalcSpec=0;
250  X->iReStart=0;
251  X->iFlgMPI=0;
252  X->iFlgScaLAPACK=0;
253 #ifdef _MAGMA
254  X->iNGPU=2;
255 #else
256  X->iNGPU=0;
257 #endif
258  /*=======================================================================*/
259  fp = fopenMPI(defname, "r");
260  if(fp==NULL) return ReadDefFileError(defname);
261  /* read Parameters from calcmod.def*/
262  while( fgetsMPI(ctmpLine, D_CharTmpReadDef+D_CharKWDMAX, fp)!=NULL ){
263  if( (iret=GetKWWithIdx(ctmpLine, ctmp, &itmp)) !=0){
264  if(iret==1) continue;
265  return(-1);
266  }
267  if(CheckWords(ctmp, "CalcType")==0){
268  X->iCalcType=itmp;
269  }
270  else if(CheckWords(ctmp, "FlgFiniteTemperature")==0){
271  X->iFlgFiniteTemperature = itmp;
272  }
273  else if(CheckWords(ctmp, "CalcModel")==0){
274  X->iCalcModel=itmp;
275  }
276  else if(CheckWords(ctmp, "OutputMode")==0){
277  X->iOutputMode=itmp;
278  }
279  else if(CheckWords(ctmp, "CalcEigenVec")==0){
280  X->iCalcEigenVec=itmp;
281  }
282  else if(CheckWords(ctmp, "InitialVecType")==0){
283  X->iInitialVecType=itmp;
284  }
285  else if(CheckWords(ctmp, "OutputEigenVec")==0 || CheckWords(ctmp, "OEV")==0){
286  X->iOutputEigenVec=itmp;
287  }
288  else if(CheckWords(ctmp, "InputEigenVec")==0 || CheckWords(ctmp, "IEV")==0){
289  X->iInputEigenVec=itmp;
290  }
291  else if(CheckWords(ctmp, "OutputHam")==0){
292  X->iOutputHam=itmp;
293  }
294  else if(CheckWords(ctmp, "InputHam")==0){
295  X->iInputHam=itmp;
296  }
297  else if(CheckWords(ctmp, "CalcSpec")==0 || CheckWords(ctmp, "CalcSpectrum")==0){
298  X->iFlgCalcSpec=itmp;
299  }
300  else if(CheckWords(ctmp, "ReStart")==0){
301  X->iReStart=itmp;
302  }
303  else if(CheckWords(ctmp, "NGPU")==0){
304  X->iNGPU=itmp;
305  }
306  else if(CheckWords(ctmp, "ScaLAPACK")==0){
307 #ifdef _SCALAPACK
308  X->iFlgScaLAPACK=itmp;
309 #endif
310  }
311  else{
312  fprintf(stdoutMPI, cErrDefFileParam, defname, ctmp);
313  return(-1);
314  }
315  }
316  fclose(fp);
317 
318  /* Check values*/
319  if(ValidateValue(X->iCalcModel, 0, NUM_CALCMODEL-1)){
320  fprintf(stdoutMPI, cErrCalcType, defname);
321  return (-1);
322  }
323  if(ValidateValue(X->iCalcType, 0, NUM_CALCTYPE-1)){
324  fprintf(stdoutMPI, cErrCalcType, defname);
325  return (-1);
326  }
327  if(ValidateValue(X->iOutputMode, 0, NUM_OUTPUTMODE-1)){
328  fprintf(stdoutMPI, cErrOutputMode, defname);
329  return (-1);
330  }
331 
332  if(ValidateValue(X->iCalcEigenVec, -1, NUM_CALCEIGENVEC-1)){
333  fprintf(stdoutMPI, cErrCalcEigenVec, defname);
334  return (-1);
335  }
336 
337  if(ValidateValue(X->iInitialVecType, 0, NUM_SETINITAILVEC-1)){
338  fprintf(stdoutMPI, cErrSetIniVec, defname);
339  return (-1);
340  }
341 
342  if(ValidateValue(X->iOutputHam, 0, NUM_OUTPUTHAM-1)){
343  fprintf(stdoutMPI, cErrOutputHam, defname);
344  return (-1);
345  }
346  if(ValidateValue(X->iInputHam, 0, NUM_INPUTHAM-1)){
347  fprintf(stdoutMPI, cErrInputHam, defname);
348  return (-1);
349  }
350  if(X->iInputHam == 1 && X->iOutputHam==1){
351  fprintf(stdoutMPI, cErrInputOutputHam, defname);
352  return (-1);
353  }
354  if(ValidateValue(X->iReStart, 0, NUM_RESTART-1)){
355  fprintf(stdoutMPI, cErrRestart, defname);
356  return (-1);
357  }
358  if(X->iNGPU < 0){
359  fprintf(stdoutMPI, cErrCUDA, defname);
360  return (-1);
361  }
362  if(ValidateValue(X->iFlgScaLAPACK, 0, 1)){
363  fprintf(stdoutMPI, cErrCUDA, defname);
364  return (-1);
365  }
366 
367  /* In the case of Full Diagonalization method(iCalcType=2)*/
368  if(X->iCalcType==2 && ValidateValue(X->iFlgFiniteTemperature, 0, 1)){
369  fprintf(stdoutMPI, cErrFiniteTemp, defname);
370  return (-1);
371  }
372 
373  if(X->iCalcType !=2 && X->iOutputHam ==TRUE) {
374  fprintf(stdoutMPI, cErrOutputHamForFullDiag, defname);
375  return (-1);
376  }
377 
378  return 0;
379 }
380 
392  const char* cFileListNameFile,
393  char cFileNameList[][D_CharTmpReadDef]
394  )
395 {
396  FILE *fplist;
397  int itmpKWidx=-1;
398  char ctmpFileName[D_FileNameMaxReadDef];
399  char ctmpKW[D_CharTmpReadDef], ctmp2[256];
400  int i;
401  for(i=0; i< D_iKWNumDef; i++){
402  strcpy(cFileNameList[i],"");
403  }
404 
405  fplist = fopenMPI(cFileListNameFile, "r");
406  if(fplist==NULL) return ReadDefFileError(cFileListNameFile);
407 
408  while(fgetsMPI(ctmp2, 256, fplist) != NULL){
409  memset(ctmpKW, '\0', strlen(ctmpKW));
410  memset(ctmpFileName, '\0', strlen(ctmpFileName));
411  sscanf(ctmp2,"%s %s\n", ctmpKW, ctmpFileName);
412 
413  if(strncmp(ctmpKW, "#", 1)==0 || *ctmp2=='\n' || (strcmp(ctmpKW, "")&&strcmp(ctmpFileName,""))==0){
414  continue;
415  }
416  else if(strcmp(ctmpKW, "")*strcmp(ctmpFileName, "")==0){
417  fprintf(stdoutMPI, cErrKW_InCorPair, cFileListNameFile);
418  fclose(fplist);
419  return(-1);
420  }
422  if( CheckKW(ctmpKW, cKWListOfFileNameList, D_iKWNumDef, &itmpKWidx)!=0 ){
423  fprintf(stdoutMPI, cErrKW, ctmpKW, cFileListNameFile);
424  fprintf(stdoutMPI, "%s", cErrKW_ShowList);
425  for(i=0; i<D_iKWNumDef;i++){
426  fprintf(stdoutMPI, "%s \n", cKWListOfFileNameList[i]);
427  }
428  fclose(fplist);
429  return(-1);
430  }
432  if(strcmp(cFileNameList[itmpKWidx], "") !=0){
433  fprintf(stdoutMPI, cErrKW_Same, cFileListNameFile);
434  fclose(fplist);
435  return(-1);
436  }
437 
439  strcpy(cFileNameList[itmpKWidx], ctmpFileName);
440  }
441  fclose(fplist);
442  return 0;
443 }
444 
457  char *xNameListFile,
458  struct DefineList *X,
459  struct BoostList *xBoost
460  )
461 {
462  FILE *fp;
463  char defname[D_FileNameMaxReadDef];
464  char ctmp[D_CharTmpReadDef], ctmp2[256];
465  int i,itmp;
466  unsigned int iline=0;
467  X->nvec=0;
468  X->iFlgSpecOmegaMax=FALSE;
469  X->iFlgSpecOmegaMin=FALSE;
470  X->iFlgSpecOmegaOrg=FALSE;
471  X->iNOmega=1000;
472  X->NCond=0;
473  X->iFlgSzConserved=FALSE;
474  X->dcOmegaOrg=0;
475  int iReadNCond=FALSE;
476  xBoost->flgBoost=FALSE;
478  NumAve=1;
479  X->Param.ExpecInterval=1;
480  cFileNameListFile = malloc(sizeof(char)*D_CharTmpReadDef*D_iKWNumDef);
481 
482  fprintf(stdoutMPI, cReadFileNamelist, xNameListFile);
483  if(GetFileName(xNameListFile, cFileNameListFile)!=0){
484  return(-1);
485  }
486 
487  /*=======================================================================*/
488  int iKWidx=0;
489  //Check the existence of Essensial Files.
490  X->READ=0;
491  X->WRITE=0;
492 
493  for(iKWidx=0; iKWidx< D_iKWNumDef; iKWidx++){
494  strcpy(defname, cFileNameListFile[iKWidx]);
495  if(strcmp(defname,"")==0){
496  switch (iKWidx){
497  case KWCalcMod:
498  case KWModPara:
499  case KWLocSpin:
500  fprintf(stdoutMPI, cErrMakeDef, cKWListOfFileNameList[iKWidx]);
501  return(-1);
502  default:
503  break;
504  }
505  }
506  }
507 
508 
509  for(iKWidx=0; iKWidx< D_iKWNumDef; iKWidx++) {
510  strcpy(defname, cFileNameListFile[iKWidx]);
511 
512  if (strcmp(defname, "") == 0) continue;
513  if(iKWidx==KWSpectrumVec){
514  continue;
515  }
516  fprintf(stdoutMPI, cReadFile, defname, cKWListOfFileNameList[iKWidx]);
517  fp = fopenMPI(defname, "r");
518  if (fp == NULL) return ReadDefFileError(defname);
519  switch (iKWidx) {
520  case KWCalcMod:
521  /* Read calcmod.def---------------------------------------*/
522  if (ReadcalcmodFile(defname, X) != 0) {
523  fclose(fp);
524  return ReadDefFileError(defname);
525  }
526  break;
527 
528  case KWModPara:
529  /* Read modpara.def---------------------------------------*/
530  //TODO: add error procedure here when parameters are not enough.
532  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp);
533  fgetsMPI(ctmp2, 256, fp);
534  sscanf(ctmp2, "%s %d\n", ctmp, &itmp); //2
535  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp); //3
536  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp); //4
537  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp); //5
539  fgetsMPI(ctmp2, 256, fp);
540  sscanf(ctmp2, "%s %s\n", ctmp, X->CDataFileHead); //6
542  fgetsMPI(ctmp2, 256, fp);
543  sscanf(ctmp2, "%s %s\n", ctmp, X->CParaFileHead); //7
545  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp); //8
546  double dtmp, dtmp2;
547  X->read_hacker = 1;
549  while (fgetsMPI(ctmp2, 256, fp) != NULL) {
550  if (*ctmp2 == '\n') continue;
551  sscanf(ctmp2, "%s %lf %lf\n", ctmp, &dtmp, &dtmp2);
552  if (CheckWords(ctmp, "Nsite") == 0) {
553  X->Nsite = (int) dtmp;
554  }
555  else if (CheckWords(ctmp, "Nup") == 0) {
556  X->Nup = (int) dtmp;
557  }
558  else if (CheckWords(ctmp, "Ndown") == 0) {
559  X->Ndown = (int) dtmp;
560  X->Total2Sz = X->Nup - X->Ndown;
561  }
562  else if (CheckWords(ctmp, "2Sz") == 0) {
563  X->Total2Sz = (int) dtmp;
564  X->iFlgSzConserved = TRUE;
565  }
566  else if (CheckWords(ctmp, "Ncond") == 0) {
567  if((int) dtmp <0) {
568  fprintf(stdoutMPI, cErrNcond, defname);
569  return (-1);
570  }
571  X->NCond = (int) dtmp;
572  iReadNCond = TRUE;
573  }
574  else if (CheckWords(ctmp, "Lanczos_max") == 0) {
575  X->Lanczos_max = (int) dtmp;
576  }
577  else if (CheckWords(ctmp, "initial_iv") == 0) {
578  X->initial_iv = (int) dtmp;
579  }
580  else if (CheckWords(ctmp, "nvec") == 0) {
581  X->nvec = (int) dtmp;
582  }
583  else if (CheckWords(ctmp, "exct") == 0) {
584  X->k_exct = (int) dtmp;
585  }
586  else if (CheckWords(ctmp, "LanczosEps") == 0) {
587  X->LanczosEps = (int) dtmp;
588  }
589  else if (CheckWords(ctmp, "LanczosTarget") == 0) {
590  X->LanczosTarget = (int) dtmp;
591  }
592  else if (CheckWords(ctmp, "LargeValue") == 0) {
593  LargeValue = dtmp;
594  }
595  else if (CheckWords(ctmp, "NumAve") == 0) {
596  NumAve = (int) dtmp;
597  }
598  else if(strcmp(ctmp, "TimeSlice")==0){
599  X->Param.TimeSlice=dtmp;
600  }
601  else if(strcmp(ctmp, "ExpandCoef")==0){
602  X->Param.ExpandCoef=(int)dtmp;
603  }
604  else if(strcmp(ctmp, "OutputInterval")==0){
605  X->Param.OutputInterval=(int)dtmp;
606  }
607  else if (CheckWords(ctmp, "ExpecInterval") == 0) {
608  X->Param.ExpecInterval = (int) dtmp;
609  }
610  else if(strcmp(ctmp, "Tinit")==0){
611  X->Param.Tinit=dtmp;
612  }
613  else if (CheckWords(ctmp, "CalcHS") == 0) {
614  X->read_hacker = (int) dtmp;
615  }
616  else if(CheckWords(ctmp, "OmegaMax")==0){
617  X->dcOmegaMax=dtmp+dtmp2*I;
618  X->iFlgSpecOmegaMax=TRUE;
619  }
620  else if(CheckWords(ctmp, "OmegaMin")==0){
621  X->dcOmegaMin =dtmp+dtmp2*I;
622  X->iFlgSpecOmegaMin=TRUE;
623  }
624  else if(CheckWords(ctmp, "OmegaIm")==0){
625  X->dcOmegaOrg +=dtmp*I;
626  X->iFlgSpecOmegaOrg=TRUE;
627  }
628  else if(CheckWords(ctmp, "OmegaOrg")==0){
629  X->dcOmegaOrg +=dtmp+dtmp2*I;
630  X->iFlgSpecOmegaOrg=TRUE;
631  }
632  else if(CheckWords(ctmp, "NOmega")==0){
633  X->iNOmega=(int)dtmp;
634  }
635  else if(CheckWords(ctmp, "TargetTPQRand")==0) {
636  X->irand=(int)dtmp;
637  }
638  else {
639  return (-1);
640  }
641  }
642  break;
643 
644  case KWLocSpin:
645  // Read locspn.def
646  X->iFlgGeneralSpin = FALSE;
647  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp);
648  fgetsMPI(ctmp2, 256, fp);
649  sscanf(ctmp2, "%s %d\n", ctmp, &(X->NLocSpn));
650  break;
651  case KWTrans:
652  // Read transfer.def
653  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp);
654  fgetsMPI(ctmp2, 256, fp);
655  sscanf(ctmp2, "%s %d\n", ctmp, &(X->NTransfer));
656  break;
657  case KWCoulombIntra:
658  /* Read coulombintra.def----------------------------------*/
659  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp);
660  fgetsMPI(ctmp2, 256, fp);
661  sscanf(ctmp2, "%s %d\n", ctmp, &(X->NCoulombIntra));
662  break;
663  case KWCoulombInter:
664  /* Read coulombinter.def----------------------------------*/
665  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp);
666  fgetsMPI(ctmp2, 256, fp);
667  sscanf(ctmp2, "%s %d\n", ctmp, &(X->NCoulombInter));
668  break;
669  case KWHund:
670  /* Read hund.def------------------------------------------*/
671  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp);
672  fgetsMPI(ctmp2, 256, fp);
673  sscanf(ctmp2, "%s %d\n", ctmp, &(X->NHundCoupling));
674  break;
675  case KWPairHop:
676  /* Read pairhop.def---------------------------------------*/
677  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp);
678  fgetsMPI(ctmp2, 256, fp);
679  sscanf(ctmp2, "%s %d\n", ctmp, &(X->NPairHopping));
680  X->NPairHopping*=2;
681  break;
682  case KWExchange:
683  /* Read exchange.def--------------------------------------*/
684  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp);
685  fgetsMPI(ctmp2, 256, fp);
686  sscanf(ctmp2, "%s %d\n", ctmp, &(X->NExchangeCoupling));
687  break;
688  case KWIsing:
689  /* Read ising.def--------------------------------------*/
690  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp);
691  fgetsMPI(ctmp2, 256, fp);
692  sscanf(ctmp2, "%s %d\n", ctmp, &(X->NIsingCoupling));
693  break;
694  case KWPairLift:
695  /* Read exchange.def--------------------------------------*/
696  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp);
697  fgetsMPI(ctmp2, 256, fp);
698  sscanf(ctmp2, "%s %d\n", ctmp, &(X->NPairLiftCoupling));
699  break;
700  case KWInterAll:
701  /* Read InterAll.def--------------------------------------*/
702  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp);
703  fgetsMPI(ctmp2, 256, fp);
704  sscanf(ctmp2, "%s %d\n", ctmp, &(X->NInterAll));
705  break;
706  case KWOneBodyG:
707  /* Read cisajs.def----------------------------------------*/
708  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp);
709  fgetsMPI(ctmp2, 256, fp);
710  sscanf(ctmp2, "%s %d\n", ctmp, &(X->NCisAjt));
711  break;
712  case KWTwoBodyG:
713  /* Read cisajscktaltdc.def--------------------------------*/
714  fgetsMPI(ctmp, sizeof(ctmp) / sizeof(char), fp);
715  fgetsMPI(ctmp2, 256, fp);
716  sscanf(ctmp2, "%s %d\n", ctmp, &(X->NCisAjtCkuAlvDC));
717  break;
718  case KWLaser:
719  /* Read laser.def--------------------------------*/
720  fgetsMPI(ctmp, sizeof(ctmp)/sizeof(char), fp);
721  fgetsMPI(ctmp2, 256, fp);
722  sscanf(ctmp2,"%s %d\n", ctmp, &(X->NLaser));
723  break;
724 
725  case KWTEOneBody:
726  if(X->iCalcType != TimeEvolution) break;
727  /* Read TEOnebody.def--------------------------------*/
728  fgetsMPI(ctmp, sizeof(ctmp)/sizeof(char), fp);
729  fgetsMPI(ctmp2, 256, fp);
730  sscanf(ctmp2,"%s %d\n", ctmp, &(X->NTETimeSteps));
731  fgetsMPI(ctmp2, 256, fp);
732  fgetsMPI(ctmp2, 256, fp);
733  fgetsMPI(ctmp2, 256, fp);
734  int iTETransMax=0;
735  if(X->NTETimeSteps>0) {
736  while (fgetsMPI(ctmp2, 256, fp) != NULL) {
737  sscanf(ctmp2, "%lf %d \n", &dtmp, &itmp);
738  for (i = 0; i < itmp; ++i) {
739  fgetsMPI(ctmp2, 256, fp);
740  }
741  if(iTETransMax < itmp) iTETransMax=itmp;
742  }
743  }
744  X->NTETransferMax=iTETransMax;
745  break;
746 
747  case KWTETwoBody:
748  if(X->iCalcType != TimeEvolution) break;
749  /* Read TETwobody.def--------------------------------*/
750  fgetsMPI(ctmp, sizeof(ctmp)/sizeof(char), fp);
751  fgetsMPI(ctmp2, 256, fp);
752  sscanf(ctmp2,"%s %d\n", ctmp, &(X->NTETimeSteps));
753  fgetsMPI(ctmp2, 256, fp);
754  fgetsMPI(ctmp2, 256, fp);
755  fgetsMPI(ctmp2, 256, fp);
756  int iTEInterAllMax=0;
757  if(X->NTETimeSteps>0) {
758  while (fgetsMPI(ctmp2, 256, fp) != NULL) {
759  sscanf(ctmp2, "%lf %d \n", &dtmp, &itmp);
760  for (i = 0; i < itmp; ++i) {
761  fgetsMPI(ctmp2, 256, fp);
762  }
763  if(iTEInterAllMax < itmp) iTEInterAllMax=itmp;
764  }
765  }
766  X->NTEInterAllMax=iTEInterAllMax;
767  break;
768 
769 
770  case KWBoost:
771  /* Read boost.def--------------------------------*/
772  xBoost->NumarrayJ = 0;
773  xBoost->W0 = 0;
774  xBoost->R0 = 0;
775  xBoost->num_pivot = 0;
776  xBoost->ishift_nspin = 0;
777  xBoost->flgBoost = TRUE;
778  //first line is skipped
779  fgetsMPI(ctmp2, 256, fp);
780  //read numarrayJ
781  fgetsMPI(ctmp2, 256, fp);
782  sscanf(ctmp2, "%d\n", &(xBoost->NumarrayJ));
783  //skipp arrayJ
784  for (iline = 0; iline < xBoost->NumarrayJ * 3; iline++) {
785  fgetsMPI(ctmp2, 256, fp);
786  }
787  //read W0 R0 num_pivot ishift_nspin
788  fgetsMPI(ctmp2, 256, fp);
789  sscanf(ctmp2, "%ld %ld %ld %ld\n", &(xBoost->W0), &(xBoost->R0), &(xBoost->num_pivot),
790  &(xBoost->ishift_nspin));
791 
792  break;
793 
794  case KWSingleExcitation:
795  /* Read singleexcitation.def----------------------------------------*/
796  fgetsMPI(ctmp, sizeof(ctmp)/sizeof(char), fp);
797  fgetsMPI(ctmp2, 256, fp);
798  sscanf(ctmp2,"%s %d\n", ctmp, &(X->NSingleExcitationOperator));
799  break;
800 
801  case KWPairExcitation:
802  /* Read pairexcitation.def----------------------------------------*/
803  fgetsMPI(ctmp, sizeof(ctmp)/sizeof(char), fp);
804  fgetsMPI(ctmp2, 256, fp);
805  sscanf(ctmp2,"%s %d\n", ctmp, &(X->NPairExcitationOperator));
806  break;
807 
808  default:
809  fprintf(stdoutMPI, "%s", cErrIncorrectDef);
810  fclose(fp);
811  return (-1);
812  }
813  /*=======================================================================*/
814  fclose(fp);
815  }
816 
817  //Sz, Ncond
818  switch(X->iCalcModel){
819  case Spin:
820  case Hubbard:
821  case Kondo:
822  case SpinlessFermion:
823 
824  if(iReadNCond==TRUE){
825  if(X->iCalcModel==Spin){
826  fprintf(stdoutMPI, "For Spin, Ncond should not be defined.\n");
827  return(-1);
828  }
829  else{
830  if(X->iFlgSzConserved==TRUE){
831  if(X->iCalcModel==SpinlessFermion){
832  fprintf(stdoutMPI, " Warning: For Spinless fermion, 2Sz should not be defined.\n");
833  X->Ne=X->NCond;
834  X->Nup=X->NCond;
835  X->Ndown=0;
836  break;
837  }
838  X->Nup=X->NLocSpn+X->NCond+X->Total2Sz;
839  X->Ndown=X->NLocSpn+X->NCond-X->Total2Sz;
840  X->Nup/=2;
841  X->Ndown/=2;
842  }
843  else{
844  if(X->iCalcModel == Hubbard){
845  X->Ne=X->NCond;
846  if(X->Ne <1){
847  fprintf(stdoutMPI, "Ncond is incorrect.\n");
848  return(-1);
849  }
850  X->iCalcModel=HubbardNConserved;
851  }
852  else if(X->iCalcModel ==SpinlessFermion){
853  X->Ne=X->NCond;
854  X->Nup=X->NCond;
855  X->Ndown=0;
856  }
857  else{
858  fprintf(stdoutMPI, " 2Sz is not defined.\n");
859  return(-1);
860  }
861  }
862  }
863  }
864  else if(iReadNCond == FALSE && X->iFlgSzConserved==TRUE){
865  if(X->iCalcModel != Spin){
866  fprintf(stdoutMPI, " NCond is not defined.\n");
867  return(-1);
868  }
869  X->Nup=X->NLocSpn+X->Total2Sz;
870  X->Ndown=X->NLocSpn-X->Total2Sz;
871  X->Nup /= 2;
872  X->Ndown /= 2;
873  }
874  else{
875  if(X->Nup==0 && X->Ndown==0){
876  if(X->iCalcModel == Spin){
877  fprintf(stdoutMPI, " 2Sz is not defined.\n");
878  return(-1);
879  }
880  else{
881  fprintf(stdoutMPI, " NCond is not defined.\n");
882  return(-1);
883  }
884  }
885  }
886 
887  if(X->iCalcModel == Spin){
888  X->Ne=X->Nup;
889  }
890  else{
891  if(X->Ne==0) {
892  X->Ne = X->Nup + X->Ndown;
893  }
894  if(X->NLocSpn>X->Ne){
895  fprintf(stdoutMPI, "%s", cErrNLoc);
896  fprintf(stdoutMPI, "NLocalSpin=%d, Ne=%d\n", X->NLocSpn, X->Ne);
897  return(-1);
898  }
899  }
900  break;
901  case SpinGC:
902  case KondoGC:
903  case HubbardGC:
904  case SpinlessFermionGC:
905  if(iReadNCond == TRUE || X->iFlgSzConserved ==TRUE){
906  fprintf(stdoutMPI, "\n Warning: For GC, both Ncond and 2Sz should not be defined.\n");
907  //return(-1);
908  }
909  break;
910  default:
911  break;
912  }
913 
914  /* Check values (Positive)*/
915  if(X->Nsite<=0) {// Nsite must be positve
916  fprintf(stdoutMPI, cErrNsite, defname);
917  return (-1);
918  }
919  if(X->Lanczos_max<=0) {// Lanczos_max must be positive
920  fprintf(stdoutMPI, cErrLanczos_max, defname);
921  return (-1);
922  }
923  if(X->LanczosEps<=0) {// Lanczos_eps must be positive
924  fprintf(stdoutMPI, cErrLanczos_eps, defname);
925  return (-1);
926  }
927  if(NumAve<=0) { // Average number must be positive
928  fprintf(stdoutMPI, cErrNumAve, defname);
929  return (-1);
930  }
931  if(X->Param.ExpecInterval<=0){// Interval to calculate expected values must be positive
932  fprintf(stdoutMPI, cErrExpecInterval, defname);
933  return (-1);
934  }
935  if(X->nvec==0){
936  X->nvec=X->Lanczos_max;
937  }
938 
939  if(X->nvec < X->k_exct){
940  X->nvec=X->k_exct;
941  }
942  if (X->LanczosTarget < X->k_exct) X->LanczosTarget = X->k_exct;
943 
944  if(ValidateValue(X->k_exct, 1, X->nvec)) {
945  fprintf(stdoutMPI, cErrLanczosExct, defname, X->nvec);
946  return (-1);
947  }
948 
949  if( X->k_exct>X->LanczosTarget ){
950  fprintf(stdoutMPI, cErrLanczosTarget, defname, X->LanczosTarget, X->k_exct);
951  return (-1);
952  }
953 
954 
955  X->fidx = 0;
956  X->NeMPI=X->Ne;
957  X->NupMPI=X->Nup;
958  X->NdownMPI=X->Ndown;
959  X->NupOrg=X->Nup;
960  X->NdownOrg=X->Ndown;
961  return 0;
962 }
963 
977  struct DefineList *X,
978  struct BoostList *xBoost
979  )
980 {
981  FILE *fp;
982  char defname[D_FileNameMaxReadDef];
983  char ctmp[D_CharTmpReadDef], ctmp2[256];
984 
985  unsigned int i,j, idx, itype;
986  int xitmp[8];
987  int iKWidx=0;
988  int iboolLoc=0;
989  int isite1, isite2, isite3, isite4;
990  int isigma1, isigma2, isigma3, isigma4;
991  double dvalue_re, dvalue_im;
992  double dArrayValue_re[3];
993  int icnt_diagonal=0;
994  int ieps_CheckImag0=-12;
995  eps_CheckImag0=pow(10.0, ieps_CheckImag0);
996  unsigned int iline=0;
997  int ilineIn=0;
998  int ilineIn2=0;
999  int itmp=0;
1000  int icnt_trans=0;
1001  int iflg_trans=0;
1002  int icnt_interall=0;
1003  int iflg_interall=0;
1004 
1005  unsigned int iloop=0;
1006 
1007  for(iKWidx=KWLocSpin; iKWidx< D_iKWNumDef; iKWidx++){
1008  strcpy(defname, cFileNameListFile[iKWidx]);
1009  if(strcmp(defname,"")==0 || iKWidx==KWSpectrumVec) continue;
1010  fprintf(stdoutMPI, cReadFileNamelist, defname);
1011  fp = fopenMPI(defname, "r");
1012  if(fp==NULL) return ReadDefFileError(defname);
1013  if(iKWidx != KWBoost){
1014  for(i=0;i<IgnoreLinesInDef;i++) fgetsMPI(ctmp, sizeof(ctmp)/sizeof(char), fp);
1015  }
1016 
1017  idx=0;
1018  /*=======================================================================*/
1019  switch(iKWidx){
1020  case KWLocSpin:
1021  /* Read locspn.def----------------------------------------*/
1022  while( fgetsMPI(ctmp2, 256, fp) != NULL){
1023  if(idx==X->Nsite){
1024  fclose(fp);
1025  return ReadDefFileError(defname);
1026  }
1027 
1028  sscanf(ctmp2, "%d %d\n", &(xitmp[0]), &(xitmp[1]) );
1029  X->LocSpn[xitmp[0]] = xitmp[1];
1030  X->SiteToBit[xitmp[0]]=(X->LocSpn[xitmp[0]]+1);//2S+1
1031  if(CheckSite(xitmp[0], X->Nsite) !=0){
1032  fclose(fp);
1033  return ReadDefFileError(defname);
1034  }
1035  idx++;
1036  }
1037  if(CheckLocSpin(X)==FALSE){
1038  fclose(fp);
1039  return ReadDefFileError(defname);
1040  }
1041 
1042  break;
1043 
1044  case KWTrans:
1045  /* transfer.def--------------------------------------*/
1046  if(X->NTransfer>0){
1047  icnt_trans=0;
1048  while( fgetsMPI(ctmp2, 256, fp) != NULL )
1049  {
1050  if(idx==X->NTransfer){
1051  fclose(fp);
1052  return ReadDefFileError(defname);
1053  }
1054 
1055  sscanf(ctmp2, "%d %d %d %d %lf %lf\n",
1056  &isite1,
1057  &isigma1,
1058  &isite2,
1059  &isigma2,
1060  &dvalue_re,
1061  &dvalue_im
1062  );
1063 
1064  if(CheckPairSite(isite1, isite2,X->Nsite) !=0){
1065  fclose(fp);
1066  return ReadDefFileError(defname);
1067  }
1068 
1069  if(isite1==isite2 && isigma1==isigma2){
1070  if(fabs(dvalue_im)> eps_CheckImag0){
1071  //NonHermite
1072  fprintf(stdoutMPI, cErrNonHermiteTrans, isite1, isigma1, isite2, isigma2, dvalue_re, dvalue_im);
1073  fclose(fp);
1074  return ReadDefFileError(defname);
1075  }
1076  }
1077 
1078  if(X->iCalcModel==Spin){
1079  if(isite1 != isite2){
1080  iboolLoc=1;
1081  fprintf(stdoutMPI, cWarningIncorrectFormatForSpin2, isite1, isite2);
1082  }
1083  }
1084  else if(X->iCalcModel==Kondo){
1085  if(X->LocSpn[isite1]!=ITINERANT || X->LocSpn[isite2] !=ITINERANT){
1086  if(isite1 != isite2){
1087  iboolLoc=1;
1088  fprintf(stdoutMPI, cErrIncorrectFormatForKondoTrans, isite1, isite2);
1089  }
1090  }
1091  }
1092  else if(X->iCalcModel==SpinlessFermion || X->iCalcModel==SpinlessFermionGC){
1093  if(isigma1 != 0 || isigma2 !=0){
1094  //Not allowed
1095  fprintf(stderr, cErrNonHermiteTrans, isite1, isigma1, isite2, isigma2, dvalue_re, dvalue_im);
1096  fclose(fp);
1097  return ReadDefFileError(defname);
1098  }
1099  }
1100 
1101  iflg_trans=0;
1102  for( i=0; i < icnt_trans; i++){
1103  if(isite1 ==X->GeneralTransfer[i][0] && isite2 == X->GeneralTransfer[i][2]
1104  && isigma1 == X->GeneralTransfer[i][1] && isigma2 == X->GeneralTransfer[i][3])
1105  {
1106  X->ParaGeneralTransfer[i] += dvalue_re+dvalue_im*I;
1107  iflg_trans=1;
1108  continue;
1109  }
1110  }
1111 
1112  if(iflg_trans == 0){
1113  X->GeneralTransfer[icnt_trans][0]=isite1;
1114  X->GeneralTransfer[icnt_trans][1]=isigma1;
1115  X->GeneralTransfer[icnt_trans][2]=isite2;
1116  X->GeneralTransfer[icnt_trans][3]=isigma2;
1117  X->ParaGeneralTransfer[icnt_trans] = dvalue_re+dvalue_im*I;
1118  icnt_trans++;
1119  }
1120  idx++;
1121  }
1122 
1123  if(iboolLoc ==1){
1124  fclose(fp);
1125  return(-1);
1126  }
1127  }
1128 
1129  X->NTransfer = icnt_trans;
1130 
1132  fclose(fp);
1133  return(-1);
1134  }
1135 
1136  if(CheckTransferHermite(X) !=0){
1137  fprintf(stdoutMPI, "%s", cErrNonHermiteTransForAll);
1138  fclose(fp);
1139  return(-1);
1140  }
1141  break;
1142 
1143  case KWCoulombIntra:
1144  /*coulombintra.def----------------------------------*/
1145  if(X->NCoulombIntra>0){
1146  while(fgetsMPI(ctmp2, 256, fp) != NULL){
1147  if(idx==X->NCoulombIntra){
1148  fclose(fp);
1149  return ReadDefFileError(defname);
1150  }
1151  sscanf(ctmp2, "%d %lf\n",
1152  &(X->CoulombIntra[idx][0]),
1153  &(X->ParaCoulombIntra[idx])
1154  );
1155 
1156  if(CheckSite(X->CoulombIntra[idx][0], X->Nsite) !=0){
1157  fclose(fp);
1158  return ReadDefFileError(defname);
1159  }
1160  idx++;
1161  }
1162  }
1163  break;
1164 
1165  case KWCoulombInter:
1166  /*coulombinter.def----------------------------------*/
1167  if(X->NCoulombInter>0){
1168  while(fgetsMPI(ctmp2, 256, fp) != NULL){
1169  if(idx==X->NCoulombInter){
1170  fclose(fp);
1171  return ReadDefFileError(defname);
1172  }
1173 
1174  sscanf(ctmp2, "%d %d %lf\n",
1175  &(X->CoulombInter[idx][0]),
1176  &(X->CoulombInter[idx][1]),
1177  &(X->ParaCoulombInter[idx])
1178  );
1179 
1180  if(CheckPairSite(X->CoulombInter[idx][0], X->CoulombInter[idx][1],X->Nsite) !=0){
1181  fclose(fp);
1182  return ReadDefFileError(defname);
1183  }
1184 
1185  idx++;
1186  }
1187  }
1188  break;
1189 
1190  case KWHund:
1191  /*hund.def------------------------------------------*/
1192  if(X->NHundCoupling>0){
1193  while(fgetsMPI(ctmp2,256,fp) != NULL)
1194  {
1195  if(idx==X->NHundCoupling){
1196  fclose(fp);
1197  return ReadDefFileError(defname);
1198  }
1199 
1200  sscanf(ctmp2, "%d %d %lf\n",
1201  &(X->HundCoupling[idx][0]),
1202  &(X->HundCoupling[idx][1]),
1203  &(X->ParaHundCoupling[idx])
1204  );
1205 
1206  if(CheckPairSite(X->HundCoupling[idx][0], X->HundCoupling[idx][1],X->Nsite) !=0){
1207  fclose(fp);
1208  return ReadDefFileError(defname);
1209  }
1210 
1211  idx++;
1212  }
1213  }
1214  break;
1215  case KWPairHop:
1216  /*pairhop.def---------------------------------------*/
1217  if(X->iCalcModel == Spin || X->iCalcModel == SpinGC){
1218  fprintf(stdoutMPI, "PairHop is not active in Spin and SpinGC.\n");
1219  return(-1);
1220  }
1221 
1222  if(X->NPairHopping>0){
1223  while(fgetsMPI(ctmp2, 256, fp) != NULL){
1224  if(idx==X->NPairHopping/2){
1225  fclose(fp);
1226  return ReadDefFileError(defname);
1227  }
1228  sscanf(ctmp2, "%d %d %lf\n",
1229  &(X->PairHopping[2*idx][0]),
1230  &(X->PairHopping[2*idx][1]),
1231  &(X->ParaPairHopping[2*idx])
1232  );
1233 
1234  if(CheckPairSite(X->PairHopping[2*idx][0], X->PairHopping[2*idx][1],X->Nsite) !=0){
1235  fclose(fp);
1236  return ReadDefFileError(defname);
1237  }
1238  X->PairHopping[2*idx+1][0]=X->PairHopping[2*idx][1];
1239  X->PairHopping[2*idx+1][1]=X->PairHopping[2*idx][0];
1240  X->ParaPairHopping[2*idx+1]=X->ParaPairHopping[2*idx];
1241  idx++;
1242  }
1243  }
1244  break;
1245 
1246  case KWExchange:
1247  /*exchange.def--------------------------------------*/
1248  if(X->NExchangeCoupling>0){
1249  while(fgetsMPI(ctmp2,256,fp) != NULL){
1250  if(idx==X->NExchangeCoupling){
1251  fclose(fp);
1252  return ReadDefFileError(defname);
1253  }
1254 
1255  sscanf(ctmp2, "%d %d %lf\n",
1256  &(X->ExchangeCoupling[idx][0]),
1257  &(X->ExchangeCoupling[idx][1]),
1258  &(X->ParaExchangeCoupling[idx])
1259  );
1260 
1261  if(CheckPairSite(X->ExchangeCoupling[idx][0], X->ExchangeCoupling[idx][1],X->Nsite) !=0){
1262  fclose(fp);
1263  return ReadDefFileError(defname);
1264  }
1265 
1266  idx++;
1267  }
1268  }
1269  break;
1270 
1271  case KWIsing:
1272  /*ising.def--------------------------------------*/
1273  if(X->NIsingCoupling>0){
1274  while(fgetsMPI(ctmp2,256,fp) != NULL){
1275  if(idx==X->NIsingCoupling){
1276  fclose(fp);
1277  return ReadDefFileError(defname);
1278  }
1279 
1280  sscanf(ctmp2, "%d %d %lf\n",
1281  &isite1,
1282  &isite2,
1283  &dvalue_re
1284  );
1285 
1286  if(CheckPairSite(isite1,isite2,X->Nsite) !=0){
1287  fclose(fp);
1288  return ReadDefFileError(defname);
1289  }
1290 
1291  //input into exchange couplings
1292  X->HundCoupling[X->NHundCoupling+idx][0]=isite1;
1293  X->HundCoupling[X->NHundCoupling+idx][1]=isite2;
1294  X->ParaHundCoupling[X->NHundCoupling+idx]= -dvalue_re/2.0;
1295  //input into inter Coulomb
1296  X->CoulombInter[X->NCoulombInter+idx][0]=isite1;
1297  X->CoulombInter[X->NCoulombInter+idx][1]=isite2;
1298  X->ParaCoulombInter[X->NCoulombInter+idx]=-dvalue_re/4.0;
1299  idx++;
1300  }
1301  }
1302  break;
1303 
1304  case KWPairLift:
1305  /*pairlift.def--------------------------------------*/
1306  if(X->NPairLiftCoupling>0){
1307  if(X->iCalcModel != SpinGC){
1308  fprintf(stdoutMPI, "PairLift is active only in SpinGC.\n");
1309  return(-1);
1310  }
1311  while(fgetsMPI(ctmp2,256,fp) != NULL)
1312  {
1313  if(idx==X->NPairLiftCoupling){
1314  fclose(fp);
1315  return ReadDefFileError(defname);
1316  }
1317 
1318  sscanf(ctmp2, "%d %d %lf\n",
1319  &(X->PairLiftCoupling[idx][0]),
1320  &(X->PairLiftCoupling[idx][1]),
1321  &(X->ParaPairLiftCoupling[idx])
1322  );
1323 
1324  if(CheckPairSite(X->PairLiftCoupling[idx][0], X->PairLiftCoupling[idx][1],X->Nsite) !=0){
1325  fclose(fp);
1326  return ReadDefFileError(defname);
1327  }
1328 
1329  idx++;
1330  }
1331  }
1332  break;
1333 
1334  case KWInterAll:
1335  /*interall.def---------------------------------------*/
1336  X->NInterAll_Diagonal=0;
1337  X->NInterAll_OffDiagonal=0;
1338  if(X->NInterAll>0) {
1339  icnt_interall =0;
1340  icnt_diagonal=0;
1341  while (fgetsMPI(ctmp2, 256, fp) != NULL) {
1342  if (idx == X->NInterAll) {
1343  fclose(fp);
1344  return ReadDefFileError(defname);
1345  }
1346  sscanf(ctmp2, "%d %d %d %d %d %d %d %d %lf %lf\n",
1347  &isite1,
1348  &isigma1,
1349  &isite2,
1350  &isigma2,
1351  &isite3,
1352  &isigma3,
1353  &isite4,
1354  &isigma4,
1355  &dvalue_re,
1356  &dvalue_im
1357  );
1358 
1359  if (CheckInterAllCondition(X->iCalcModel, X->Nsite, X->iFlgGeneralSpin, X->LocSpn,
1360  isite1, isigma1, isite2, isigma2,
1361  isite3, isigma3, isite4, isigma4) != 0) {
1362  fclose(fp);
1363  return ReadDefFileError(defname);
1364  }
1365 
1366  if (InputInterAllInfo(&icnt_interall,
1367  X->InterAll,
1368  X->ParaInterAll,
1369  isite1, isigma1,
1370  isite2, isigma2,
1371  isite3, isigma3,
1372  isite4, isigma4,
1373  dvalue_re, dvalue_im
1374  ) != 0) {
1375  icnt_diagonal += 1;
1376  }
1377  idx++;
1378  }
1379  }
1380 
1381  X->NInterAll = icnt_interall;
1382  X->NInterAll_Diagonal=icnt_diagonal;
1383  X->NInterAll_OffDiagonal = X->NInterAll-X->NInterAll_Diagonal;
1384 
1385 /*
1386  setmem_IntAll_Diagonal(
1387  X->InterAll_OffDiagonal, X->ParaInterAll_OffDiagonal,
1388  X->InterAll_Diagonal, X->ParaInterAll_Diagonal, NInterAllSet);
1389 */
1391  X->InterAll, X->ParaInterAll, X->NInterAll,
1392  X->InterAll_Diagonal, X->ParaInterAll_Diagonal,
1393  X->InterAll_OffDiagonal, X->ParaInterAll_OffDiagonal,
1394  X->EDChemi, X->EDSpinChemi, X->EDParaChemi, &X->EDNChemi,
1395  X->iCalcModel
1396  )!=0){
1397  fclose(fp);
1398  return(-1);
1399  }
1400 
1402  X->InterAll, X->ParaInterAll,
1403  X->InterAll_OffDiagonal, X->ParaInterAll_OffDiagonal,
1404  X->NInterAll_OffDiagonal, X->iCalcModel
1405  )!=0) {
1406  fprintf(stdoutMPI, "%s", cErrNonHermiteInterAllForAll);
1407  fclose(fp);
1408  return (-1);
1409  }
1410  break;
1411 
1412  case KWOneBodyG:
1413  /*cisajs.def----------------------------------------*/
1414  if(X->NCisAjt>0){
1415  while(fgetsMPI(ctmp2, 256, fp) != NULL){
1416  if(idx==X->NCisAjt){
1417  fclose(fp);
1418  return ReadDefFileError(defname);
1419  }
1420  sscanf(ctmp2, "%d %d %d %d\n",
1421  &isite1,
1422  &isigma1,
1423  &isite2,
1424  &isigma2);
1425 
1426  if(X->iCalcModel == Spin){
1427  if(isite1 != isite2){
1428  fprintf(stdoutMPI, cWarningIncorrectFormatForSpin2, isite1, isite2);
1429  X->NCisAjt--;
1430  continue;
1431  }
1432  }
1433 
1434  X->CisAjt[ idx ][0] = isite1;
1435  X->CisAjt[ idx ][1] = isigma1;
1436  X->CisAjt[ idx ][2] = isite2;
1437  X->CisAjt[ idx ][3] = isigma2;
1438 
1439  if(CheckPairSite(isite1, isite2,X->Nsite) !=0){
1440  fclose(fp);
1441  return ReadDefFileError(defname);
1442  }
1443 
1444  idx++;
1445  }
1446  }
1447  break;
1448 
1449  case KWTwoBodyG:
1450  /*cisajscktaltdc.def--------------------------------*/
1451  if(X->NCisAjtCkuAlvDC>0){
1452  while(fgetsMPI(ctmp2, 256, fp) != NULL){
1453  if(idx==X->NCisAjtCkuAlvDC){
1454  fclose(fp);
1455  return ReadDefFileError(defname);
1456  }
1457 
1458  sscanf(ctmp2, "%d %d %d %d %d %d %d %d\n",
1459  &isite1,
1460  &isigma1,
1461  &isite2,
1462  &isigma2,
1463  &isite3,
1464  &isigma3,
1465  &isite4,
1466  &isigma4
1467  );
1468 
1469  if(X->iCalcModel == Spin || X->iCalcModel == SpinGC){
1470  if(CheckFormatForSpinInt(isite1, isite2, isite3, isite4)!=0){
1471  exitMPI(-1);
1472  //X->NCisAjtCkuAlvDC--;
1473  //continue;
1474  }
1475  }
1476 
1477 
1478  X->CisAjtCkuAlvDC[idx][0] = isite1;
1479  X->CisAjtCkuAlvDC[idx][1] = isigma1;
1480  X->CisAjtCkuAlvDC[idx][2] = isite2;
1481  X->CisAjtCkuAlvDC[idx][3] = isigma2;
1482  X->CisAjtCkuAlvDC[idx][4] = isite3;
1483  X->CisAjtCkuAlvDC[idx][5] = isigma3;
1484  X->CisAjtCkuAlvDC[idx][6] = isite4;
1485  X->CisAjtCkuAlvDC[idx][7] = isigma4;
1486 
1487  if(CheckQuadSite(isite1, isite2, isite3, isite4,X->Nsite) !=0){
1488  fclose(fp);
1489  return ReadDefFileError(defname);
1490  }
1491  idx++;
1492  }
1493  }
1494  break;
1495 
1496  case KWLaser:
1497  //printf("KWLaser\n");
1498  /*laser.def----------------------------------*/
1499  if(X->NLaser>0){
1500  //printf("Read Start\n");
1501  while(fgetsMPI(ctmp2, 256, fp) != NULL){
1502  sscanf(ctmp2, "%s %lf\n", &(ctmp), &(X->ParaLaser[idx]));
1503  //printf("[%d]:%f\n",idx,X->ParaLaser[idx]);
1504  idx++;
1505  }
1506  if(idx!=X->NLaser){
1507  fclose(fp);
1508  return ReadDefFileError(defname);
1509  }
1510  }
1511  break;
1512 
1513  case KWTEOneBody:
1514  if(X->NTETimeSteps>0){
1515  idx=0;
1516  while(fgetsMPI(ctmp2, 256, fp) != NULL){
1517  sscanf(ctmp2, "%lf %d\n", &(X->TETime[idx]), &(X->NTETransfer[idx]));
1518  for(i=0; i<X->NTETransfer[idx]; ++i ){
1519  fgetsMPI(ctmp2, 256, fp);
1520  sscanf(ctmp2, "%d %d %d %d %lf %lf\n",
1521  &isite1,
1522  &isigma1,
1523  &isite2,
1524  &isigma2,
1525  &dvalue_re,
1526  &dvalue_im
1527  );
1528  X->TETransfer[idx][i][0]= isite1;
1529  X->TETransfer[idx][i][1]= isigma1;
1530  X->TETransfer[idx][i][2]= isite2;
1531  X->TETransfer[idx][i][3] = isigma2;
1532  X->ParaTETransfer[idx][i]=dvalue_re+dvalue_im*I;
1533  }
1534  //check Transfer Hermite
1535  if(CheckTETransferHermite(X, X->NTETransfer[idx], idx)!=0){
1536  fclose(fp);
1537  return ReadDefFileError(defname);
1538  }
1539  idx++;
1540  }
1541  if(idx!=X->NTETimeSteps){
1542  fclose(fp);
1543  return ReadDefFileError(defname);
1544  }
1545  }
1546  break;
1547 
1548  case KWTETwoBody:
1549  if(X->NTETimeSteps>0){
1550  idx=0;
1551  while(fgetsMPI(ctmp2, 256, fp) != NULL) {
1552  sscanf(ctmp2, "%lf %d\n", &(X->TETime[idx]), &(X->NTEInterAll[idx]));
1553  icnt_interall =0;
1554  icnt_diagonal=0;
1555  for (i = 0; i < X->NTEInterAll[idx]; ++i) {
1556  fgetsMPI(ctmp2, 256, fp);
1557  sscanf(ctmp2, "%d %d %d %d %d %d %d %d %lf %lf\n",
1558  &isite1,
1559  &isigma1,
1560  &isite2,
1561  &isigma2,
1562  &isite3,
1563  &isigma3,
1564  &isite4,
1565  &isigma4,
1566  &dvalue_re,
1567  &dvalue_im
1568  );
1569  if (CheckInterAllCondition(X->iCalcModel, X->Nsite, X->iFlgGeneralSpin, X->LocSpn,
1570  isite1, isigma1, isite2, isigma2,
1571  isite3, isigma3, isite4, isigma4) != 0) {
1572  fclose(fp);
1573  return ReadDefFileError(defname);
1574  }
1575  if (InputInterAllInfo(&icnt_interall,
1576  X->TEInterAll[idx],
1577  X->ParaTEInterAll[idx],
1578  isite1, isigma1,
1579  isite2, isigma2,
1580  isite3, isigma3,
1581  isite4, isigma4,
1582  dvalue_re, dvalue_im
1583  ) != 0) {
1584  icnt_diagonal += 1;
1585  }
1586  }
1587 
1588  X->NTEInterAll[idx] = icnt_interall;
1589  X->NTEInterAllDiagonal[idx] = icnt_diagonal;
1590  X->NTEInterAllOffDiagonal[idx] = icnt_interall - icnt_diagonal;
1591  //Diagonal -> OffDiagonal -> search pair -> hermite
1592  if (GetDiagonalInterAll(X->TEInterAll[idx], X->ParaTEInterAll[idx], X->NTEInterAll[idx], X->TEInterAllDiagonal[idx], X->ParaTEInterAllDiagonal[idx],
1593  X->TEInterAllOffDiagonal[idx], X->ParaTEInterAllOffDiagonal[idx], X->TEChemi[idx], X->SpinTEChemi[idx], X->ParaTEChemi[idx], &X->NTEChemi[idx], X->iCalcModel) != 0)
1594  {
1595  fclose(fp);
1596  return (-1);
1597  }
1598 
1600  X->TEInterAll[idx], X->ParaTEInterAll[idx],
1601  X->TEInterAllOffDiagonal[idx], X->ParaTEInterAllOffDiagonal[idx],
1602  X->NTEInterAllOffDiagonal[idx], X->iCalcModel
1603  )!=0) {
1604  fprintf(stdoutMPI, "%s", cErrNonHermiteInterAllForAll);
1605  fclose(fp);
1606  return (-1);
1607  }
1608  idx++;
1609  }
1610 
1611  if(idx!=X->NTETimeSteps){
1612  fclose(fp);
1613  return ReadDefFileError(defname);
1614  }
1615  }
1616  break;
1617 
1618  case KWBoost:
1619  /* boost.def--------------------------------*/
1620  //input magnetic field
1621  fgetsMPI(ctmp2, 256, fp);
1622  sscanf(ctmp2, "%lf %lf %lf\n",
1623  &dArrayValue_re[0],
1624  &dArrayValue_re[1],
1625  &dArrayValue_re[2]);
1626  for(iline=0; iline<3; iline++){
1627  xBoost->vecB[iline]= dArrayValue_re[iline];
1628  }
1629 
1630  //this line is skipped;
1631  fgetsMPI(ctmp2, 256, fp);
1632 
1633  //input arrayJ
1634  if(xBoost->NumarrayJ>0){
1635  for(iline=0; iline<xBoost->NumarrayJ; iline++){
1636  for(ilineIn=0; ilineIn<3; ilineIn++){
1637  fgetsMPI(ctmp2, 256, fp);
1638  sscanf(ctmp2, "%lf %lf %lf\n",
1639  &dArrayValue_re[0],
1640  &dArrayValue_re[1],
1641  &dArrayValue_re[2]);
1642  for(ilineIn2=0; ilineIn2<3; ilineIn2++){
1643  xBoost->arrayJ[iline][ilineIn][ilineIn2]= dArrayValue_re[ilineIn2];
1644  }
1645  }
1646  }
1647  }
1648 
1649  //this line is skipped;
1650  fgetsMPI(ctmp2, 256, fp);
1651 
1652  //read list_6spin_star
1653  if(xBoost->num_pivot>0){
1654  for(iline=0; iline<xBoost->num_pivot; iline++){
1655  //input
1656  fgetsMPI(ctmp2, 256, fp);
1657  sscanf(ctmp2, "%d %d %d %d %d %d %d\n",
1658  &xBoost->list_6spin_star[iline][0],
1659  &xBoost->list_6spin_star[iline][1],
1660  &xBoost->list_6spin_star[iline][2],
1661  &xBoost->list_6spin_star[iline][3],
1662  &xBoost->list_6spin_star[iline][4],
1663  &xBoost->list_6spin_star[iline][5],
1664  &xBoost->list_6spin_star[iline][6]
1665  );
1666  //copy
1667  for(iloop=0; iloop<xBoost->R0; iloop++){
1668  for(itmp=0; itmp<7; itmp++){
1669  xBoost->list_6spin_star[iloop*xBoost->num_pivot+iline][itmp]=xBoost->list_6spin_star[iline][itmp];
1670  }
1671  }
1672  }
1673  }
1674 
1675  //read list_6spin_pair
1676  if(xBoost->num_pivot>0){
1677  for(iline=0; iline<xBoost->num_pivot; iline++){
1678  //input
1679  for(ilineIn2=0; ilineIn2<xBoost->list_6spin_star[iline][0]; ilineIn2++){
1680  fgetsMPI(ctmp2, 256, fp);
1681  sscanf(ctmp2, "%d %d %d %d %d %d %d\n",
1682  &xBoost->list_6spin_pair[iline][0][ilineIn2],
1683  &xBoost->list_6spin_pair[iline][1][ilineIn2],
1684  &xBoost->list_6spin_pair[iline][2][ilineIn2],
1685  &xBoost->list_6spin_pair[iline][3][ilineIn2],
1686  &xBoost->list_6spin_pair[iline][4][ilineIn2],
1687  &xBoost->list_6spin_pair[iline][5][ilineIn2],
1688  &xBoost->list_6spin_pair[iline][6][ilineIn2]
1689  );
1690 
1691  //copy
1692  for(iloop=0; iloop<xBoost->R0; iloop++){
1693  for(itmp=0; itmp<7; itmp++){
1694  xBoost->list_6spin_pair[iloop*xBoost->num_pivot+iline][itmp][ilineIn2]=xBoost->list_6spin_pair[iline][itmp][ilineIn2];
1695  }
1696  }
1697  }
1698  }
1699 
1700  }
1701 
1702  break;
1703 
1704  case KWSingleExcitation:
1705  /*singleexcitation.def----------------------------------------*/
1706  if(X->NSingleExcitationOperator>0) {
1707  if(X->iCalcModel == Spin || X->iCalcModel == SpinGC) {
1708  fprintf(stderr, "SingleExcitation is not allowed for spin system.\n");
1709  fclose(fp);
1710  return ReadDefFileError(defname);
1711  }
1712  while (fgetsMPI(ctmp2, 256, fp) != NULL) {
1713  sscanf(ctmp2, "%d %d %d %lf %lf\n",
1714  &isite1,
1715  &isigma1,
1716  &itype,
1717  &dvalue_re,
1718  &dvalue_im
1719  );
1720 
1721  if (CheckSite(isite1, X->Nsite) != 0) {
1722  fclose(fp);
1723  return ReadDefFileError(defname);
1724  }
1725 
1726  X->SingleExcitationOperator[idx][0] = isite1;
1727  X->SingleExcitationOperator[idx][1] = isigma1;
1728  X->SingleExcitationOperator[idx][2] = itype;
1729  X->ParaSingleExcitationOperator[idx] = dvalue_re + I * dvalue_im;
1730  idx++;
1731  }
1732  if (idx != X->NSingleExcitationOperator) {
1733  fclose(fp);
1734  return ReadDefFileError(defname);
1735  }
1736  }
1737  break;
1738 
1739  case KWPairExcitation:
1740  /*pairexcitation.def----------------------------------------*/
1741  if(X->NPairExcitationOperator>0) {
1742  while (fgetsMPI(ctmp2, 256, fp) != NULL) {
1743  sscanf(ctmp2, "%d %d %d %d %d %lf %lf\n",
1744  &isite1,
1745  &isigma1,
1746  &isite2,
1747  &isigma2,
1748  &itype,
1749  &dvalue_re,
1750  &dvalue_im
1751  );
1752  if (CheckPairSite(isite1, isite2, X->Nsite) != 0) {
1753  fclose(fp);
1754  return ReadDefFileError(defname);
1755  }
1756 
1757  if(itype==1){
1758  X->PairExcitationOperator[idx][0] = isite1;
1759  X->PairExcitationOperator[idx][1] = isigma1;
1760  X->PairExcitationOperator[idx][2] = isite2;
1761  X->PairExcitationOperator[idx][3] = isigma2;
1762  X->PairExcitationOperator[idx][4] = itype;
1763  X->ParaPairExcitationOperator[idx] = dvalue_re + I * dvalue_im;
1764  }
1765  else{
1766  X->PairExcitationOperator[idx][0] = isite2;
1767  X->PairExcitationOperator[idx][1] = isigma2;
1768  X->PairExcitationOperator[idx][2] = isite1;
1769  X->PairExcitationOperator[idx][3] = isigma1;
1770  X->PairExcitationOperator[idx][4] = itype;
1771  X->ParaPairExcitationOperator[idx] = -(dvalue_re + I * dvalue_im);
1772  }
1773 
1774  idx++;
1775  }
1776  if (idx != X->NPairExcitationOperator) {
1777  fclose(fp);
1778  return ReadDefFileError(defname);
1779  }
1780  }
1781  break;
1782 
1783  default:
1784  break;
1785  }
1786  fclose(fp);
1787 
1788  switch(iKWidx){
1789  case KWCoulombIntra:
1790  case KWCoulombInter:
1791  case KWHund:
1792  case KWPairHop:
1793  case KWExchange:
1794  case KWIsing:
1795  case KWPairLift:
1796  if(X->iFlgGeneralSpin==TRUE){
1797  fprintf(stdoutMPI, "%s", cErrIncorrectFormatInter);
1798  return(-1);
1799  }
1800  break;
1801  default:
1802  break;
1803  }
1804  }
1805 
1807  /*=======================================================================*/
1808  return 0;
1809 }
1810 
1822  const int iSite,
1823  const int iMaxNum
1824  )
1825 {
1826  if(iSite>=iMaxNum) return(-1);
1827  return 0;
1828 }
1829 
1842  const int iSite1,
1843  const int iSite2,
1844  const int iMaxNum
1845  )
1846 {
1847  if(CheckSite(iSite1, iMaxNum)!=0){
1848  return(-1);
1849  }
1850  if(CheckSite(iSite2, iMaxNum)!=0){
1851  return(-1);
1852  }
1853  return 0;
1854 }
1855 
1870  const int iSite1,
1871  const int iSite2,
1872  const int iSite3,
1873  const int iSite4,
1874  const int iMaxNum
1875  )
1876 {
1877  if(CheckPairSite(iSite1, iSite2, iMaxNum)!=0){
1878  return(-1);
1879  }
1880  if(CheckPairSite(iSite3, iSite4, iMaxNum)!=0){
1881  return(-1);
1882  }
1883  return 0;
1884 }
1885 
1899  struct DefineList *X
1900  )
1901 {
1902  unsigned int i,j;
1903  int isite1, isite2;
1904  int isigma1, isigma2;
1905  int itmpsite1, itmpsite2;
1906  int itmpsigma1, itmpsigma2;
1907  int itmperrsite1, itmperrsite2;
1908  int itmperrsigma1, itmperrsigma2;
1909  double complex dcerrTrans;
1910  int icheckHermiteCount=FALSE;
1911  int iCount=0;
1912 
1913  double complex ddiff_trans;
1914  unsigned int itmpIdx, icntHermite, icntchemi;
1915  icntHermite=0;
1916  icntchemi=0;
1917 
1918  for(i=0; i<X->NTransfer; i++){
1919  isite1=X->GeneralTransfer[i][0];
1920  isigma1=X->GeneralTransfer[i][1];
1921  isite2=X->GeneralTransfer[i][2];
1922  isigma2=X->GeneralTransfer[i][3];
1923  icheckHermiteCount=FALSE;
1924  // fprintf(stdoutMPI, "Debug: isite1=%d, isigma1=%d, isite2=%d, isigma2=%d, reTrans=%lf, imTrans = %lf\n",
1925  // isite1, isigma1, isite2, isigma2, creal(X->ParaGeneralTransfer[i]), cimag((X->ParaGeneralTransfer[i])));
1926  for(j=0; j<X->NTransfer; j++){
1927  itmpsite1=X->GeneralTransfer[j][0];
1928  itmpsigma1=X->GeneralTransfer[j][1];
1929  itmpsite2=X->GeneralTransfer[j][2];
1930  itmpsigma2=X->GeneralTransfer[j][3];
1931  if(isite1 == itmpsite2 && isite2 == itmpsite1){
1932  if(isigma1 == itmpsigma2 && isigma2 == itmpsigma1){
1933 
1934  ddiff_trans = X->ParaGeneralTransfer[i]-conj(X->ParaGeneralTransfer[j]);
1935  if(cabs(ddiff_trans) > eps_CheckImag0 ){
1936  itmperrsite1=itmpsite1;
1937  itmperrsigma1=itmpsigma1;
1938  itmperrsite2=itmpsite2;
1939  itmperrsigma2=itmpsigma2;
1940  dcerrTrans=X->ParaGeneralTransfer[j];
1941  fprintf(stdoutMPI, cErrNonHermiteTrans, isite1, isigma1, isite2, isigma2, creal(X->ParaGeneralTransfer[i]), cimag(X->ParaGeneralTransfer[i]));
1942  fprintf(stdoutMPI, cErrNonHermiteTrans, itmperrsite1, itmperrsigma1, itmperrsite2, itmperrsigma2, creal(dcerrTrans), cimag(dcerrTrans));
1943  iCount++;
1944  }
1945  else{
1946  if (icheckHermiteCount == FALSE) {
1947  if(i<=j){
1948  if(2*icntHermite >= X->NTransfer){
1949  fprintf(stderr, "Elements of Transfers are incorrect.\n");
1950  return(-1);
1951  }
1952  if(isite1 !=isite2 || isigma1 !=isigma2){
1953  for(itmpIdx=0; itmpIdx<4; itmpIdx++){
1954  X->EDGeneralTransfer[2*icntHermite][itmpIdx]=X->GeneralTransfer[i][itmpIdx];
1955  X->EDGeneralTransfer[2*icntHermite+1][itmpIdx]=X->GeneralTransfer[j][itmpIdx];
1956  }
1957  X->EDParaGeneralTransfer[2*icntHermite]=X->ParaGeneralTransfer[i];
1958  X->EDParaGeneralTransfer[2*icntHermite+1]=X->ParaGeneralTransfer[j];
1959  icntHermite++;
1960  }
1961  else{
1962  X->EDChemi[icntchemi] = X->GeneralTransfer[i][0];
1963  X->EDSpinChemi[icntchemi] = X->GeneralTransfer[i][1];
1964  X->EDParaChemi[icntchemi] = creal(X->ParaGeneralTransfer[i]);
1965  icntchemi+=1;
1966  }
1967  }
1968  icheckHermiteCount = TRUE;
1969  }
1970  }
1971  }
1972 
1973  }
1974  }
1975 
1976  //if counterpart for satisfying hermite conjugate does not exist.
1977  if(icheckHermiteCount == FALSE){
1978  fprintf(stdoutMPI, cErrNonHermiteTrans, isite1, isigma1, isite2, isigma2, creal(X->ParaGeneralTransfer[i]), cimag(X->ParaGeneralTransfer[i]));
1979  iCount++;
1980  }
1981  }
1982 
1983  if(iCount !=0){
1984  return -1;
1985  }
1986  X->EDNTransfer=2*icntHermite;
1987  X->EDNChemi=icntchemi;
1988 
1989  //To realize ido-san's result
1990  for(i=0; i<X->EDNTransfer; i++){
1991  for(itmpIdx=0; itmpIdx<4; itmpIdx++){
1992  X->GeneralTransfer[i][itmpIdx]=X->EDGeneralTransfer[i][itmpIdx];
1993  }
1994  X->ParaGeneralTransfer[i]=X->EDParaGeneralTransfer[i];
1995  }
1996 
1997 
1998  return 0;
1999 }
2000 
2001 
2019  (
2020  int **InterAll,
2021  double complex* ParaInterAll,
2022  int **InterAllOffDiagonal,
2023  double complex*ParaInterAllOffDiagonal,
2024  const int NInterAllOffDiagonal,
2025  const int iCalcModel
2026  ) {
2027  unsigned int i, j, icntincorrect, itmpret;
2028  int isite1, isite2, isite3, isite4;
2029  int isigma1, isigma2, isigma3, isigma4;
2030  int itmpsite1, itmpsite2, itmpsite3, itmpsite4;
2031  int itmpsigma1, itmpsigma2, itmpsigma3, itmpsigma4;
2032  unsigned int itmpIdx, icntHermite;
2033  int icheckHermiteCount = FALSE;
2034  double complex ddiff_intall;
2035  icntincorrect = 0;
2036  icntHermite = 0;
2037  for (i = 0; i < NInterAllOffDiagonal; i++) {
2038  itmpret = 0;
2039  isite1 = InterAllOffDiagonal[i][0];
2040  isigma1 = InterAllOffDiagonal[i][1];
2041  isite2 = InterAllOffDiagonal[i][2];
2042  isigma2 = InterAllOffDiagonal[i][3];
2043  isite3 = InterAllOffDiagonal[i][4];
2044  isigma3 = InterAllOffDiagonal[i][5];
2045  isite4 = InterAllOffDiagonal[i][6];
2046  isigma4 = InterAllOffDiagonal[i][7];
2047  icheckHermiteCount = FALSE;
2048 
2049  for (j = 0; j < NInterAllOffDiagonal; j++) {
2050  itmpsite1 = InterAllOffDiagonal[j][0];
2051  itmpsigma1 = InterAllOffDiagonal[j][1];
2052  itmpsite2 = InterAllOffDiagonal[j][2];
2053  itmpsigma2 = InterAllOffDiagonal[j][3];
2054  itmpsite3 = InterAllOffDiagonal[j][4];
2055  itmpsigma3 = InterAllOffDiagonal[j][5];
2056  itmpsite4 = InterAllOffDiagonal[j][6];
2057  itmpsigma4 = InterAllOffDiagonal[j][7];
2058 
2059  if (isite1 == itmpsite4 && isite2 == itmpsite3 && isite3 == itmpsite2 && isite4 == itmpsite1) {
2060  if (isigma1 == itmpsigma4 && isigma2 == itmpsigma3 && isigma3 == itmpsigma2 && isigma4 == itmpsigma1) {
2061  ddiff_intall = cabs(ParaInterAllOffDiagonal[i] - conj(ParaInterAllOffDiagonal[j]));
2062 
2063  if (cabs(ddiff_intall) < eps_CheckImag0) {
2064  itmpret = 1;
2065  if (icheckHermiteCount == FALSE) {
2066  icheckHermiteCount = TRUE; //for not double counting
2067  if (i <= j) {
2068  if (2 * icntHermite >= NInterAllOffDiagonal) {
2069  fprintf(stdoutMPI, "Elements of InterAll are incorrect.\n");
2070  return (-1);
2071  }
2072 
2073  for (itmpIdx = 0; itmpIdx < 8; itmpIdx++) {
2074  InterAll[2 * icntHermite][itmpIdx] = InterAllOffDiagonal[i][itmpIdx];
2075  InterAll[2 * icntHermite + 1][itmpIdx] = InterAllOffDiagonal[j][itmpIdx];
2076  }
2077 
2078  ParaInterAll[2 * icntHermite] = ParaInterAllOffDiagonal[i];
2079  ParaInterAll[2 * icntHermite + 1] = ParaInterAllOffDiagonal[j];
2080  icntHermite++;
2081  }
2082  break;
2083  }
2084  }
2085  }
2086  } else if (isite1 == itmpsite2 && isite2 == itmpsite1 && isite3 == itmpsite4 &&
2087  isite4 == itmpsite3) { //for spin and Kondo
2088  if (iCalcModel == Kondo || iCalcModel == KondoGC || iCalcModel == Spin || iCalcModel == SpinGC) {
2089  if (isigma1 == itmpsigma2 && isigma2 == itmpsigma1 && isigma3 == itmpsigma4 && isigma4 == itmpsigma3) {
2090  ddiff_intall = ParaInterAllOffDiagonal[i] - conj(ParaInterAllOffDiagonal[j]);
2091  if (cabs(ddiff_intall) < eps_CheckImag0) {
2092  itmpret = 1;
2093  if (icheckHermiteCount == FALSE) {
2094  icheckHermiteCount = TRUE; // for not double-counting
2095  if (i <= j) {
2096  if (2 * icntHermite >= NInterAllOffDiagonal) {
2097  fprintf(stdoutMPI, "Elements of InterAll are incorrect.\n");
2098  return (-1);
2099  }
2100  for (itmpIdx = 0; itmpIdx < 8; itmpIdx++) {
2101  InterAll[2 * icntHermite][itmpIdx] = InterAllOffDiagonal[i][itmpIdx];
2102  }
2103  for (itmpIdx = 0; itmpIdx < 4; itmpIdx++) {
2104  InterAll[2 * icntHermite + 1][2 * itmpIdx] = InterAllOffDiagonal[i][6 -
2105  2 *
2106  itmpIdx];
2107  InterAll[2 * icntHermite + 1][2 * itmpIdx + 1] = InterAllOffDiagonal[i][7 - 2 *
2108  itmpIdx];
2109 
2110  }
2111  ParaInterAll[2 * icntHermite] = ParaInterAllOffDiagonal[i];
2112  ParaInterAll[2 * icntHermite + 1] = ParaInterAllOffDiagonal[j];
2113  icntHermite++;
2114  }
2115  break;
2116  }
2117  }
2118  }
2119  }
2120  }
2121  }
2122  //if counterpart for satisfying hermite conjugate does not exist.
2123  if (itmpret != 1) {
2124  fprintf(stdoutMPI, cErrNonHermiteInterAll, isite1, isigma1, isite2, isigma2, isite3, isigma3, isite4, isigma4,
2125  creal(ParaInterAllOffDiagonal[i]), cimag(ParaInterAllOffDiagonal[i]));
2126  icntincorrect++;
2127  }
2128  }
2129 
2130  if (icntincorrect != 0) {
2131  return (-1);
2132  }
2133 
2134  for (i = 0; i < NInterAllOffDiagonal; i++) {
2135  for (itmpIdx = 0; itmpIdx < 8; itmpIdx++) {
2136  InterAllOffDiagonal[i][itmpIdx] = InterAll[i][itmpIdx];
2137  }
2138  ParaInterAllOffDiagonal[i] = ParaInterAll[i];
2139  }
2140 
2141  return 0;
2142 }
2143 
2162  (
2163  int **InterAll,
2164  complex double *ParaInterAll,
2165  const int NInterAll,
2166  int **InterAllDiagonal,
2167  double *ParaInterAllDiagonal,
2168  int **InterAllOffDiagonal,
2169  complex double *ParaInterAllOffDiagonal,
2170  int *Chemi,
2171  int *SpinChemi,
2172  double *ParaChemi,
2173  unsigned int *NChemi,
2174  const int iCalcModel
2175  )
2176 {
2177  unsigned int i,icnt_diagonal, icnt_offdiagonal, tmp_i;
2178  int isite1, isite2, isite3, isite4;
2179  int isigma1, isigma2, isigma3, isigma4;
2180  int iret=0;
2181  icnt_diagonal=0;
2182  icnt_offdiagonal=0;
2183 
2184  for(i=0; i<NInterAll; i++){
2185  isite1=InterAll[i][0];
2186  isigma1=InterAll[i][1];
2187  isite2=InterAll[i][2];
2188  isigma2=InterAll[i][3];
2189  isite3=InterAll[i][4];
2190  isigma3=InterAll[i][5];
2191  isite4=InterAll[i][6];
2192  isigma4=InterAll[i][7];
2193 
2194  //Get Diagonal term
2195  if(isite1 == isite2 && isite3 == isite4 &&
2196  isigma1 == isigma2 && isigma3 == isigma4)
2197  {
2198  InterAllDiagonal[icnt_diagonal][0]=isite1;
2199  InterAllDiagonal[icnt_diagonal][1]=isigma1;
2200  InterAllDiagonal[icnt_diagonal][2]=isite3;
2201  InterAllDiagonal[icnt_diagonal][3]=isigma3;
2202  ParaInterAllDiagonal[icnt_diagonal] = creal(ParaInterAll[i]);
2203  icnt_diagonal++;
2204  continue;
2205  }
2206  else if(isite1 == isite4 && isite2 ==isite3 &&
2207  isigma1 == isigma4 && isigma2 ==isigma3)
2208  {
2209  InterAllDiagonal[icnt_diagonal][0]=isite1;
2210  InterAllDiagonal[icnt_diagonal][1]=isigma1;
2211  InterAllDiagonal[icnt_diagonal][2]=isite2;
2212  InterAllDiagonal[icnt_diagonal][3]=isigma2;
2213  ParaInterAllDiagonal[icnt_diagonal] = -creal(ParaInterAll[i]);
2214  Chemi[*NChemi] = isite1;
2215  SpinChemi[*NChemi] = isigma1;
2216  //transfer integral has minus sign for default setting
2217  ParaChemi[*NChemi] = -creal(ParaInterAll[i]);
2218  icnt_diagonal++;
2219  *NChemi +=1;
2220  continue;
2221  }
2222  else{
2223  //Get Off-Diagonal term
2224  switch(iCalcModel){
2225  case Hubbard:
2226  case HubbardNConserved:
2227  case Kondo:
2228  case KondoGC:
2229  case HubbardGC:
2230  if(isigma1 == isigma2 && isigma3 == isigma4){
2231  for(tmp_i=0; tmp_i<8; tmp_i++){
2232  InterAllOffDiagonal[icnt_offdiagonal][tmp_i]=InterAll[i][tmp_i];
2233  }
2234  ParaInterAllOffDiagonal[icnt_offdiagonal] = ParaInterAll[i];
2235  }
2236  else if(isigma1==isigma4 && isigma2 == isigma3){
2237  InterAllOffDiagonal[icnt_offdiagonal][0]=isite1;
2238  InterAllOffDiagonal[icnt_offdiagonal][1]=isigma1;
2239  InterAllOffDiagonal[icnt_offdiagonal][2]=isite4;
2240  InterAllOffDiagonal[icnt_offdiagonal][3]=isigma1;
2241  InterAllOffDiagonal[icnt_offdiagonal][4]=isite3;
2242  InterAllOffDiagonal[icnt_offdiagonal][5]=isigma2;
2243  InterAllOffDiagonal[icnt_offdiagonal][6]=isite2;
2244  InterAllOffDiagonal[icnt_offdiagonal][7]=isigma2;
2245  ParaInterAllOffDiagonal[icnt_offdiagonal] = -ParaInterAll[i];
2246  }
2247  else{
2248  // Sz symmetry is assumed
2249  if(iCalcModel==Hubbard || iCalcModel==Kondo){
2251  isite1,
2252  isigma1,
2253  isite2,
2254  isigma2,
2255  isite3,
2256  isigma3,
2257  isite4,
2258  isigma4,
2259  creal(ParaInterAll[i]),
2260  cimag(ParaInterAll[i])
2261  );
2262  iret=-1;
2263  }
2264  else{
2265  for(tmp_i=0; tmp_i<8; tmp_i++){
2266  InterAllOffDiagonal[icnt_offdiagonal][tmp_i]=InterAll[i][tmp_i];
2267  }
2268  ParaInterAllOffDiagonal[icnt_offdiagonal] = ParaInterAll[i];
2269  }
2270  }
2271  break;
2272  case Spin:
2273  case SpinGC:
2274  if(isite1 == isite2 && isite3 == isite4){
2275  for(tmp_i=0; tmp_i<8; tmp_i++){
2276  InterAllOffDiagonal[icnt_offdiagonal][tmp_i]=InterAll[i][tmp_i];
2277  }
2278  ParaInterAllOffDiagonal[icnt_offdiagonal] =ParaInterAll[i];
2279  }
2280  break;
2281  default:
2282  return(-1);
2283  }
2284  if(iret != -1){
2285  icnt_offdiagonal++;
2286  }
2287  }
2288 
2289  if(iret !=0){
2290  return(-1);
2291  }
2292  }
2293 
2294  return 0;
2295 }
2296 
2297 
2311 int JudgeDefType
2313  const int argc,
2314  char *argv[],
2315  int *mode
2316  )
2317 {
2318  int ver_maj =
2319 #include "version_major.h"
2320 ;
2321  int ver_min =
2322 #include "version_miner.h"
2323 ;
2324  int ver_pat =
2325 #include "version_patch.h"
2326 ;
2327 
2328  if(argc == 3 &&
2329  (CheckWords(argv[1], "-e") == 0 ||
2330  CheckWords(argv[1], "--Expert") == 0)){
2331  *mode=EXPERT_MODE;
2332  }
2333  else if (argc == 3 &&
2334  (CheckWords(argv[1], "-s") ==0 ||
2335  CheckWords(argv[1], "--Standard") == 0 )){
2336  *mode=STANDARD_MODE;
2337  }
2338  else if (argc == 3 &&
2339  (CheckWords(argv[1], "-sdry") == 0 ||
2340  CheckWords(argv[1], "-s-dry") == 0)
2341  ){
2342  *mode = STANDARD_DRY_MODE;
2343  }
2344  else if (argc >= 2 &&
2345  (CheckWords(argv[1], "-v") == 0
2346  || CheckWords(argv[1], "--version") == 0)
2347  ) {
2348  fprintf(stdoutMPI, "\nHPhi version %d.%d.%d \n\n", ver_maj, ver_min, ver_pat);
2349  exit(-1);
2350  }
2351  else{
2352  /*fprintf(stdoutMPI, cErrArgv, argv[1]);*/
2353  fprintf(stdoutMPI, "\n[Usage] \n");
2354  fprintf(stdoutMPI, "* Expert mode \n");
2355  fprintf(stdoutMPI, " $ HPhi -e {namelist_file} \n");
2356  fprintf(stdoutMPI, "* Standard mode \n");
2357  fprintf(stdoutMPI, " $ HPhi -s {input_file} \n");
2358  fprintf(stdoutMPI, "* Standard DRY mode \n");
2359  fprintf(stdoutMPI, " $ HPhi -sdry {input_file} \n");
2360  fprintf(stdoutMPI, " In this mode, Hphi stops after it generats expert input files. \n");
2361  fprintf(stdoutMPI, "* Print the version \n");
2362  fprintf(stdoutMPI, " $ HPhi -v \n\n");
2363  exit(-1);
2364  }
2365 
2366  return 0;
2367 }
2368 
2385  const int site1,
2386  const int site2,
2387  const int site3,
2388  const int site4
2389  ){
2390  if(site1==site2 && site3==site4){
2391  return 0;
2392  }
2393 
2394  fprintf(stdoutMPI, cWarningIncorrectFormatForSpin, site1, site2, site3, site4);
2395  return(-1);
2396 
2397 }
2398 
2399 
2412  (
2413  const int isite1, const int isite2,
2414  const int isite3, const int isite4,
2415  int* iLocInfo
2416  )
2417 {
2418  if (iLocInfo[isite1] != ITINERANT || iLocInfo[isite2] != ITINERANT) {
2419  if (isite1 != isite2) {
2420  fprintf(stdoutMPI, cErrIncorrectFormatForKondoInt, isite1, isite2, isite3, isite4);
2421  return -1;
2422  }
2423  }
2424  if (iLocInfo[isite3] != ITINERANT || iLocInfo[isite4] != ITINERANT) {
2425  if (isite3 != isite4) {
2426  fprintf(stdoutMPI, cErrIncorrectFormatForKondoInt, isite1, isite2, isite3, isite4);
2427  return -1;
2428  }
2429  }
2430  return 0;
2431 }
2432 
2443  struct DefineList *X
2444  )
2445 {
2446  //In future, convergence facator can be set by a def file.
2447  int neps = -8;
2448  int nepsCG =-8;
2449  int nbisec =-14;
2450  int nEnergy = -12;
2451  int nShiftBeta=8;
2452  int nepsvec12=-14;
2453  eps=pow(10.0, neps);
2454  eps_CG=pow(10.0, nepsCG);
2455  eps_Lanczos = pow(10,-X->LanczosEps);
2456  eps_Energy = pow(10.0, nEnergy);
2457 }
2458 
2470 int CheckLocSpin
2472  struct DefineList *X
2473  )
2474 {
2475 
2476  unsigned int i=0;
2477  switch(X->iCalcModel){
2478  case Hubbard:
2479  case HubbardNConserved:
2480  case HubbardGC:
2481  case SpinlessFermion:
2482  case SpinlessFermionGC:
2483  for(i=0; i<X->Nsite; i++){
2484  if(X->LocSpn[i]!=ITINERANT){
2485  return FALSE;
2486  }
2487  }
2488  break;
2489 
2490  case Kondo:
2491  case KondoGC:
2492  for(i=0; i<X->Nsite; i++){
2493  if(X->LocSpn[i]>LOCSPIN){
2494  X->iFlgGeneralSpin=TRUE;
2495  }
2496  else if(X->LocSpn[i]<ITINERANT){
2497  return FALSE;
2498  }
2499  }
2500  break;
2501 
2502  case Spin:
2503  case SpinGC:
2504  for(i=0; i<X->Nsite; i++){
2505  if(X->LocSpn[i]>LOCSPIN){
2506  X->iFlgGeneralSpin=TRUE;
2507  }
2508  else if(X->LocSpn[i]<LOCSPIN){
2509  return FALSE;
2510  }
2511  }
2512  break;
2513 
2514  default:
2515  return FALSE;
2516  //break;
2517  }
2518 
2519  if(CheckTotal2Sz(X) != TRUE){
2520  return FALSE;
2521  }
2522  return TRUE;
2523 }
2524 
2536  struct DefineList *X
2537  )
2538 {
2539  X->NHundCoupling += X->NIsingCoupling;
2540  X->NCoulombInter += X->NIsingCoupling;
2541 }
2542 
2553  struct DefineList *X
2554  )
2555 {
2556  X->NTransfer=0;
2557  X->NCoulombIntra=0;
2558  X->NCoulombInter=0;
2559  X->NIsingCoupling=0;
2560  X->NPairLiftCoupling=0;
2561  X->NInterAll=0;
2562  X->NCisAjt=0;
2563  X->NCisAjtCkuAlvDC=0;
2564  X->NSingleExcitationOperator=0;
2565  X->NPairExcitationOperator=0;
2566  //[s] Time Evolution
2567  X->NTETimeSteps=0;
2568  X->NLaser=0;
2569  X->NTEInterAll=0;
2570  X->NTETransfer=0;
2571  //[e] Time Evolution
2572 
2573 }
2574 
2575 
2594  const int isite1, const int isigma1,
2595  const int isite2, const int isigma2,
2596  const int isite3, const int isigma3,
2597  const int isite4, const int isigma4,
2598  int* iLocInfo
2599  )
2600 {
2601  if( isigma1 > iLocInfo[isite1] || isigma2 >iLocInfo[isite2]
2602  ||isigma3 > iLocInfo[isite3] || isigma4 >iLocInfo[isite4]){
2603  fprintf(stdoutMPI, "%s", cErrIncorrectSpinIndexForInter);
2604  return FALSE;
2605  }
2606  return TRUE;
2607 }
2608 
2621  struct DefineList *X
2622  )
2623 {
2624  unsigned int i=0;
2625  int isite1, isite2;
2626  int isigma1, isigma2;
2627  if(X->iFlgGeneralSpin==TRUE){
2628  for(i=0; i<X->NTransfer; i++){
2629  isite1 =X->GeneralTransfer[i][0];
2630  isigma1=X->GeneralTransfer[i][1];
2631  isite2 =X->GeneralTransfer[i][2];
2632  isigma2=X->GeneralTransfer[i][3];
2633  if(isigma1 > X->LocSpn[isite1] || isigma2 >X->LocSpn[isite2]){
2634  fprintf(stdoutMPI, "%s", cErrIncorrectSpinIndexForTrans);
2635  return FALSE;
2636  }
2637  }
2638  }
2639  return TRUE;
2640 }
2641 
2652 int CheckTotal2Sz
2654  struct DefineList *X
2655  )
2656 {
2657  if(X->iFlgSzConserved==TRUE && X->iFlgGeneralSpin==FALSE){
2658  int tmp_Nup=X->NLocSpn+X->NCond+X->Total2Sz;
2659  int tmp_Ndown=X->NLocSpn+X->NCond-X->Total2Sz;
2660  if(tmp_Nup%2 != 0 && tmp_Ndown%2 !=0){
2661  printf("Nup=%d, Ndown=%d\n",X->Nup,X->Ndown);
2662  fprintf(stdoutMPI, "2Sz is incorrect.\n");
2663  return FALSE;
2664  }
2665  }
2666  return TRUE;
2667 }
2668 
2681  const char* ctmp,
2682  const char* cKeyWord
2683  )
2684 {
2685  unsigned int i=0;
2686 
2687  char ctmp_small[256]={0};
2688  char cKW_small[256]={0};
2689  unsigned int n;
2690  n=strlen(cKeyWord);
2691  strncpy(cKW_small, cKeyWord, n);
2692 
2693  for(i=0; i<n; i++){
2694  cKW_small[i]=tolower(cKW_small[i]);
2695  }
2696  n=strlen(ctmp);
2697  strncpy(ctmp_small, ctmp, n);
2698  for(i=0; i<n; i++){
2699  ctmp_small[i]=tolower(ctmp_small[i]);
2700  }
2701  if(n<strlen(cKW_small)) n=strlen(cKW_small);
2702  return(strncmp(ctmp_small, cKW_small, n));
2703 }
2704 
2712  int iKWidx,
2713  char **FileName
2714 ){
2715  if(cFileNameListFile == NULL){
2716  return -1;
2717  }
2718  *FileName=cFileNameListFile[iKWidx];
2719  return 0;
2720 }
2721 
2722 
2743  int iCalcModel,
2744  int Nsite,
2745  int iFlgGeneralSpin,
2746  int *iLocInfo,
2747  int isite1, int isigma1,
2748  int isite2, int isigma2,
2749  int isite3, int isigma3,
2750  int isite4, int isigma4
2751 ){
2752  if(CheckQuadSite(isite1, isite2, isite3, isite4, Nsite) !=0){
2753  fprintf(stderr, "%s", "Error: Site index of InterAll is incorrect.\n");
2754  return(-1);
2755  }
2756 
2757  if(iCalcModel == Spin || iCalcModel ==SpinGC){
2758  if(CheckFormatForSpinInt(isite1, isite2, isite3, isite4)!=0){
2759  fprintf(stderr, "%s", "Error: Spin index of InterAll is incorrect.\n");
2760  return(-1);
2761  }
2762  }
2763  else if(iCalcModel == SpinlessFermion || iCalcModel==SpinlessFermionGC){
2764  if(isigma1 !=0 || isigma2 != 0 || isigma3 != 0 || isigma4 !=0){
2765  fprintf(stderr, "%s", "Error: Spin index of InterAll is incorrect.\n");
2766  return -1;
2767  }
2768  }
2769  else if(iCalcModel == Kondo){
2770  if(CheckFormatForKondoInt(isite1, isite2, isite3, isite4, iLocInfo)!=0){
2771  return -1;
2772  }
2773  }
2774 
2775  if(iFlgGeneralSpin ==TRUE) {
2776  if(CheckGeneralSpinIndexForInterAll(isite1, isigma1, isite2, isigma2, isite3, isigma3, isite4, isigma4, iLocInfo)!=TRUE){
2777  return -1;
2778  }
2779  }
2780  return 0;
2781 }
2782 
2801  int *icnt_interall,
2802  int **iInterAllInfo,
2803  double complex *cInterAllValue,
2804  int isite1, int isigma1,
2805  int isite2, int isigma2,
2806  int isite3, int isigma3,
2807  int isite4, int isigma4,
2808  double dvalue_re, double dvalue_im
2809 ) {
2810  int i = 0;
2811  int iflg_interall = 0;
2812  //Collect and sum same components of InterAll interactions
2813  for (i = 0; i < *icnt_interall; i++) {
2814  if (isite1 == iInterAllInfo[i][0] && isite2 == iInterAllInfo[i][2] &&
2815  isite3 == iInterAllInfo[i][4] && isite4 == iInterAllInfo[i][6] &&
2816  isigma1 == iInterAllInfo[i][1] && isigma2 == iInterAllInfo[i][3] &&
2817  isigma3 == iInterAllInfo[i][5] && isigma4 == iInterAllInfo[i][7]) {
2818  cInterAllValue[i] += dvalue_re + dvalue_im * I;
2819  iflg_interall = 1;
2820  return 0;
2821  }
2822  }
2823 
2824  //Input all InterAll interactions
2825  if (iflg_interall == 0) {
2826  iInterAllInfo[*icnt_interall][0] = isite1;
2827  iInterAllInfo[*icnt_interall][1] = isigma1;
2828  iInterAllInfo[*icnt_interall][2] = isite2;
2829  iInterAllInfo[*icnt_interall][3] = isigma2;
2830  iInterAllInfo[*icnt_interall][4] = isite3;
2831  iInterAllInfo[*icnt_interall][5] = isigma3;
2832  iInterAllInfo[*icnt_interall][6] = isite4;
2833  iInterAllInfo[*icnt_interall][7] = isigma4;
2834  cInterAllValue[*icnt_interall] = dvalue_re + I * dvalue_im;
2835  *icnt_interall+=1;
2836  //Check Diagonal part or not
2837  if (isite1 == isite2 && isite3 == isite4 &&
2838  isigma1 == isigma2 && isigma3 == isigma4) { //normal diagonal part
2839  return 1;
2840  } else if (isite1 == isite4 && isite2 == isite3 &&
2841  isigma1 == isigma4 && isigma2 == isigma3) { //hund term
2842  return 1;
2843  }
2844  }
2845  return 0;
2846 }
2847 
2848 
2849 
2859  (
2860  struct DefineList *X,
2861  const int NTETransfer,
2862  const int idx
2863  )
2864 {
2865  unsigned int i,j;
2866  int isite1, isite2;
2867  int isigma1, isigma2;
2868  int itmpsite1, itmpsite2;
2869  int itmpsigma1, itmpsigma2;
2870  int itmperrsite1, itmperrsite2;
2871  int itmperrsigma1, itmperrsigma2;
2872  double complex dcerrTrans;
2873  int icheckHermiteCount;
2874  int iCount=0;
2875 
2876  double complex ddiff_trans;
2877  unsigned int itmpIdx, icntHermite, icntchemi;
2878  icntHermite=0;
2879  icntchemi=0;
2880 
2881  int** tmp_TETransfer = (int**)malloc((NTETransfer)*sizeof(int*));
2882  for(i =0; i<NTETransfer; i++ ){
2883  tmp_TETransfer[i] = (int*)malloc((4*sizeof(int)));
2884  }
2885  double complex*tmp_paraTETransfer = (double complex*)malloc((NTETransfer)*sizeof(double complex));
2886 
2887  //copy
2888  for(i=0; i<NTETransfer; i++){
2889  for(j=0; j<4; j++){
2890  tmp_TETransfer[i][j]=X->TETransfer[idx][i][j];
2891  X->TETransfer[idx][i][j]=0;
2892  }
2893  tmp_paraTETransfer[i] = X->ParaTETransfer[idx][i];
2894  X->ParaTETransfer[idx][i]=0.0;
2895  }
2896 
2897  for(i=0; i<NTETransfer; i++){
2898  isite1=tmp_TETransfer[i][0];
2899  isigma1=tmp_TETransfer[i][1];
2900  isite2=tmp_TETransfer[i][2];
2901  isigma2=tmp_TETransfer[i][3];
2902  icheckHermiteCount=FALSE;
2903  for(j=0; j<NTETransfer; j++){
2904  itmpsite1=tmp_TETransfer[j][0];
2905  itmpsigma1=tmp_TETransfer[j][1];
2906  itmpsite2=tmp_TETransfer[j][2];
2907  itmpsigma2=tmp_TETransfer[j][3];
2908  if(isite1 == itmpsite2 && isite2 == itmpsite1){
2909  if(isigma1 == itmpsigma2 && isigma2 == itmpsigma1){
2910 
2911  ddiff_trans = tmp_paraTETransfer[i]-conj(tmp_paraTETransfer[j]);
2912  if(cabs(ddiff_trans) > eps_CheckImag0 ){
2913  itmperrsite1=itmpsite1;
2914  itmperrsigma1=itmpsigma1;
2915  itmperrsite2=itmpsite2;
2916  itmperrsigma2=itmpsigma2;
2917  dcerrTrans=tmp_paraTETransfer[j];
2918  fprintf(stdoutMPI, cErrNonHermiteTrans, isite1, isigma1, isite2, isigma2, creal(tmp_paraTETransfer[i]), cimag(tmp_paraTETransfer[i]));
2919  fprintf(stdoutMPI, cErrNonHermiteTrans, itmperrsite1, itmperrsigma1, itmperrsite2, itmperrsigma2, creal(dcerrTrans), cimag(dcerrTrans));
2920  iCount++;
2921  }
2922  else{
2923  if (icheckHermiteCount == FALSE) {
2924  if(i<=j){
2925  if(2*icntHermite >= NTETransfer){
2926  fprintf(stderr, "Elements of Transfers are incorrect.\n");
2927  return(-1);
2928  }
2929  if(isite1 !=isite2 || isigma1 !=isigma2){
2930  for(itmpIdx=0; itmpIdx<4; itmpIdx++){
2931  X->TETransfer[idx][2*icntHermite][itmpIdx]=tmp_TETransfer[i][itmpIdx];
2932  X->TETransfer[idx][2*icntHermite+1][itmpIdx]=tmp_TETransfer[j][itmpIdx];
2933  }
2934  X->ParaTETransfer[idx][2*icntHermite]=tmp_paraTETransfer[i];
2935  X->ParaTETransfer[idx][2*icntHermite+1]=tmp_paraTETransfer[j];
2936  icntHermite++;
2937  }
2938  else{
2939  X->TETransferDiagonal[idx][icntchemi][0] = tmp_TETransfer[i][0];
2940  X->TETransferDiagonal[idx][icntchemi][1] = tmp_TETransfer[i][1];
2941  X->ParaTETransferDiagonal[idx][icntchemi] = creal(tmp_paraTETransfer[i]);
2942  icntchemi+=1;
2943  }
2944  }
2945  icheckHermiteCount = TRUE;
2946  }
2947  }
2948  }
2949  }
2950  }
2951 
2952  //if counterpart for satisfying hermite conjugate does not exist.
2953  if(icheckHermiteCount == FALSE){
2954  fprintf(stdoutMPI, cErrNonHermiteTrans, isite1, isigma1, isite2, isigma2, creal(tmp_paraTETransfer[i]), cimag(tmp_paraTETransfer[i]));
2955  iCount++;
2956  //fprintf(stdoutMPI, cErrNonHermiteTrans, itmperrsite1, itmperrsigma1, itmperrsite2, itmperrsigma2, creal(dcerrTrans), cimag(dcerrTrans));
2957  //return(-1);
2958  }
2959  }
2960 
2961  if(iCount !=0){
2962  return -1;
2963  }
2964 
2965  X->NTETransfer[idx]=2*icntHermite;
2966  X->NTETransferDiagonal[idx]=icntchemi;
2967 
2968  free(tmp_TETransfer);
2969  free(tmp_paraTETransfer);
2970  return 0;
2971 }
3104 
3111 
3114 
3117 
const char * cReadFile
Definition: LogMessage.c:28
char * cErrDefFileFormat
Definition: ErrorMessage.c:40
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
int CheckFormatForKondoInt(const int isite1, const int isite2, const int isite3, const int isite4, int *iLocInfo)
function of checking format of Kondo interactions
Definition: readdef.c:2412
char * cErrIncorrectFormatForKondoTrans
Definition: ErrorMessage.c:80
long unsigned int W0
Definition: struct.h:396
char * cErrIncorrectSpinIndexForInter
Definition: ErrorMessage.c:87
int CheckInterAllHermite(int **InterAll, double complex *ParaInterAll, int **InterAllOffDiagonal, double complex *ParaInterAllOffDiagonal, const int NInterAllOffDiagonal, const int iCalcModel)
function of checking hermite conditions about interall interactions
Definition: readdef.c:2019
int JudgeDefType(const int argc, char *argv[], int *mode)
function of judging a type of define files.
Definition: readdef.c:2312
int ** list_6spin_star
Definition: struct.h:402
void SetConvergenceFactor(struct DefineList *X)
function to set convergence factors
Definition: readdef.c:2442
char * cErrReadDefFile
Error Message in readdef.c.
Definition: ErrorMessage.c:39
int GetKWWithIdx(char *ctmpLine, char *ctmp, int *itmp)
Function of Getting keyword and it&#39;s variable from characters.
Definition: readdef.c:187
int GetFileName(const char *cFileListNameFile, char cFileNameList[][D_CharTmpReadDef])
Function of Fitting FileName.
Definition: readdef.c:391
char * cErrNsite
Definition: ErrorMessage.c:58
#define ITINERANT
Definition: global.h:31
void InitializeInteractionNum(struct DefineList *X)
function of initializing interactions
Definition: readdef.c:2552
int CheckLocSpin(struct DefineList *X)
function of checking indexies of localized spin
Definition: readdef.c:2471
char * cErrLanczosExct
Definition: ErrorMessage.c:63
int NumAve
Definition: global.h:65
char * cErrMakeDef
Definition: ErrorMessage.c:72
int CheckInterAllCondition(int iCalcModel, int Nsite, int iFlgGeneralSpin, int *iLocSpin, int isite1, int isigma1, int isite2, int isigma2, int isite3, int isigma3, int isite4, int isigma4)
Check InterAll condition.
Definition: readdef.c:2742
char * cErrIncorrectDef
Definition: ErrorMessage.c:73
char * cErrLanczos_max
Definition: ErrorMessage.c:61
long unsigned int num_pivot
Definition: struct.h:397
int ReadDefFileIdxPara(struct DefineList *X, struct BoostList *xBoost)
function of reading def files to get keyword index
Definition: readdef.c:976
int ReadDefFileError(const char *defname)
Error Function of reading def files.
Definition: readdef.c:112
int *** list_6spin_pair
Definition: struct.h:403
char * cErrKW_ShowList
Definition: ErrorMessage.c:68
char * cErrIncorrectFormatInter
Definition: ErrorMessage.c:86
int ReadDefFileNInt(char *xNameListFile, struct DefineList *X, struct BoostList *xBoost)
Function of reading information about "ModPara" file and total number of parameters from other def fi...
Definition: readdef.c:456
#define TRUE
Definition: global.h:26
int ValidateValue(const int icheckValue, const int ilowestValue, const int iHighestValue)
Function of Validating value.
Definition: readdef.c:130
int InputInterAllInfo(int *icnt_interall, int **iInterAllInfo, double complex *cInterAllValue, int isite1, int isigma1, int isite2, int isigma2, int isite3, int isigma3, int isite4, int isigma4, double re_value, double im_value)
Input InterAll Interactions (Operators of the same kinds are grouped together).
Definition: readdef.c:2800
int CheckKW(const char *cKW, char cKWList[][D_CharTmpReadDef], int iSizeOfKWidx, int *iKWidx)
Function of Checking keyword in NameList file.
Definition: readdef.c:154
char * cErrKW_Same
Definition: ErrorMessage.c:69
#define LOCSPIN
Definition: global.h:32
static char(* cFileNameListFile)[D_CharTmpReadDef]
Definition: readdef.c:75
int CheckQuadSite(const int iSite1, const int iSite2, const int iSite3, const int iSite4, const int iMaxNum)
Check Site Number for a quad -> (siteA, siteB, siteC, siteD).
Definition: readdef.c:1869
double eps_CG
Definition: global.h:153
char * cWarningIncorrectFormatForSpin
Definition: ErrorMessage.c:84
char * cErrLanczos_eps
Definition: ErrorMessage.c:62
int CheckTransferHermite(struct DefineList *X)
Check Hermite for Transfer integrals.
Definition: readdef.c:1898
int CheckFormatForSpinInt(const int site1, const int site2, const int site3, const int site4)
function of checking format of spin interactions
Definition: readdef.c:2384
char * cErrSetIniVec
Definition: ErrorMessage.c:52
int CheckSpinIndexForTrans(struct DefineList *X)
function of checking spin index for transfers
Definition: readdef.c:2620
char * cErrDefFileParam
Definition: ErrorMessage.c:42
int CheckWords(const char *ctmp, const char *cKeyWord)
function of checking whether ctmp is same as cKeyWord or not
Definition: readdef.c:2680
double eps_Energy
Definition: global.h:155
char * cErrNumAve
Definition: ErrorMessage.c:59
double complex *** arrayJ
Definition: struct.h:400
double eps
Definition: global.h:152
char * cErrCalcEigenVec
Definition: ErrorMessage.c:45
char * cErrNonConservedInterAll
Definition: ErrorMessage.c:77
int ReadcalcmodFile(const char *defname, struct DefineList *X)
Function of Reading calcmod file.
Definition: readdef.c:230
char * cErrOutputHamForFullDiag
Definition: ErrorMessage.c:47
char * cErrIncorrectFormatForKondoInt
Definition: ErrorMessage.c:79
char * cErrNonHermiteInterAll
Definition: ErrorMessage.c:76
char * cErrKW
Definition: ErrorMessage.c:67
long unsigned int R0
Definition: struct.h:395
char * cErrFiniteTemp
Definition: ErrorMessage.c:51
char * cErrNonHermiteTrans
Definition: ErrorMessage.c:74
char * cErrRestart
Definition: ErrorMessage.c:53
long unsigned int ishift_nspin
Definition: struct.h:398
char * cErrIncorrectSpinIndexForTrans
Definition: ErrorMessage.c:88
#define FALSE
Definition: global.h:25
int D_iKWNumDef
Definition: readdef.c:70
char * cErrInputHam
Definition: ErrorMessage.c:48
char * cErrCalcType
Definition: ErrorMessage.c:43
double complex vecB[3]
Definition: struct.h:401
int CheckTETransferHermite(struct DefineList *X, const int NTETrans, const int idx)
Check Hermite for TETransfer integrals.
Definition: readdef.c:2859
char * cWarningIncorrectFormatForSpin2
Definition: ErrorMessage.c:85
char * cErrOutputMode
Definition: ErrorMessage.c:44
char * cErrInputOutputHam
Definition: ErrorMessage.c:49
int GetDiagonalInterAll(int **InterAll, complex double *ParaInterAll, const int NInterAll, int **InterAllDiagonal, double *ParaInterAllDiagonal, int **InterAllOffDiagonal, complex double *ParaInterAllOffDiagonal, int *Chemi, int *SpinChemi, double *ParaChemi, unsigned int *NChemi, const int iCalcModel)
function of getting diagonal components
Definition: readdef.c:2162
char * cErrNonHermiteInterAllForAll
Definition: ErrorMessage.c:78
double eps_CheckImag0
Definition: global.h:156
char * cErrNonHermiteTransForAll
Definition: ErrorMessage.c:75
char * cErrExpecInterval
Definition: ErrorMessage.c:60
char * cErrOutputHam
Definition: ErrorMessage.c:46
For Boost.
Definition: struct.h:393
char * cErrLanczosTarget
Definition: ErrorMessage.c:64
char * cErrKW_InCorPair
Definition: ErrorMessage.c:70
struct EDMainCalStruct X
Definition: struct.h:431
int CheckSite(const int iSite, const int iMaxNum)
Check Site Number.
Definition: readdef.c:1821
static char cKWListOfFileNameList[][D_CharTmpReadDef]
Definition: readdef.c:46
int CheckGeneralSpinIndexForInterAll(const int isite1, const int isigma1, const int isite2, const int isigma2, const int isite3, const int isigma3, const int isite4, const int isigma4, int *iLocInfo)
function of checking spin index for all interactions
Definition: readdef.c:2593
char * cErrCUDA
Definition: ErrorMessage.c:54
int flgBoost
Flag whether use CMA algorithm.
Definition: struct.h:394
char * cErrNcond
Definition: ErrorMessage.c:57
char * cErrNLoc
Definition: ErrorMessage.c:41
void ResetInteractionNum(struct DefineList *X)
function of resetting number of interactions
Definition: readdef.c:2535
const char * cReadFileNamelist
Definition: LogMessage.c:27
int GetFileNameByKW(int iKWidx, char **FileName)
function of getting file name labeled by the keyword
Definition: readdef.c:2711
unsigned int NumarrayJ
Definition: struct.h:399
int CheckPairSite(const int iSite1, const int iSite2, const int iMaxNum)
Check Site Number for a pair -> (siteA, siteB).
Definition: readdef.c:1841
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
FILE * fopenMPI(const char *FileName, const char *mode)
MPI file I/O (open) wrapper. Only the root node (myrank = 0) should be open/read/write (small) parame...
Definition: wrapperMPI.c:105
double LargeValue
Definition: global.h:64
double eps_Lanczos
Definition: global.h:154
int CheckTotal2Sz(struct DefineList *X)
function of checking an input data of total2Sz
Definition: readdef.c:2653
Definision of system (Hamiltonian) etc.
Definition: struct.h:41
FILE * stdoutMPI
File pointer to the standard output defined in InitializeMPI()
Definition: global.h:164