HΦ  3.2.0
CalcByTEM.c File Reference

File to define functions to calculate expected values by Time evolution method. More...

#include "Common.h"
#include "readdef.h"
#include "FirstMultiply.h"
#include "Multiply.h"
#include "diagonalcalc.h"
#include "expec_energy_flct.h"
#include "expec_cisajs.h"
#include "expec_cisajscktaltdc.h"
#include "CalcByTEM.h"
#include "FileIO.h"
#include "wrapperMPI.h"
#include "HPhiTrans.h"
+ Include dependency graph for CalcByTEM.c:

Go to the source code of this file.

Functions

void MakeTEDTransfer (struct BindStruct *X, const int timeidx)
 Set transfer integrals at timeidx-th time. More...
 
void MakeTEDInterAll (struct BindStruct *X, const int timeidx)
 Set interall interactions at timeidx-th time. More...
 
int CalcByTEM (const int ExpecInterval, struct EDMainCalStruct *X)
 main function of time evolution calculation More...
 

Detailed Description

File to define functions to calculate expected values by Time evolution method.

Definition in file CalcByTEM.c.

Function Documentation

◆ CalcByTEM()

int CalcByTEM ( const int  ExpecInterval,
struct EDMainCalStruct X 
)

main function of time evolution calculation

Parameters
ExpecIntervalinterval to output expected values
Xstruct to get information of calculations.
Returns
0 normally finished
-1 unnormally finished
Author
Kota Ido (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)

Definition at line 53 of file CalcByTEM.c.

References EDMainCalStruct::Bind, c_InputEigenVectorStart, DefineList::CDataFileHead, cFileNameFlct, cFileNameNorm, cFileNameOutputEigen, cFileNameSS, cFileNameTEStep, cFileNameTimeKeep, BindStruct::Check, childfopenALL(), childfopenMPI(), cLogFlct, cLogNorm, cLogSS, cLogTEM_End, cLogTEM_Start, cLogTEStep, cTEStep, D_FileNameMax, BindStruct::Def, PhysList::doublon, PhysList::doublon2, DefineList::EDNTransfer, PhysList::energy, exitMPI(), expec_cisajs(), expec_cisajscktaltdc(), expec_energy_flct(), FALSE, GetFileNameByKW(), global_norm, CheckList::idim_max, DefineList::iInputEigenVec, DefineList::iOutputEigenVec, DefineList::iReStart, DefineList::istep, DefineList::Lanczos_max, MakeTEDInterAll(), MakeTEDTransfer(), MultiplyForTEM(), myrank, DefineList::NInterAll_OffDiagonal, DefineList::NLaser, DefineList::NTEInterAllMax, DefineList::NTETimeSteps, DefineList::NTETransferMax, PhysList::num, PhysList::num2, ParamList::OutputInterval, DefineList::Param, BindStruct::Phys, DefineList::St, stdoutMPI, step_i, step_spin, PhysList::Sz, PhysList::Sz2, DefineList::TETime, TimeKeeper(), TimeKeeperWithStep(), ParamList::TimeSlice, ParamList::Tinit, TransferWithPeierls(), TRUE, v1, and PhysList::var.

Referenced by main().

56  {
57  size_t byte_size;
58  char *defname;
59  char sdt[D_FileNameMax];
60  char sdt_phys[D_FileNameMax];
61  char sdt_norm[D_FileNameMax];
62  char sdt_flct[D_FileNameMax];
63  int rand_i=0;
64  int step_initial = 0;
65  long int i_max = 0;
66  FILE *fp;
67  double Time = X->Bind.Def.Param.Tinit;
68  double dt = ((X->Bind.Def.NLaser==0)? 0.0: X->Bind.Def.Param.TimeSlice);
69 
71  fprintf(stdoutMPI, "Error: NTETimeSteps must be larger than Lanczos_max.\n");
72  return -1;
73  }
74  step_spin = ExpecInterval;
75  X->Bind.Def.St = 0;
76  fprintf(stdoutMPI, "%s", cLogTEM_Start);
77  if (X->Bind.Def.iInputEigenVec == FALSE) {
78  fprintf(stderr, "Error: A file of Inputvector is not inputted.\n");
79  return -1;
80  } else {
81  //input v1
82  fprintf(stdoutMPI, "%s","An Initial Vector is inputted.\n");
84  GetFileNameByKW(KWSpectrumVec, &defname);
85  strcat(defname, "_rank_%d.dat");
86  sprintf(sdt, defname, myrank);
87  childfopenALL(sdt, "rb", &fp);
88  if (fp == NULL) {
89  fprintf(stderr, "Error: A file of Inputvector does not exist.\n");
90  fclose(fp);
91  exitMPI(-1);
92  }
93  byte_size = fread(&step_initial, sizeof(int), 1, fp);
94  byte_size = fread(&i_max, sizeof(long int), 1, fp);
95  if (i_max != X->Bind.Check.idim_max) {
96  fprintf(stderr, "Error: A file of Inputvector is incorrect.\n");
97  fclose(fp);
98  exitMPI(-1);
99  }
100  fread(v1, sizeof(complex double), X->Bind.Check.idim_max + 1, fp);
101  fclose(fp);
102  if (X->Bind.Def.iReStart == RESTART_NOT || X->Bind.Def.iReStart == RESTART_OUT) {
103  step_initial = 0;
104  }
105  }
106 
107  sprintf(sdt_phys, "%s", cFileNameSS);
108  if (childfopenMPI(sdt_phys, "w", &fp) != 0) {
109  return -1;
110  }
111  fprintf(fp, "%s",cLogSS);
112  fclose(fp);
113 
114  sprintf(sdt_norm, "%s", cFileNameNorm);
115  if (childfopenMPI(sdt_norm, "w", &fp) != 0) {
116  return -1;
117  }
118  fprintf(fp, "%s",cLogNorm);
119  fclose(fp);
120 
121  sprintf(sdt_flct, "%s", cFileNameFlct);
122  if (childfopenMPI(sdt_flct, "w", &fp) != 0) {
123  return -1;
124  }
125  fprintf(fp, "%s",cLogFlct);
126  fclose(fp);
127 
128 
129  int iInterAllOffDiagonal_org = X->Bind.Def.NInterAll_OffDiagonal;
130  int iTransfer_org = X->Bind.Def.EDNTransfer;
131  for (step_i = step_initial; step_i < X->Bind.Def.Lanczos_max; step_i++) {
132  X->Bind.Def.istep = step_i;
133 
134  //Reset total number of interactions (changed in MakeTED***function.)
135  X->Bind.Def.EDNTransfer = iTransfer_org;
136  X->Bind.Def.NInterAll_OffDiagonal = iInterAllOffDiagonal_org;
137 
138  if (step_i % (X->Bind.Def.Lanczos_max / 10) == 0) {
140  }
141 
142  if(X->Bind.Def.NLaser !=0) {
143  TransferWithPeierls(&(X->Bind), Time);
144  }
145  else {
146  // common procedure
147  Time = X->Bind.Def.TETime[step_i];
148  if (step_i == 0) dt = 0.0;
149  else {
150  dt = X->Bind.Def.TETime[step_i] - X->Bind.Def.TETime[step_i - 1];
151  }
152  X->Bind.Def.Param.TimeSlice = dt;
153 
154  // Set interactions
155  if(X->Bind.Def.NTETransferMax != 0 && X->Bind.Def.NTEInterAllMax!=0){
156  fprintf(stdoutMPI,
157  "Error: Time Evoluation mode does not support TEOneBody and TETwoBody interactions at the same time. \n");
158  return -1;
159  }
160  else if (X->Bind.Def.NTETransferMax > 0) { //One-Body type
161  MakeTEDTransfer(&(X->Bind), step_i);
162  }else if (X->Bind.Def.NTEInterAllMax > 0) { //Two-Body type
163  MakeTEDInterAll(&(X->Bind), step_i);
164  }
165  //[e] Yoshimi
166  }
167 
168  if(step_i == step_initial){
170  }
171  else {
173  }
174  MultiplyForTEM(&(X->Bind));
175  //Add Diagonal Parts
176  //Multiply Diagonal
177  expec_energy_flct(&(X->Bind));
178 
179  if(X->Bind.Def.NLaser >0 ) Time+=dt;
180  if (childfopenMPI(sdt_phys, "a", &fp) != 0) {
181  return -1;
182  }
183  fprintf(fp, "%.16lf %.16lf %.16lf %.16lf %.16lf %d\n", Time, X->Bind.Phys.energy, X->Bind.Phys.var,
184  X->Bind.Phys.doublon, X->Bind.Phys.num, step_i);
185  fclose(fp);
186 
187  if (childfopenMPI(sdt_norm, "a", &fp) != 0) {
188  return -1;
189  }
190  fprintf(fp, "%.16lf %.16lf %d\n", Time, global_norm, step_i);
191  fclose(fp);
192 
193  if (childfopenMPI(sdt_flct, "a", &fp) != 0) {
194  return -1;
195  }
196  fprintf(fp, "%.16lf %.16lf %.16lf %.16lf %.16lf %.16lf %.16lf %d\n", Time,X->Bind.Phys.num,X->Bind.Phys.num2, X->Bind.Phys.doublon,X->Bind.Phys.doublon2, X->Bind.Phys.Sz,X->Bind.Phys.Sz2,step_i);
197  fclose(fp);
198 
199 
200 
201  if (step_i % step_spin == 0) {
202  expec_cisajs(&(X->Bind), v1);
203  expec_cisajscktaltdc(&(X->Bind), v1);
204  }
205  if (X->Bind.Def.iOutputEigenVec == TRUE) {
206  if (step_i % X->Bind.Def.Param.OutputInterval == 0) {
208  if (childfopenALL(sdt, "wb", &fp) != 0) {
209  fclose(fp);
210  exitMPI(-1);
211  }
212  fwrite(&step_i, sizeof(step_i), 1, fp);
213  fwrite(&X->Bind.Check.idim_max, sizeof(long int), 1, fp);
214  fwrite(v1, sizeof(complex double), X->Bind.Check.idim_max + 1, fp);
215  fclose(fp);
216  }
217  }
218  }
219  if (X->Bind.Def.iOutputEigenVec == TRUE) {
220  sprintf(sdt, cFileNameOutputEigen, X->Bind.Def.CDataFileHead, rand_i, myrank);
221  if (childfopenALL(sdt, "wb", &fp) != 0) {
222  fclose(fp);
223  exitMPI(-1);
224  }
225  fwrite(&step_i, sizeof(step_i), 1, fp);
226  fwrite(&X->Bind.Check.idim_max, sizeof(long int), 1, fp);
227  fwrite(v1, sizeof(complex double), X->Bind.Check.idim_max + 1, fp);
228  fclose(fp);
229  }
230 
231  fprintf(stdoutMPI, "%s",cLogTEM_End);
232  return 0;
233 }
unsigned int NTETimeSteps
Definition: struct.h:246
int step_spin
Definition: global.h:69
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:410
int expec_cisajscktaltdc(struct BindStruct *X, double complex *vec)
Parent function to calculate two-body green&#39;s functions.
int St
0 or 1, but it affects nothing.
Definition: struct.h:80
unsigned int NTEInterAllMax
Definition: struct.h:270
unsigned long int idim_max
The dimension of the Hilbert space of this process.
Definition: struct.h:303
int iReStart
Definition: struct.h:221
const char * cFileNameTEStep
Definition: global.c:68
unsigned int NLaser
Definition: struct.h:250
int TransferWithPeierls(struct BindStruct *X, const double time)
Function of getting transfer with peierls.
Definition: HPhiTrans.c:92
double complex * v1
Definition: global.h:35
const char * cLogSS
Definition: LogMessage.c:158
#define D_FileNameMax
Definition: global.h:23
int iOutputEigenVec
ASwitch for outputting an eigenvector. 0: no output, 1:output.
Definition: struct.h:202
double * TETime
Definition: struct.h:247
const char * cFileNameFlct
Definition: global.c:71
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 * cFileNameOutputEigen
Definition: global.c:49
struct PhysList Phys
Physical quantities.
Definition: struct.h:413
const char * cLogTEStep
Definition: LogMessage.c:157
double Sz2
Expectation value of the Square of total Sz.
Definition: struct.h:361
#define TRUE
Definition: global.h:26
double global_norm
Definition: global.h:67
unsigned int NInterAll_OffDiagonal
Number of interall term (off-diagonal)
Definition: struct.h:165
void MakeTEDInterAll(struct BindStruct *X, const int timeidx)
Set interall interactions at timeidx-th time.
Definition: CalcByTEM.c:261
void MakeTEDTransfer(struct BindStruct *X, const int timeidx)
Set transfer integrals at timeidx-th time.
Definition: CalcByTEM.c:238
double num
Expectation value of the Number of electrons.
Definition: struct.h:358
const char * c_InputEigenVectorStart
Definition: LogMessage.c:60
int expec_cisajs(struct BindStruct *X, double complex *vec)
function of calculation for one body green&#39;s function
Definition: expec_cisajs.c:69
double TimeSlice
Definition: struct.h:33
const char * cTEStep
Definition: LogMessage.c:83
struct ParamList Param
Definition: struct.h:241
const char * cLogTEM_Start
double var
Expectation value of the Energy variance.
Definition: struct.h:363
int MultiplyForTEM(struct BindStruct *X)
Function of multiplying Hamiltonian for Time Evolution.
Definition: Multiply.c:80
#define FALSE
Definition: global.h:25
const char * cLogFlct
Definition: LogMessage.c:160
const char * cFileNameNorm
Definition: global.c:70
int istep
Index of TPQ step ???
Definition: struct.h:78
double Tinit
Definition: struct.h:32
const char * cFileNameTimeKeep
Definition: global.c:23
const char * cFileNameSS
Definition: global.c:69
int OutputInterval
Definition: struct.h:34
unsigned int EDNTransfer
Number of transfer integrals for calculation.
Definition: struct.h:105
double doublon
Expectation value of the Doublon.
Definition: struct.h:356
double Sz
Expectation value of the Total Sz.
Definition: struct.h:360
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:163
int expec_energy_flct(struct BindStruct *X)
Parent function to calculate expected values of energy and physical quantities.
const char * cLogNorm
Definition: LogMessage.c:159
int GetFileNameByKW(int iKWidx, char **FileName)
function of getting file name labeled by the keyword
Definition: readdef.c:2715
double energy
Expectation value of the total energy.
Definition: struct.h:355
unsigned int NTETransferMax
Definition: struct.h:255
int iInputEigenVec
Switch for reading an eigenvector. 0: no input, 1:input.
Definition: struct.h:203
struct CheckList Check
Size of the Hilbert space.
Definition: struct.h:411
double num2
Expectation value of the quare of the number of electrons.
Definition: struct.h:359
int step_i
Definition: global.h:66
char * CDataFileHead
Read from Calcmod in readdef.h. Header of output file such as Green&#39;s function.
Definition: struct.h:42
int TimeKeeperWithStep(struct BindStruct *X, const char *cFileName, const char *cTimeKeeper_Message, const char *cWriteType, const int istep)
Functions for writing a time log.
Definition: log.c:78
const char * cLogTEM_End
struct BindStruct Bind
Binded struct.
Definition: struct.h:420
double doublon2
Expectation value of the Square of doublon.
Definition: struct.h:357
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
unsigned int Lanczos_max
Maximum number of iterations.
Definition: struct.h:74
FILE * stdoutMPI
File pointer to the standard output defined in InitializeMPI()
Definition: global.h:165
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ MakeTEDInterAll()

void MakeTEDInterAll ( struct BindStruct X,
const int  timeidx 
)

Set interall interactions at timeidx-th time.

Parameters
Xstruct for getting information of interall interactions
timeidxindex of time

Definition at line 261 of file CalcByTEM.c.

References BindStruct::Def, DefineList::InterAll_OffDiagonal, DefineList::NInterAll_OffDiagonal, DefineList::NTEInterAllMax, DefineList::NTEInterAllOffDiagonal, DefineList::ParaInterAll_OffDiagonal, DefineList::ParaTEInterAllOffDiagonal, and DefineList::TEInterAllOffDiagonal.

Referenced by CalcByTEM().

261  {
262  int i, j;
263  //Clear values
264  for (i = 0; i < X->Def.NTEInterAllMax; i++) {
265  for (j = 0; j < 8; j++) {
267  }
269  }
270 
271  //Input values
272  for (i = 0; i < X->Def.NTEInterAllOffDiagonal[timeidx]; i++) {
273  for (j = 0; j < 8; j++) {
275  }
277  }
279 }
unsigned int * NTEInterAllOffDiagonal
Definition: struct.h:273
struct DefineList Def
Definision of system (Hamiltonian) etc.
Definition: struct.h:410
unsigned int NTEInterAllMax
Definition: struct.h:270
double complex * ParaInterAll_OffDiagonal
[DefineList::NInterAll_OffDiagonal] Coupling constant of off-diagonal inter-all term. malloc in setmem_def().
Definition: struct.h:170
unsigned int NInterAll_OffDiagonal
Number of interall term (off-diagonal)
Definition: struct.h:165
int ** InterAll_OffDiagonal
[DefineList::NinterAll_OffDiagonal][8] Interacted quartet
Definition: struct.h:161
double complex ** ParaTEInterAllOffDiagonal
Definition: struct.h:288
int *** TEInterAllOffDiagonal
Definition: struct.h:281
+ Here is the caller graph for this function:

◆ MakeTEDTransfer()

void MakeTEDTransfer ( struct BindStruct X,
const int  timeidx 
)

Set transfer integrals at timeidx-th time.

Parameters
Xstruct for getting information of transfer integrals
timeidxindex of time

Definition at line 238 of file CalcByTEM.c.

References BindStruct::Def, DefineList::EDGeneralTransfer, DefineList::EDNTransfer, DefineList::EDParaGeneralTransfer, DefineList::NTETransfer, DefineList::NTETransferMax, DefineList::ParaTETransfer, and DefineList::TETransfer.

Referenced by CalcByTEM().

238  {
239  int i,j;
240  //Clear values
241  for(i=0; i<X->Def.NTETransferMax ;i++) {
242  for(j =0; j<4; j++) {
243  X->Def.EDGeneralTransfer[i + X->Def.EDNTransfer][j] = 0;
244  }
246  }
247 
248  //Input values
249  for(i=0; i<X->Def.NTETransfer[timeidx] ;i++){
250  for(j =0; j<4; j++) {
251  X->Def.EDGeneralTransfer[i + X->Def.EDNTransfer][j] = X->Def.TETransfer[timeidx][i][j];
252  }
253  X->Def.EDParaGeneralTransfer[i+X->Def.EDNTransfer]=X->Def.ParaTETransfer[timeidx][i];
254  }
255  X->Def.EDNTransfer += X->Def.NTETransfer[timeidx];
256 }
struct DefineList Def
Definision of system (Hamiltonian) etc.
Definition: struct.h:410
int ** EDGeneralTransfer
Index of transfer integrals for calculation. malloc in setmem_def(). Data Format [DefineList::NTransf...
Definition: struct.h:110
double complex ** ParaTETransfer
Definition: struct.h:264
unsigned int * NTETransfer
Definition: struct.h:256
double complex * EDParaGeneralTransfer
Value of general transfer integrals by a def file. malloc in setmem_def(). Data Format [DefineList::N...
Definition: struct.h:116
int *** TETransfer
Definition: struct.h:260
unsigned int EDNTransfer
Number of transfer integrals for calculation.
Definition: struct.h:105
unsigned int NTETransferMax
Definition: struct.h:255
+ Here is the caller graph for this function: