HPhi++  3.1.0
StdFace_main.cpp
Go to the documentation of this file.
1 /*
2 HPhi-mVMC-StdFace - Common input generator
3 Copyright (C) 2015 The University of Tokyo
4 
5 This program is free software: you can redistribute it and/or modify
6 it under the terms of the GNU General Public License as published by
7 the Free Software Foundation, either version 3 of the License, or
8 (at your option) any later version.
9 
10 This program is distributed in the hope that it will be useful,
11 but WITHOUT ANY WARRANTY; without even the implied warranty of
12 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 GNU General Public License for more details.
14 
15 You should have received a copy of the GNU General Public License
16 along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
35 #include <cstdlib>
36 #include <cstdio>
37 #include <cstring>
38 #include <cctype>
39 #include <cmath>
40 #include "StdFace_vals.hpp"
41 #include "StdFace_ModelUtil.hpp"
42 #include <complex>
43 #include <iostream>
44 
45 #if defined(_HPhi)
46 
50 static void StdFace_LargeValue(struct StdIntList *StdI) {
51  int ktrans, kintr;
52  double LargeValue0;
53 
54  LargeValue0 = 0.0;
55  for (ktrans = 0; ktrans < StdI->ntrans; ktrans++) {
56  LargeValue0 += std::abs(StdI->trans[ktrans]);
57  }
58  for (kintr = 0; kintr < StdI->nintr; kintr++) {
59  LargeValue0 += std::abs(StdI->intr[kintr]);
60  }
61  for (kintr = 0; kintr < StdI->NCintra; kintr++) {
62  LargeValue0 += std::abs(StdI->Cintra[kintr]);
63  }
64  for (kintr = 0; kintr < StdI->NCinter; kintr++) {
65  LargeValue0 += std::abs(StdI->Cinter[kintr]);
66  }
67  for (kintr = 0; kintr < StdI->NEx; kintr++) {
68  LargeValue0 += 2.0 * std::abs(StdI->Ex[kintr]);
69  }
70  for (kintr = 0; kintr < StdI->NPairLift; kintr++) {
71  LargeValue0 += 2.0 * std::abs(StdI->PairLift[kintr]);
72  }
73  for (kintr = 0; kintr < StdI->NHund; kintr++) {
74  LargeValue0 += 2.0 * std::abs(StdI->Hund[kintr]);
75  }
76  LargeValue0 /= (double)StdI->nsite;
77  StdFace_PrintVal_d("LargeValue", &StdI->LargeValue, LargeValue0);
78 }/*static void StdFace_LargeValue*/
83 static void PrintCalcMod(struct StdIntList *StdI)
84 {
85  FILE *fp;
86  int iCalcType, iCalcModel, iRestart, iCalcSpec,
87  iCalcEigenvec, iInitialVecTpye, InputEigenVec, OutputEigenVec;
88  /*
89  First, check all parameters and exit if invalid parameters
90  */
91  fprintf(stdout, "\n @ CalcMod\n\n");
92  /*
93  Method
94  */
95  iCalcEigenvec = 0;
96  if (strcmp(StdI->method, "****") == 0){
97  fprintf(stdout, "ERROR ! Method is NOT specified !\n");
98  StdFace_exit(-1);
99  }
100  else if (strcmp(StdI->method, "lanczos") == 0) iCalcType = 0;
101  else if (strcmp(StdI->method, "lanczosenergy") == 0) {
102  iCalcType = 0;
103  iCalcEigenvec = 1;
104  }
105  else if (strcmp(StdI->method, "tpq") == 0) iCalcType = 1;
106  else if (strcmp(StdI->method, "fulldiag") == 0 ) iCalcType = 2;
107  else if (strcmp(StdI->method, "cg") == 0) iCalcType = 3;
108  else if (strcmp(StdI->method, "timeevolution") == 0) iCalcType = 4;
109  else{
110  fprintf(stdout, "\n ERROR ! Unsupported Solver : %s\n", StdI->method);
111  StdFace_exit(-1);
112  }/*if (strcmp(StdI->method, METHODS) != 0*/
113  if (iCalcType != 4) StdI->PumpBody = 0;
114  /*
115  Model
116  */
117  if (strcmp(StdI->model, "hubbard") == 0) {
118  if (StdI->lGC == 0)iCalcModel = 0;
119  else iCalcModel = 3;
120  }/*if (strcmp(StdI->model, "hubbard") == 0)*/
121  else if (strcmp(StdI->model, "spin") == 0) {
122  if (StdI->lGC == 0)iCalcModel = 1;
123  else iCalcModel = 4;
124  }/*if (strcmp(StdI->model, "spin") == 0)*/
125  else if (strcmp(StdI->model, "kondo") == 0) {
126  if (StdI->lGC == 0)iCalcModel = 2;
127  else iCalcModel = 5;
128  }/*if (strcmp(StdI->model, "kondo") == 0)*/
129  /*
130  Restart
131  */
132  if (strcmp(StdI->Restart, "****") == 0) {
133  strcpy(StdI->Restart, "none\0");
134  fprintf(stdout, " Restart = none ###### DEFAULT VALUE IS USED ######\n");
135  iRestart = 0;
136  }/*if (strcmp(StdI->Restart, "****") == 0)*/
137  else {
138  fprintf(stdout, " Restart = %s\n", StdI->Restart);
139  if (strcmp(StdI->Restart, "none") == 0) iRestart = 0;
140  else if (strcmp(StdI->Restart, "restart_out") == 0 ||
141  strcmp(StdI->Restart, "save") == 0) iRestart = 1;
142  else if (strcmp(StdI->Restart, "restartsave") == 0 ||
143  strcmp(StdI->Restart, "restart") == 0) iRestart = 2;
144  else if (strcmp(StdI->Restart, "restart_in") == 0) iRestart = 3;
145  else {
146  fprintf(stdout, "\n ERROR ! Restart Mode : %s\n", StdI->Restart);
147  StdFace_exit(-1);
148  }
149  }/*if (strcmp(StdI->Restart, "****") != 0)*/
150  /*
151  InitialVecType
152  */
153  if (strcmp(StdI->InitialVecType, "****") == 0) {
154  strcpy(StdI->InitialVecType, "c\0");
155  fprintf(stdout, " InitialVecType = c ###### DEFAULT VALUE IS USED ######\n");
156  iInitialVecTpye = 0;
157  }/*if (strcmp(StdI->InitialVecType, "****") == 0)*/
158  else {
159  fprintf(stdout, " InitialVecType = %s\n", StdI->InitialVecType);
160  if (strcmp(StdI->InitialVecType, "c") == 0) iInitialVecTpye = 0;
161  else if (strcmp(StdI->InitialVecType, "r") == 0) iInitialVecTpye = 1;
162  else {
163  fprintf(stdout, "\n ERROR ! Restart Mode : %s\n", StdI->Restart);
164  StdFace_exit(-1);
165  }
166  }/*if (strcmp(StdI->InitialVecType, "****") != 0)*/
167  /*
168  EigenVecIO
169  */
170  InputEigenVec = 0;
171  OutputEigenVec = 0;
172  if (strcmp(StdI->EigenVecIO, "****") == 0) {
173  strcpy(StdI->EigenVecIO, "none\0");
174  fprintf(stdout, " EigenVecIO = none ###### DEFAULT VALUE IS USED ######\n");
175  }/*if (strcmp(StdI->EigenVecIO, "****") == 0)*/
176  else {
177  fprintf(stdout, " EigenVecIO = %s\n", StdI->EigenVecIO);
178  if (strcmp(StdI->EigenVecIO, "none") == 0) InputEigenVec = 0;
179  else if (strcmp(StdI->EigenVecIO, "in") == 0) InputEigenVec = 1;
180  else if (strcmp(StdI->EigenVecIO, "out") == 0) OutputEigenVec = 1;
181  else if (strcmp(StdI->EigenVecIO, "inout") == 0) {
182  InputEigenVec = 1;
183  OutputEigenVec = 1;
184  }/*if (strcmp(StdI->EigenVecIO, "inout") == 0)*/
185  else {
186  fprintf(stdout, "\n ERROR ! EigenVecIO Mode : %s\n", StdI->Restart);
187  StdFace_exit(-1);
188  }
189  }/*if (strcmp(StdI->EigenVecIO, "****") != 0)*/
190  if (strcmp(StdI->method, "timeevolution") == 0) InputEigenVec = 1;
191  /*
192  CalcSpec
193  */
194  if (strcmp(StdI->CalcSpec, "****") == 0) {
195  strcpy(StdI->CalcSpec, "none\0");
196  fprintf(stdout, " CalcSpec = none ###### DEFAULT VALUE IS USED ######\n");
197  iCalcSpec = 0;
198  }/*if (strcmp(StdI->CalcSpec, "****") == 0)*/
199  else {
200  fprintf(stdout, " CalcSpec = %s\n", StdI->CalcSpec);
201  if (strcmp(StdI->CalcSpec, "none") == 0) iCalcSpec = 0;
202  else if (strcmp(StdI->CalcSpec, "normal") == 0) iCalcSpec = 1;
203  else if (strcmp(StdI->CalcSpec, "noiteration") == 0) iCalcSpec = 2;
204  else if (strcmp(StdI->CalcSpec, "restart_out") == 0) iCalcSpec = 3;
205  else if (strcmp(StdI->CalcSpec, "restart_in") == 0) iCalcSpec = 4;
206  else if (strcmp(StdI->CalcSpec, "restartsave") == 0 ||
207  strcmp(StdI->CalcSpec, "restart") == 0) iCalcSpec = 5;
208  else if (strcmp(StdI->CalcSpec, "scratch") == 0) iCalcSpec = 6;
209  else {
210  fprintf(stdout, "\n ERROR ! CalcSpec : %s\n", StdI->CalcSpec);
211  StdFace_exit(-1);
212  }
213  }/*if (strcmp(StdI->CalcSpec, "****") != 0)*/
214 
215  fp = fopen("calcmod.def", "w");
216  fprintf(fp, "#CalcType = 0:Lanczos, 1:TPQCalc, 2:FullDiag, 3:CG, 4:Time-evolution\n");
217  fprintf(fp, "#CalcModel = 0:Hubbard, 1:Spin, 2:Kondo, 3:HubbardGC, 4:SpinGC, 5:KondoGC\n");
218  fprintf(fp, "#Restart = 0:None, 1:Save, 2:Restart&Save, 3:Restart\n");
219  fprintf(fp, "#CalcSpec = 0:None, 1:Normal, 2:No H*Phi, 3:Save, 4:Restart, 5:Restart&Save, 6:Scratch\n");
220  fprintf(fp, "CalcType %3d\n", iCalcType);
221  fprintf(fp, "CalcModel %3d\n", iCalcModel);
222  fprintf(fp, "ReStart %3d\n", iRestart);
223  fprintf(fp, "CalcSpec %3d\n", iCalcSpec);
224  fprintf(fp, "CalcEigenVec %3d\n", iCalcEigenvec);
225  fprintf(fp, "InitialVecType %3d\n", iInitialVecTpye);
226  fprintf(fp, "InputEigenVec %3d\n", InputEigenVec);
227  fprintf(fp, "OutputEigenVec %3d\n", OutputEigenVec);
228  fflush(fp);
229  fclose(fp);
230  fprintf(stdout, " calcmod.def is written.\n\n");
231 }/*static void PrintCalcMod*/
236 static void PrintExcitation(struct StdIntList *StdI) {
237  FILE *fp;
238  int NumOp, **spin, isite, ispin, icell, itau, iEx, lR;
239  double *coef, Cphase, S, Sz;
240  double *fourier_r, *fourier_i;
241 
242  if (strcmp(StdI->model, "spin") == 0 && StdI->S2 > 1) {
243  coef = (double *)malloc(sizeof(double) * (StdI->S2 + 1));
244  spin = (int **)malloc(sizeof(int*) * (StdI->S2 + 1));
245  for (ispin = 0; ispin < StdI->S2 + 1; ispin++) spin[ispin] = (int *)malloc(sizeof(int) * 2);
246  }
247  else {
248  coef = (double *)malloc(sizeof(double) * 2);
249  spin = (int **)malloc(sizeof(int*) * 2);
250  for (ispin = 0; ispin < 2; ispin++) spin[ispin] = (int *)malloc(sizeof(int) * 2);
251  }
252 
253  fourier_r = (double *)malloc(sizeof(double) * StdI->nsite);
254  fourier_i = (double *)malloc(sizeof(double) * StdI->nsite);
255 
256  fprintf(stdout, "\n @ Spectrum\n\n");
257 
258  StdFace_PrintVal_d("SpectrumQW", &StdI->SpectrumQ[0], 0.0);
259  StdFace_PrintVal_d("SpectrumQL", &StdI->SpectrumQ[1], 0.0);
260  StdFace_PrintVal_d("SpectrumQH", &StdI->SpectrumQ[2], 0.0);
261 
262  if (strcmp(StdI->SpectrumType, "****") == 0) {
263  strcpy(StdI->SpectrumType, "szsz\0");
264  fprintf(stdout, " SpectrumType = szsz ###### DEFAULT VALUE IS USED ######\n");
265  if (strcmp(StdI->model, "spin") == 0) {
266  NumOp = StdI->S2 + 1;
267  for (ispin = 0; ispin <= StdI->S2; ispin++) {
268  Sz = (double)ispin - (double)StdI->S2 * 0.5;
269  coef[ispin] = Sz;
270  spin[ispin][0] = ispin;
271  spin[ispin][1] = ispin;
272  }
273  }
274  else {
275  NumOp = 2;
276  coef[0] = 0.5;
277  coef[1] = -0.5;
278  spin[0][0] = 0;
279  spin[0][1] = 0;
280  spin[1][0] = 1;
281  spin[1][1] = 1;
282  }
283  StdI->SpectrumBody = 2;
284  lR = 0;
285  }
286  else {
287  fprintf(stdout, " SpectrumType = %s\n", StdI->SpectrumType);
288  if (strcmp(StdI->SpectrumType, "szsz") == 0 ||
289  strcmp(StdI->SpectrumType, "szsz_r") == 0) {
290  if (strcmp(StdI->model, "spin") == 0) {
291  NumOp = StdI->S2 + 1;
292  for (ispin = 0; ispin <= StdI->S2; ispin++) {
293  Sz = (double)ispin - (double)StdI->S2 * 0.5;
294  coef[ispin] = Sz;
295  spin[ispin][0] = ispin;
296  spin[ispin][1] = ispin;
297  }
298  }
299  else {
300  NumOp = 2;
301  coef[0] = 0.5;
302  coef[1] = -0.5;
303  spin[0][0] = 0;
304  spin[0][1] = 0;
305  spin[1][0] = 1;
306  spin[1][1] = 1;
307  }
308  if (strcmp(StdI->SpectrumType, "szsz") == 0) lR = 0;
309  else lR = 1;
310  StdI->SpectrumBody = 2;
311  }
312  else if (strcmp(StdI->SpectrumType, "s+s-") == 0 ||
313  strcmp(StdI->SpectrumType, "s+s-_r") == 0) {
314  if (strcmp(StdI->model, "spin") == 0 && StdI->S2 > 1) {
315  NumOp = StdI->S2;
316  S = (double)StdI->S2 * 0.5;
317  for (ispin = 0; ispin < StdI->S2; ispin++) {
318  Sz = (double)ispin - (double)StdI->S2 * 0.5;
319  coef[ispin] = sqrt(S*(S + 1.0) - Sz*(Sz + 1.0));
320  spin[ispin][0] = ispin;
321  spin[ispin][1] = ispin + 1;
322  }
323  }
324  else {
325  NumOp = 1;
326  coef[0] = 1.0;
327  spin[0][0] = 0;
328  spin[0][1] = 1;
329  }
330  if (strcmp(StdI->SpectrumType, "s+s-") == 0) lR = 0;
331  else lR = 1;
332  StdI->SpectrumBody = 2;
333  }
334  else if (strcmp(StdI->SpectrumType, "density") == 0 ||
335  strcmp(StdI->SpectrumType, "density_r") == 0) {
336  NumOp = 2;
337  coef[0] = 1.0;
338  coef[1] = 1.0;
339  spin[0][0] = 0;
340  spin[0][1] = 0;
341  spin[1][0] = 1;
342  spin[1][1] = 1;
343  if (strcmp(StdI->SpectrumType, "density") == 0) lR = 0;
344  else lR = 1;
345  StdI->SpectrumBody = 2;
346  }
347  else if (strcmp(StdI->SpectrumType, "up") == 0 ||
348  strcmp(StdI->SpectrumType, "up_r") == 0) {
349  NumOp = 1;
350  coef[0] = 1.0;
351  spin[0][0] = 0;
352  if (strcmp(StdI->SpectrumType, "up") == 0) lR = 0;
353  else lR = 1;
354  StdI->SpectrumBody = 1;
355  }
356  else if (strcmp(StdI->SpectrumType, "down") == 0 ||
357  strcmp(StdI->SpectrumType, "down_r") == 0) {
358  NumOp = 1;
359  coef[0] = 1.0;
360  spin[0][0] = 1;
361  if (strcmp(StdI->SpectrumType, "down") == 0) lR = 0;
362  else lR = 1;
363  StdI->SpectrumBody = 1;
364  }
365  else {
366  fprintf(stdout, "\n ERROR ! SpectrumType : %s\n", StdI->SpectrumType);
367  StdFace_exit(-1);
368  }
369  }
370 
371  isite = 0;
372  for (icell = 0; icell < StdI->NCell; icell++) {
373  for (itau = 0; itau < StdI->NsiteUC; itau++) {
374  Cphase = (StdI->Cell[icell][0] + StdI->tau[itau][0])*StdI->SpectrumQ[0]
375  + (StdI->Cell[icell][1] + StdI->tau[itau][1])*StdI->SpectrumQ[1]
376  + (StdI->Cell[icell][2] + StdI->tau[itau][2])*StdI->SpectrumQ[2];
377  fourier_r[isite] = cos(2.0*StdI->pi*Cphase);
378  fourier_i[isite] = sin(2.0*StdI->pi*Cphase);
379  isite += 1;
380  }
381  }
382  if (strcmp(StdI->model, "kondo") == 0) {
383  for (isite = 0; isite < StdI->nsite / 2; isite++) {
384  fourier_r[isite + StdI->nsite / 2] = fourier_r[isite];
385  fourier_i[isite + StdI->nsite / 2] = fourier_i[isite];
386  }/*for (isite = 0; isite < StdI->nsite; isite++)*/
387  }/*if (strcmp(StdI->model, "kondo") == 0)*/
388 
389  if (StdI->SpectrumBody == 1) {
390  fp = fopen("single.def", "w");
391  fprintf(fp, "=============================================\n");
392  if (lR == 0) fprintf(fp, "NSingle %d\n", 2);
393  else fprintf(fp, "NSingle %d\n", 1+ StdI->nsite);
394  fprintf(fp, "=============================================\n");
395  fprintf(fp, "============== Single Excitation ============\n");
396  fprintf(fp, "=============================================\n");
397  if (lR == 0) {
398  if (strcmp(StdI->model, "kondo") == 0) {
399  for (iEx = 0; iEx < 2; iEx++) {
400  fprintf(fp, "%d\n", StdI->nsite / 2 * NumOp);
401  for (isite = StdI->nsite / 2; isite < StdI->nsite; isite++) {
402  fprintf(fp, "%d %d 0 %25.15f %25.15f\n", isite, spin[0][0],
403  fourier_r[isite] * coef[0], fourier_i[isite] * coef[0]);
404  }/*for (isite = 0; isite < StdI->nsite; isite++)*/
405  }/*for (iEx = 0; iEx < 2; iEx++)*/
406  }/*if (strcmp(StdI->model, "kondo") == 0)*/
407  else {
408  for (iEx = 0; iEx < 2; iEx++) {
409  fprintf(fp, "%d\n", StdI->nsite * NumOp);
410  for (isite = 0; isite < StdI->nsite; isite++) {
411  fprintf(fp, "%d %d 0 %25.15f %25.15f\n", isite, spin[0][0],
412  fourier_r[isite] * coef[0], fourier_i[isite] * coef[0]);
413  }/*for (isite = 0; isite < StdI->nsite; isite++)*/
414  }/*for (iEx = 0; iEx < 2; iEx++)*/
415  }
416  }/*if (lR == 0)*/
417  else {
418  if (strcmp(StdI->model, "kondo") == 0) {
419  fprintf(fp, "%d\n", NumOp);
420  fprintf(fp, "%d %d 0 %25.15f 0.0\n", StdI->nsite / 2, spin[0][0], coef[0]);
421  for (isite = StdI->nsite / 2; isite < StdI->nsite; isite++) {
422  fprintf(fp, "%d\n", NumOp);
423  fprintf(fp, "%d %d 0 %25.15f 0.0\n", isite, spin[0][0], coef[0]);
424  }
425  }
426  else {
427  fprintf(fp, "%d\n", NumOp);
428  fprintf(fp, "%d %d 0 %25.15f 0.0\n", 0, spin[0][0], coef[0]);
429  for (isite = 0; isite < StdI->nsite; isite++) {
430  fprintf(fp, "%d\n", NumOp);
431  fprintf(fp, "%d %d 0 %25.15f 0.0\n", isite, spin[0][0], coef[0]);
432  }
433  }
434  }/*if (lR != 0)*/
435  fprintf(stdout, " single.def is written.\n\n");
436  }/*if (StdI->SpectrumBody == 1)*/
437  else {
438  fp = fopen("pair.def", "w");
439  fprintf(fp, "=============================================\n");
440  if (lR == 0) fprintf(fp, "NPair %d\n", 2);
441  else fprintf(fp, "NSingle %d\n", 1 + StdI->nsite);
442  fprintf(fp, "=============================================\n");
443  fprintf(fp, "=============== Pair Excitation =============\n");
444  fprintf(fp, "=============================================\n");
445  if (lR == 0) {
446  for (iEx = 0; iEx < 2; iEx++) {
447  fprintf(fp, "%d\n", StdI->nsite * NumOp);
448  for (isite = 0; isite < StdI->nsite; isite++) {
449  for (ispin = 0; ispin < NumOp; ispin++) {
450  fprintf(fp, "%d %d %d %d 1 %25.15f %25.15f\n",
451  isite, spin[ispin][0], isite, spin[ispin][1],
452  fourier_r[isite] * coef[ispin], fourier_i[isite] * coef[ispin]);
453  }
454  }
455  }/*for (iEx = 0; iEx < 2; iEx++)*/
456  }/*if (lR == 0)*/
457  else {
458  fprintf(fp, "%d\n", NumOp);
459  for (ispin = 0; ispin < NumOp; ispin++) {
460  fprintf(fp, "%d %d %d %d 1 %25.15f 0.0\n",
461  0, spin[ispin][0], 0, spin[ispin][1], coef[ispin]);
462  }
463  for (isite = 0; isite < StdI->nsite; isite++) {
464  fprintf(fp, "%d\n", NumOp);
465  for (ispin = 0; ispin < NumOp; ispin++) {
466  fprintf(fp, "%d %d %d %d 1 %25.15f 0.0\n",
467  isite, spin[ispin][0], isite, spin[ispin][1], coef[ispin]);
468  }
469  }
470  }/*if (lR != 0)*/
471  fprintf(stdout, " pair.def is written.\n\n");
472  }/*if (StdI->SpectrumBody == 2)*/
473  fflush(fp);
474  fclose(fp);
475 
476  free(fourier_r);
477  free(fourier_i);
478  if (strcmp(StdI->model, "spin") == 0)
479  for (ispin = 0; ispin < StdI->S2 + 1; ispin++) free(spin[ispin]);
480  else
481  for (ispin = 0; ispin < 2; ispin++) free(spin[ispin]);
482  free(spin);
483  free(coef);
484 
485 }/*static void PrintExcitation()*/
486 /*
487 @brief Compute vectorpotential
488 */
489 static void VectorPotential(struct StdIntList *StdI) {
490  FILE *fp;
491  int it, ii;
492  double time;
493  double **Et;
494 
495  fprintf(stdout, "\n @ Time-evolution\n\n");
496 
497  StdFace_PrintVal_d("VecPotW", &StdI->VecPot[0], 0.0);
498  StdFace_PrintVal_d("VecPotL", &StdI->VecPot[1], 0.0);
499  StdFace_PrintVal_d("VecPotH", &StdI->VecPot[2], 0.0);
500  StdFace_PrintVal_i("Lanczos_max", &StdI->Lanczos_max, 1000);
501  StdFace_PrintVal_d("dt", &StdI->dt, 0.1);
502  StdFace_PrintVal_d("freq", &StdI->freq, 0.1);
503  StdFace_PrintVal_d("tshift", &StdI->tshift, 0.0);
504  StdFace_PrintVal_d("tdump", &StdI->tdump, 0.1);
505  StdFace_PrintVal_d("Uquench", &StdI->Uquench, 0.0);
506  StdFace_PrintVal_i("ExpandCoef", &StdI->ExpandCoef, 10);
507  StdI->At = (double **)malloc(sizeof(double*) * StdI->Lanczos_max);
508  Et = (double **)malloc(sizeof(double*) * StdI->Lanczos_max);
509  for (it = 0; it < StdI->Lanczos_max; it++) {
510  StdI->At[it] = (double *)malloc(sizeof(double) * 3);
511  Et[it] = (double *)malloc(sizeof(double) * 3);
512  }
513 
514  if (strcmp(StdI->PumpType, "****") == 0) {
515  strcpy(StdI->PumpType, "quench\0");
516  fprintf(stdout, " PumpType = quench ###### DEFAULT VALUE IS USED ######\n");
517  StdI->PumpBody = 2;
518  }/*if (strcmp(StdI->PumpType, "****")*/
519  else {
520  fprintf(stdout, " PumpType = %s\n", StdI->PumpType);
521  if (strcmp(StdI->PumpType, "quench") == 0) {
522  StdI->PumpBody = 2;
523  }/*if (strcmp(StdI->PumpType, "quench")*/
524  else if (strcmp(StdI->PumpType, "pulselaser") == 0) {
525  for (it = 0; it < StdI->Lanczos_max; it++) {
526  time = StdI->dt*(double)it;
527  for (ii = 0; ii < 3; ii++) {
528  StdI->At[it][ii] = StdI->VecPot[ii] * cos(StdI->freq*(time - StdI->tshift))
529  * exp(-0.5* (time - StdI->tshift)*(time - StdI->tshift) / (StdI->tdump*StdI->tdump));
530  Et[it][ii] = -StdI->VecPot[ii]
531  * (
532  (StdI->tshift - time) / (StdI->tdump*StdI->tdump) * cos(StdI->freq*(time - StdI->tshift))
533  - StdI->freq* sin(StdI->freq*(time - StdI->tshift))
534  )
535  * exp(-0.5* (time - StdI->tshift)*(time - StdI->tshift) / (StdI->tdump*StdI->tdump));
536  }
537  }/*for (it = 0; it < StdI->Lanczos_max; it++)*/
538  StdI->PumpBody = 1;
539  }/*if (strcmp(StdI->PumpType, "pulselaser") == 0)*/
540  else if (strcmp(StdI->PumpType, "aclaser") == 0) {
541  for (it = 0; it < StdI->Lanczos_max; it++) {
542  time = StdI->dt*(double)it;
543  for (ii = 0; ii < 3; ii++) {
544  StdI->At[it][ii] = StdI->VecPot[ii] * sin(StdI->freq*(time - StdI->tshift));
545  Et[it][ii] = StdI->VecPot[ii] * cos(StdI->freq*(time - StdI->tshift)) * StdI->freq;
546  }
547  }/*for (it = 0; it < StdI->Lanczos_max; it++)*/
548  StdI->PumpBody = 1;
549  }/*if (strcmp(StdI->PumpType, "aclaser") == 0)*/
550  else if (strcmp(StdI->PumpType, "dclaser") == 0) {
551  for (it = 0; it < StdI->Lanczos_max; it++) {
552  time = StdI->dt*(double)it;
553  for (ii = 0; ii < 3; ii++) {
554  StdI->At[it][ii] = StdI->VecPot[ii] * time;
555  Et[it][ii] = -StdI->VecPot[ii];
556  }
557  }/*for (it = 0; it < StdI->Lanczos_max; it++)*/
558  StdI->PumpBody = 1;
559  }/* if (strcmp(StdI->PumpType, "dclaser") == 0)*/
560  else {
561  fprintf(stdout, "\n ERROR ! PumpType : %s\n", StdI->PumpType);
562  StdFace_exit(-1);
563  }
564  }/*if (! strcmp(StdI->PumpType, "****"))*/
565 
566  if (StdI->PumpBody == 1) {
567  fp = fopen("potential.dat", "w");
568  fprintf(fp, "# Time A_W A_L A_H E_W E_L E_H\n");
569  for (it = 0; it < StdI->Lanczos_max; it++) {
570  time = StdI->dt*(double)it;
571  fprintf(fp, "%f %f %f %f %f %f %f\n",
572  time, StdI->At[it][0], StdI->At[it][1], StdI->At[it][2], Et[it][0], Et[it][1], Et[it][2]);
573  }
574  fflush(fp);
575  fclose(fp);
576  }/*if (StdI->PumpBody == 1)*/
577 
578  for (it = 0; it < StdI->Lanczos_max; it++) free(Et[it]);
579  free(Et);
580 }/*static void VectorPotential(struct StdIntList *StdI)*/
585 static void PrintPump(struct StdIntList *StdI) {
586  FILE *fp;
587  int it, isite, ipump, jpump, npump0;
588 
589  if (StdI->PumpBody == 1) {
590 
591  fp = fopen("teone.def", "w");
592  fprintf(fp, "=============================================\n");
593  fprintf(fp, "AllTimeStep %d\n", StdI->Lanczos_max);
594  fprintf(fp, "=============================================\n");
595  fprintf(fp, "========= OneBody Time Evolution ==========\n");
596  fprintf(fp, "=============================================\n");
597  for (it = 0; it < StdI->Lanczos_max; it++) {
598  /*
599  Sum equivalent pumping
600  */
601  for (ipump = 0; ipump < StdI->npump[it]; ipump++) {
602  for (jpump = ipump + 1; jpump < StdI->npump[it]; jpump++) {
603  if (StdI->pumpindx[it][ipump][0] == StdI->pumpindx[it][jpump][0]
604  && StdI->pumpindx[it][ipump][1] == StdI->pumpindx[it][jpump][1]
605  && StdI->pumpindx[it][ipump][2] == StdI->pumpindx[it][jpump][2]
606  && StdI->pumpindx[it][ipump][3] == StdI->pumpindx[it][jpump][3]) {
607  StdI->pump[it][ipump] = StdI->pump[it][ipump] + StdI->pump[it][jpump];
608  StdI->pump[it][jpump] = 0.0;
609  }
610  }/*for (ktrans = jtrans + 1; ktrans < StdI->ntrans; ktrans++)*/
611  }/*for (jtrans = 0; jtrans < StdI->ntrans; jtrans++)*/
612  /*
613  Count the number of finite pumping
614  */
615  npump0 = 0;
616  for (ipump = 0; ipump < StdI->npump[it]; ipump++)
617  if (std::abs(StdI->pump[it][ipump]) > 0.000001) npump0 += 1;
618 
619  fprintf(fp, "%f %d\n", StdI->dt*(double)it, npump0);
620  for (ipump = 0; ipump < StdI->npump[it]; ipump++) {
621 
622  if (std::abs(StdI->pump[it][ipump]) <= 0.000001) continue;
623 
624  fprintf(fp, "%5d %5d %5d %5d %25.15f %25.15f\n",
625  StdI->pumpindx[it][ipump][0], StdI->pumpindx[it][ipump][1],
626  StdI->pumpindx[it][ipump][2], StdI->pumpindx[it][ipump][3],
627  real(StdI->pump[it][ipump]), imag(StdI->pump[it][ipump]));
628  }/*for (itrans = 0; itrans < StdI->ntrans; itrans++)*/
629  }/*for (it = 0; it < StdI->Lanczos_max; it++)*/
630  fprintf(stdout, " teone.def is written.\n\n");
631  }
632  else {
633  fp = fopen("tetwo.def", "w");
634  fprintf(fp, "=============================================\n");
635  fprintf(fp, "AllTimeStep %d\n", StdI->Lanczos_max);
636  fprintf(fp, "=============================================\n");
637  fprintf(fp, "========== TwoBody Time Evolution ===========\n");
638  fprintf(fp, "=============================================\n");
639  for (it = 0; it < StdI->Lanczos_max; it++) {
640  fprintf(fp, "%f %d\n", StdI->dt*(double)it, StdI->nsite);
641  for (isite = 0; isite < StdI->nsite; isite++) {
642  fprintf(fp, "%5d %5d %5d %5d %5d %5d %5d %5d %25.15f %25.15f\n",
643  isite, 0, isite, 0, isite, 1, isite, 1, StdI->Uquench, 0.0);
644  }/*for (isite = 0; isite < StdI->nsite; isite++)*/
645  }/*for (it = 0; it < StdI->Lanczos_max; it++)*/
646  fprintf(stdout, " tetwo.def is written.\n\n");
647  }
648  fflush(fp);
649  fclose(fp);
650 }/*tatic void PrintPump*/
651 #elif defined(_mVMC)
652 
656 static void PrintOrb(struct StdIntList *StdI) {
657  FILE *fp;
658  int isite, jsite, iOrb;
659 
660  fp = fopen("orbitalidx.def", "w");
661  fprintf(fp, "=============================================\n");
662  fprintf(fp, "NOrbitalIdx %10d\n", StdI->NOrb);
663  fprintf(fp, "ComplexType %10d\n", StdI->ComplexType);
664  fprintf(fp, "=============================================\n");
665  fprintf(fp, "=============================================\n");
666 
667  for (isite = 0; isite < StdI->nsite; isite++) {
668  for (jsite = 0; jsite < StdI->nsite; jsite++) {
669  if (StdI->AntiPeriod[0] == 1 || StdI->AntiPeriod[1] == 1 || StdI->AntiPeriod[2] == 1) {
670  fprintf(fp, "%5d %5d %5d %5d\n", isite, jsite, StdI->Orb[isite][jsite], StdI->AntiOrb[isite][jsite]);
671  }
672  else {
673  fprintf(fp, "%5d %5d %5d\n", isite, jsite, StdI->Orb[isite][jsite]);
674  }
675  }/*for (jsite = 0; jsite < isite; jsite++)*/
676  }/*for (isite = 0; isite < StdI->nsite; isite++)*/
677 
678  for (iOrb = 0; iOrb < StdI->NOrb; iOrb++)
679  fprintf(fp, "%5d %5d\n", iOrb, 1);
680 
681  fflush(fp);
682  fclose(fp);
683  fprintf(stdout, " orbitalidx.def is written.\n");
684 
685  for (isite = 0; isite < StdI->nsite; isite++) free(StdI->Orb[isite]);
686  free(StdI->Orb);
687 }/*void PrintOrb*/
692 static void PrintOrbPara(struct StdIntList *StdI) {
693  FILE *fp;
694  int isite, jsite, NOrbGC, iOrbGC, isite1, jsite1, iorb;
695  int **OrbGC, **AntiOrbGC;
699  OrbGC = (int **)malloc(sizeof(int*) * StdI->nsite);
700  AntiOrbGC = (int **)malloc(sizeof(int*) * StdI->nsite);
701  for (isite = 0; isite < StdI->nsite; isite++) {
702  OrbGC[isite] = (int *)malloc(sizeof(int) * StdI->nsite);
703  AntiOrbGC[isite] = (int *)malloc(sizeof(int) * StdI->nsite);
704  for (jsite = 0; jsite < StdI->nsite; jsite++) {
705  OrbGC[isite][jsite] = StdI->Orb[isite][jsite];
706  AntiOrbGC[isite][jsite] = StdI->AntiOrb[isite][jsite];
707  }/*for (jsite = 0; jsite < isite; jsite++)*/
708  }/*for (isite = 0; isite < StdI->nsite; isite++)*/
712  for (iorb = 0; iorb < StdI->NOrb; iorb++) {
713  for (isite = 0; isite < StdI->nsite; isite++) {
714  for (jsite = 0; jsite < StdI->nsite; jsite++) {
715  if (OrbGC[isite][jsite] == iorb) {
716  OrbGC[jsite][isite] = OrbGC[isite][jsite];
717  }
718  }/*for (jsite = 0; jsite < isite; jsite++)*/
719  }/*for (isite = 0; isite < StdI->nsite; isite++)*/
720  }/*for (iorb = 0; iorb < StdI->NOrb; iorb++)*/
721 
722  NOrbGC = 0;
723  for (isite = 0; isite < StdI->nsite; isite++) {
724  for (jsite = 0; jsite < isite; jsite++) {
725  if (OrbGC[isite][jsite] >= 0) {
726  iOrbGC = OrbGC[isite][jsite];
727  NOrbGC -= 1;
728  for (isite1 = 0; isite1 < StdI->nsite; isite1++) {
729  for (jsite1 = 0; jsite1 < StdI->nsite; jsite1++) {
730  if (OrbGC[isite1][jsite1] == iOrbGC)
731  OrbGC[isite1][jsite1] = NOrbGC;
732  }/*for (jsite1 = 0; jsite1 < StdI->nsite; jsite1++)*/
733  }/*for (isite1 = 0; isite1 < StdI->nsite; isite1++)*/
734  }/*if (OrbGC[isite][jsite] >= 0)*/
735  }/*for (jsite = 0; jsite < isite; jsite++)*/
736  }/*for (isite = 0; isite < StdI->nsite; isite++)*/
737 
738  NOrbGC = -NOrbGC;
739  for (isite = 0; isite < StdI->nsite; isite++) {
740  for (jsite = 0; jsite < StdI->nsite; jsite++) {
741  OrbGC[isite][jsite] = -1 - OrbGC[isite][jsite];
742  }/*for (jsite = 0; jsite < isite; jsite++)*/
743  }/*for (isite = 0; isite < StdI->nsite; isite++)*/
744 
745  fp = fopen("orbitalidxpara.def", "w");
746  fprintf(fp, "=============================================\n");
747  fprintf(fp, "NOrbitalIdx %10d\n", NOrbGC);
748  fprintf(fp, "ComplexType %10d\n", StdI->ComplexType);
749  fprintf(fp, "=============================================\n");
750  fprintf(fp, "=============================================\n");
751 
752  for (isite = 0; isite < StdI->nsite; isite++) {
753  for (jsite = 0; jsite < StdI->nsite; jsite++) {
754  if (isite >= jsite) continue;
755  if (StdI->AntiPeriod[0] == 1 || StdI->AntiPeriod[1] == 1 || StdI->AntiPeriod[2] == 1)
756  fprintf(fp, "%5d %5d %5d %5d\n", isite, jsite, OrbGC[isite][jsite], AntiOrbGC[isite][jsite]);
757  else
758  fprintf(fp, "%5d %5d %5d\n", isite, jsite, OrbGC[isite][jsite]);
759  }/*for (jsite = 0; jsite < isite; jsite++)*/
760  }/*for (isite = 0; isite < StdI->nsite; isite++)*/
761 
762  for (iOrbGC = 0; iOrbGC < NOrbGC; iOrbGC++)
763  fprintf(fp, "%5d %5d\n", iOrbGC, 1);
764 
765  fflush(fp);
766  fclose(fp);
767  fprintf(stdout, " orbitalidxpara.def is written.\n");
768 
769  for (isite = 0; isite < StdI->nsite; isite++) {
770  free(OrbGC[isite]);
771  free(AntiOrbGC[isite]);
772  }
773  free(OrbGC);
774  free(AntiOrbGC);
775 }/*static void PrintOrbPara*/
779 static void PrintGutzwiller(struct StdIntList *StdI)
780 {
781  FILE *fp;
782  int iCell, isite, jsite, NGutzwiller, iGutz;
783  int *Gutz;
784 
785  Gutz = (int *)malloc(sizeof(int) * StdI->nsite);
786 
787  if (std::abs(StdI->NMPTrans) == 1 || StdI->NMPTrans == StdI->NaN_i) {
788  if (strcmp(StdI->model, "hubbard") == 0) NGutzwiller = 0;
789  else NGutzwiller = -1;
790 
791  for (isite = 0; isite < StdI->nsite; isite++) Gutz[isite] = StdI->Orb[isite][isite];
792 
793  for (isite = 0; isite < StdI->nsite; isite++) {
794  /*
795  For Local spin
796  */
797  if (StdI->locspinflag[isite] != 0) {
798  Gutz[isite] = -1;
799  continue;
800  }
801 
802  if (Gutz[isite] >= 0) {
803  iGutz = Gutz[isite];
804  NGutzwiller -= 1;
805  for (jsite = 0; jsite < StdI->nsite; jsite++) {
806  if (Gutz[jsite] == iGutz)
807  Gutz[jsite] = NGutzwiller;
808  }/*for (jsite = 0; jsite < StdI->nsite; jsite++)*/
809  }/*if (Gutz[isite] >= 0)*/
810  }/*for (isite = 0; isite < StdI->nsite; isite++)*/
811 
812  NGutzwiller = -NGutzwiller;
813  for (isite = 0; isite < StdI->nsite; isite++) {
814  Gutz[isite] = -1 - Gutz[isite];
815  }/*for (isite = 0; isite < StdI->nsite; isite++)*/
816  }/*if (std::abs(StdI->NMPTrans) == 1)*/
817  else {
818  if (strcmp(StdI->model, "hubbard") == 0) NGutzwiller = StdI->NsiteUC;
819  else if (strcmp(StdI->model, "spin") == 0) NGutzwiller = 1;
820  else NGutzwiller = StdI->NsiteUC + 1;
821 
822  for (iCell = 0; iCell < StdI->NCell; iCell++) {
823  for (isite = 0; isite < StdI->NsiteUC; isite++) {
824  if (strcmp(StdI->model, "hubbard") == 0)
825  Gutz[isite + StdI->NsiteUC*iCell] = isite;
826  else if (strcmp(StdI->model, "spin") == 0)
827  Gutz[isite + StdI->NsiteUC*iCell] = 0;
828  else {
829  Gutz[isite + StdI->NsiteUC*iCell] = 0;
830  Gutz[isite + StdI->NsiteUC*(iCell + StdI->NCell)] = isite + 1;
831  }
832  }/*for (isite = 0; isite < StdI->NsiteUC; isite++)*/
833  }/*for (iCell = 0; iCell < StdI->NCell; iCell++)*/
834  }/*if (std::abs(StdI->NMPTrans) != 1)*/
835 
836  fp = fopen("gutzwilleridx.def", "w");
837  fprintf(fp, "=============================================\n");
838  fprintf(fp, "NGutzwillerIdx %10d\n", NGutzwiller);
839  fprintf(fp, "ComplexType %10d\n", 0);
840  fprintf(fp, "=============================================\n");
841  fprintf(fp, "=============================================\n");
842 
843  for (isite = 0; isite < StdI->nsite; isite++)
844  fprintf(fp, "%5d %5d\n", isite, Gutz[isite]);
845 
846  for (iGutz = 0; iGutz < NGutzwiller; iGutz++) {
847  if (strcmp(StdI->model, "hubbard") == 0 || iGutz > 0)
848  fprintf(fp, "%5d %5d\n", iGutz, 1);
849  else
850  fprintf(fp, "%5d %5d\n", iGutz, 0);
851  }/*for (iGutz = 0; iGutz < NGutzwiller; iGutz++)*/
852  fflush(fp);
853  fclose(fp);
854  fprintf(stdout, " gutzwilleridx.def is written.\n");
855 
856  free(Gutz);
857 }/*static void PrintGutzwiller*/
858 #endif
859 
864 static void StdFace_ResetVals(struct StdIntList *StdI) {
865  int i, j;
866  double NaN_d;
867  /*
868  NaN is used for not inputed variable
869  */
870  NaN_d = 0.0 / 0.0;
871  StdI->NaN_i = 2147483647;
872  StdI->pi = acos(-1.0);
873 
874  StdI->a = NaN_d;
875  for (i = 0; i < 3; i++) StdI->length[i] = NaN_d;
876  for (i = 0; i < 3; i++)
877  for (j = 0; j < 3; j++)
878  StdI->box[i][j] = StdI->NaN_i;
879  StdI->Gamma = NaN_d;
880  StdI->h = NaN_d;
881  StdI->Height = StdI->NaN_i;
882  StdI->JAll = NaN_d;
883  StdI->JpAll = NaN_d;
884  StdI->JpAll = NaN_d;
885  StdI->JppAll = NaN_d;
886  StdI->J0All = NaN_d;
887  StdI->J0pAll = NaN_d;
888  StdI->J0ppAll = NaN_d;
889  StdI->J1All = NaN_d;
890  StdI->J1pAll = NaN_d;
891  StdI->J1ppAll = NaN_d;
892  StdI->J2All = NaN_d;
893  StdI->J2pAll = NaN_d;
894  StdI->J2ppAll = NaN_d;
895  for (i = 0; i < 3; i++) {
896  for (j = 0; j < 3; j++) {
897  StdI->J[i][j] = NaN_d;
898  StdI->Jp[i][j] = NaN_d;
899  StdI->Jpp[i][j] = NaN_d;
900  StdI->J0[i][j] = NaN_d;
901  StdI->J0p[i][j] = NaN_d;
902  StdI->J0pp[i][j] = NaN_d;
903  StdI->J1[i][j] = NaN_d;
904  StdI->J1p[i][j] = NaN_d;
905  StdI->J1pp[i][j] = NaN_d;
906  StdI->J2[i][j] = NaN_d;
907  StdI->J2p[i][j] = NaN_d;
908  StdI->J2pp[i][j] = NaN_d;
909  StdI->D[i][j] = 0.0;
910  }
911  }
912  StdI->D[2][2] = NaN_d;
913  StdI->K = NaN_d;
914  StdI->L = StdI->NaN_i;
915  for (i = 0; i < 3; i++)
916  for (j = 0; j < 3; j++)
917  StdI->direct[i][j] = NaN_d;
918  StdI->mu = NaN_d;
919  StdI->S2 = StdI->NaN_i;
920  StdI->t = NaN_d;
921  StdI->tp = NaN_d;
922  StdI->tpp = NaN_d;
923  StdI->t0 = NaN_d;
924  StdI->t0p = NaN_d;
925  StdI->t0pp = NaN_d;
926  StdI->t1 = NaN_d;
927  StdI->t1p = NaN_d;
928  StdI->t1pp = NaN_d;
929  StdI->t2 = NaN_d;
930  StdI->t2p = NaN_d;
931  StdI->t2pp = NaN_d;
932  StdI->U = NaN_d;
933  StdI->V = NaN_d;
934  StdI->Vp = NaN_d;
935  StdI->Vpp = NaN_d;
936  StdI->V0 = NaN_d;
937  StdI->V0p = NaN_d;
938  StdI->V0pp = NaN_d;
939  StdI->V1 = NaN_d;
940  StdI->V1p = NaN_d;
941  StdI->V1pp = NaN_d;
942  StdI->V2 = NaN_d;
943  StdI->V2p = NaN_d;
944  StdI->V2pp = NaN_d;
945  StdI->W = StdI->NaN_i;
946  for (i = 0; i < 3; i++)StdI->phase[i] = NaN_d;
947  StdI->pi180 = StdI->pi / 180.0;
948 
949  StdI->nelec = StdI->NaN_i;
950  StdI->Sz2 = StdI->NaN_i;
951  strcpy(StdI->model, "****\0");
952  strcpy(StdI->lattice, "****\0");
953  strcpy(StdI->outputmode, "****\0");
954  strcpy(StdI->CDataFileHead, "****\0");
955  StdI->cutoff_t = NaN_d;
956  StdI->cutoff_u = NaN_d;
957  StdI->cutoff_j = NaN_d;
958  StdI->cutoff_length_t = NaN_d;
959  StdI->cutoff_length_U = NaN_d;
960  StdI->cutoff_length_J = NaN_d;
961  for (i = 0; i < 3; i++)StdI->cutoff_tR[i] = StdI->NaN_i;
962  for (i = 0; i < 3; i++)StdI->cutoff_UR[i] = StdI->NaN_i;
963  for (i = 0; i < 3; i++)StdI->cutoff_JR[i] = StdI->NaN_i;
964  StdI->double_counting = StdI->NaN_i;
965 #if defined(_HPhi)
966  StdI->LargeValue = NaN_d;
967  StdI->OmegaMax = NaN_d;
968  StdI->OmegaMin = NaN_d;
969  StdI->OmegaIm = NaN_d;
970  StdI->Nomega = StdI->NaN_i;
971  for (i = 0; i < 3; i++)StdI->SpectrumQ[i] = NaN_d;
972  strcpy(StdI->method, "****\0");
973  strcpy(StdI->Restart, "****\0");
974  strcpy(StdI->EigenVecIO, "****\0");
975  strcpy(StdI->InitialVecType, "****\0");
976  strcpy(StdI->CalcSpec, "****\0");
977  strcpy(StdI->SpectrumType, "****\0");
978  StdI->FlgTemp = 1;
979  StdI->Lanczos_max = StdI->NaN_i;
980  StdI->initial_iv = StdI->NaN_i;
981  StdI->nvec = StdI->NaN_i;
982  StdI->exct = StdI->NaN_i;
983  StdI->LanczosEps = StdI->NaN_i;
984  StdI->LanczosTarget = StdI->NaN_i;
985  StdI->NumAve = StdI->NaN_i;
986  StdI->ExpecInterval = StdI->NaN_i;
987  StdI->dt = NaN_d;
988  StdI->tdump = NaN_d;
989  StdI->tshift = NaN_d;
990  StdI->freq = NaN_d;
991  StdI->Uquench = NaN_d;
992  for (i = 0; i < 3; i++)StdI->VecPot[i] = NaN_d;;
993  strcpy(StdI->PumpType, "****\0");
994  StdI->ExpandCoef = StdI->NaN_i;
995 #elif defined(_mVMC)
996  strcpy(StdI->CParaFileHead, "****\0");
997  StdI->NVMCCalMode = StdI->NaN_i;
998  StdI->NLanczosMode = StdI->NaN_i;
999  StdI->NDataIdxStart = StdI->NaN_i;
1000  StdI->NDataQtySmp = StdI->NaN_i;
1001  StdI->NSPGaussLeg = StdI->NaN_i;
1002  StdI->NSPStot = StdI->NaN_i;
1003  StdI->NMPTrans = StdI->NaN_i;
1004  StdI->NSROptItrStep = StdI->NaN_i;
1005  StdI->NSROptItrSmp = StdI->NaN_i;
1006  StdI->DSROptRedCut = NaN_d;
1007  StdI->DSROptStaDel = NaN_d;
1008  StdI->DSROptStepDt = NaN_d;
1009  StdI->NVMCWarmUp = StdI->NaN_i;
1010  StdI->NVMCInterval = StdI->NaN_i;
1011  StdI->NVMCSample = StdI->NaN_i;
1012  StdI->NExUpdatePath = StdI->NaN_i;
1013  StdI->RndSeed = StdI->NaN_i;
1014  StdI->NSplitSize = StdI->NaN_i;
1015  StdI->NStore = StdI->NaN_i;
1016  StdI->NSRCG = StdI->NaN_i;
1017  StdI->ComplexType = StdI->NaN_i;
1018  for (i = 0; i < 3; i++)
1019  for (j = 0; j < 3; j++)
1020  StdI->boxsub[i][j] = StdI->NaN_i;
1021  StdI->Hsub = StdI->NaN_i;
1022  StdI->Lsub = StdI->NaN_i;
1023  StdI->Wsub = StdI->NaN_i;
1024 #endif
1025 }/*static void StdFace_ResetVals*/
1026 /*
1027 @brief Make all characters lower
1028 @author Mitsuaki Kawamura (The University of Tokyo)
1029 */
1030 static void Text2Lower(char *value
1031 ){
1032  char value2;
1033  int valuelen, ii;
1034 
1035  valuelen = strlen(value);
1036  for (ii = 0; ii < valuelen; ii++) {
1037  value2 = tolower(value[ii]);
1038  value[ii] = value2;
1039  }
1040 }/*static void Text2Lower*/
1045 static void TrimSpaceQuote(char *value
1046 ){
1047  char value2[256];
1048  int valuelen, valuelen2, ii;
1049 
1050  valuelen = strlen(value);
1051  valuelen2 = 0;
1052  for (ii = 0; ii < valuelen; ii++){
1053  if (value[ii] != ' ' &&
1054  value[ii] != ':' &&
1055  value[ii] != ';' &&
1056  value[ii] != '\"' &&
1057  value[ii] != '\b' &&
1058  value[ii] != '\\' &&
1059  value[ii] != '\v' &&
1060  value[ii] != '\n' &&
1061  value[ii] != '\0'){
1062  value2[valuelen2] = value[ii];
1063  valuelen2++;
1064  }
1065  }
1066 
1067  strncpy(value, value2, valuelen2);
1068  value[valuelen2] = '\0';
1069 
1070 }/*static void TrimSpaceQuote*/
1077  char *keyword,
1078  char *valuestring,
1079  char *value
1080 )
1081 {
1082  if (strcmp(value, "****") != 0){
1083  fprintf(stdout, "ERROR ! Keyword %s is duplicated ! \n", keyword);
1084  StdFace_exit(-1);
1085  }
1086  else{
1087  strcpy(value, valuestring);
1088  }
1089 }/*static void StoreWithCheckDup_s*/
1096  char *keyword,
1097  char *valuestring,
1098  char *value
1099 )
1100 {
1101  if (strcmp(value, "****") != 0) {
1102  fprintf(stdout, "ERROR ! Keyword %s is duplicated ! \n", keyword);
1103  StdFace_exit(-1);
1104  }
1105  else {
1106  strcpy(value, valuestring);
1107  Text2Lower(value);
1108  }
1109 }/*static void StoreWithCheckDup_sl*/
1116  char *keyword,
1117  char *valuestring,
1118  int *value
1119 )
1120 {
1121  int NaN_i = 2147483647;
1122 
1123  if (*value != NaN_i){
1124  fprintf(stdout, "ERROR ! Keyword %s is duplicated ! \n", keyword);
1125  StdFace_exit(-1);
1126  }
1127  else{
1128  sscanf(valuestring, "%d", value);
1129  }
1130 }/*static void StoreWithCheckDup_i*/
1137  char *keyword,
1138  char *valuestring,
1139  double *value
1140 )
1141 {
1142  if (std::isnan(*value) == 0){
1143  fprintf(stdout, "ERROR ! Keyword %s is duplicated ! \n", keyword);
1144  StdFace_exit(-1);
1145  }
1146  else{
1147  sscanf(valuestring, "%lf", value);
1148  }
1149 }/*static void StoreWithCheckDup_d*/
1156  char *keyword,
1157  char *valuestring,
1158  std::complex<double> *value
1159 )
1160 {
1161  int num;
1162  char *valuestring_r, *valuestring_i;
1163  double value_r, value_i;
1164 
1165  if (std::isnan(real(*value)) == 0) {
1166  fprintf(stdout, "ERROR ! Keyword %s is duplicated ! \n", keyword);
1167  StdFace_exit(-1);
1168  }
1169  else {
1170 
1171  if (valuestring[0] == ',') {
1172  valuestring_r = NULL;
1173  valuestring_i = strtok(valuestring, ",");
1174  }
1175  else {
1176  valuestring_r = strtok(valuestring, ",");
1177  valuestring_i = strtok(NULL, ",");
1178  }
1179 
1180  if (valuestring_r == NULL) {
1181  *value = 0.0;
1182  }
1183  else {
1184  num = sscanf(valuestring_r, "%lf", &value_r);
1185  if (num == 1) *value = value_r;
1186  else *value = 0.0;
1187  }
1188 
1189  if (valuestring_i == NULL) {
1190  *value += std::complex<double>(0.0, 0.0);
1191  }
1192  else {
1193  num = sscanf(valuestring_i, "%lf", &value_i);
1194  if (num == 1) *value += std::complex<double>(0.0, value_i);
1195  else *value += std::complex<double>(0.0, 0.0);
1196  }
1197  }
1198 }/*static void StoreWithCheckDup_c*/
1203 static void PrintLocSpin(struct StdIntList *StdI) {
1204  FILE *fp;
1205  int isite, nlocspin;
1206 
1207  nlocspin = 0;
1208  for (isite = 0; isite < StdI->nsite; isite++)
1209  if (StdI->locspinflag[isite] != 0) nlocspin = nlocspin + 1;
1210 
1211  fp = fopen("locspn.def", "w");
1212  fprintf(fp, "================================ \n");
1213  fprintf(fp, "NlocalSpin %5d \n", nlocspin);
1214  fprintf(fp, "================================ \n");
1215  fprintf(fp, "========i_0LocSpn_1IteElc ====== \n");
1216  fprintf(fp, "================================ \n");
1217 
1218  for (isite = 0; isite < StdI->nsite; isite++)
1219  fprintf(fp, "%5d %5d\n", isite, StdI->locspinflag[isite]);
1220 
1221  fflush(fp);
1222  fclose(fp);
1223  fprintf(stdout, " locspn.def is written.\n");
1224 }/*static void PrintLocSpin*/
1229 static void PrintTrans(struct StdIntList *StdI){
1230  FILE *fp;
1231  int jtrans, ktrans, ntrans0;
1232 
1233  for (jtrans = 0; jtrans < StdI->ntrans; jtrans++){
1234  for (ktrans = jtrans + 1; ktrans < StdI->ntrans; ktrans++){
1235  if (StdI->transindx[jtrans][0] == StdI->transindx[ktrans][0]
1236  && StdI->transindx[jtrans][1] == StdI->transindx[ktrans][1]
1237  && StdI->transindx[jtrans][2] == StdI->transindx[ktrans][2]
1238  && StdI->transindx[jtrans][3] == StdI->transindx[ktrans][3]){
1239  StdI->trans[jtrans] = StdI->trans[jtrans] + StdI->trans[ktrans];
1240  StdI->trans[ktrans] = 0.0;
1241  }
1242  }/*for (ktrans = jtrans + 1; ktrans < StdI->ntrans; ktrans++)*/
1243  }/*for (jtrans = 0; jtrans < StdI->ntrans; jtrans++)*/
1244 
1245  ntrans0 = 0;
1246  for (ktrans = 0; ktrans < StdI->ntrans; ktrans++){
1247  if (std::abs(StdI->trans[ktrans]) > 0.000001) ntrans0 = ntrans0 + 1;
1248  }
1249 
1250  fp = fopen("trans.def", "w");
1251  fprintf(fp, "======================== \n");
1252  fprintf(fp, "NTransfer %7d \n", ntrans0);
1253  fprintf(fp, "======================== \n");
1254  fprintf(fp, "========i_j_s_tijs====== \n");
1255  fprintf(fp, "======================== \n");
1256 
1257  ntrans0 = 0;
1258  for (ktrans = 0; ktrans < StdI->ntrans; ktrans++) {
1259  if (std::abs(StdI->trans[ktrans]) > 0.000001)
1260  fprintf(fp, "%5d %5d %5d %5d %25.15f %25.15f\n",
1261  StdI->transindx[ktrans][0], StdI->transindx[ktrans][1],
1262  StdI->transindx[ktrans][2], StdI->transindx[ktrans][3],
1263  real(StdI->trans[ktrans]), imag(StdI->trans[ktrans]));
1264  }
1265 
1266  fflush(fp);
1267  fclose(fp);
1268  fprintf(stdout, " trans.def is written.\n");
1269 }/*static void PrintTrans*/
1274 static void PrintNamelist(struct StdIntList *StdI){
1275  FILE *fp;
1276 
1277  fp = fopen("namelist.def", "w");
1278  fprintf( fp, " ModPara modpara.def\n");
1279  fprintf( fp, " LocSpin locspn.def\n");
1280  fprintf( fp, " Trans trans.def\n");
1281  if (StdI->LCintra == 1) fprintf( fp, " CoulombIntra coulombintra.def\n");
1282  if (StdI->LCinter == 1) fprintf( fp, " CoulombInter coulombinter.def\n");
1283  if (StdI->LHund == 1)fprintf( fp, " Hund hund.def\n");
1284  if (StdI->LEx == 1)fprintf( fp, " Exchange exchange.def\n");
1285  if (StdI->LPairLift == 1)fprintf(fp, " PairLift pairlift.def\n");
1286  if (StdI->LPairHopp == 1)fprintf(fp, " PairHop pairhopp.def\n");
1287  if (StdI->Lintr == 1)fprintf( fp, " InterAll interall.def\n");
1288  if (StdI->ioutputmode != 0) {
1289  fprintf( fp, " OneBodyG greenone.def\n");
1290  fprintf( fp, " TwoBodyG greentwo.def\n");
1291  }
1292 #if defined(_HPhi)
1293  fprintf( fp, " CalcMod calcmod.def\n");
1294  if(StdI->SpectrumBody == 1)
1295  fprintf( fp, "SingleExcitation single.def\n");
1296  else fprintf( fp, " PairExcitation pair.def\n");
1297  if (strcmp(StdI->method, "timeevolution") == 0) {
1298  if (StdI->PumpBody == 1)
1299  fprintf(fp, " TEOneBody teone.def\n");
1300  else if (StdI->PumpBody == 2)
1301  fprintf(fp, " TETwoBody tetwo.def\n");
1302  }/*if (strcmp(StdI->method, "timeevolution") == 0)*/
1303  fprintf( fp, " SpectrumVec %s_eigenvec_0\n",
1304  StdI->CDataFileHead);
1305  if (StdI->lBoost == 1) fprintf( fp, " Boost boost.def\n");
1306 #elif defined(_mVMC)
1307  fprintf( fp, " Gutzwiller gutzwilleridx.def\n");
1308  fprintf( fp, " Jastrow jastrowidx.def\n");
1309  fprintf( fp, " Orbital orbitalidx.def\n");
1310  if (StdI->lGC == 1 || (StdI->Sz2 != 0 && StdI->Sz2 != StdI->NaN_i))
1311  fprintf(fp, " OrbitalParallel orbitalidxpara.def\n");
1312  fprintf( fp, " TransSym qptransidx.def\n");
1313  if(strcmp(StdI->lattice, "wannier90") == 0)
1314  fprintf(fp, " Initial initial.def\n");
1315 #endif
1316 
1317  fflush(fp);
1318  fclose(fp);
1319  fprintf(stdout, " namelist.def is written.\n");
1320 }/*static void PrintNamelist*/
1325 static void PrintModPara(struct StdIntList *StdI)
1326 {
1327  FILE *fp;
1328 
1329  fp = fopen("modpara.def", "w");
1330  fprintf(fp, "--------------------\n");
1331  fprintf(fp, "Model_Parameters 0\n");
1332  fprintf(fp, "--------------------\n");
1333 #if defined(_HPhi)
1334  fprintf(fp, "HPhi_Cal_Parameters\n");
1335  fprintf(fp, "--------------------\n");
1336  fprintf(fp, "CDataFileHead %s\n", StdI->CDataFileHead);
1337  fprintf(fp, "CParaFileHead zqp\n");
1338  fprintf(fp, "--------------------\n");
1339  fprintf(fp, "Nsite %-5d\n", StdI->nsite);
1340  if (StdI->Sz2 != StdI->NaN_i) fprintf(fp, "2Sz %-5d\n", StdI->Sz2);
1341  if (StdI->nelec != StdI->NaN_i) fprintf(fp, "Ncond %-5d\n", StdI->nelec);
1342  fprintf(fp, "Lanczos_max %-5d\n", StdI->Lanczos_max);
1343  fprintf(fp, "initial_iv %-5d\n", StdI->initial_iv);
1344  if(StdI->nvec != StdI->NaN_i) fprintf(fp, "nvec %-5d\n", StdI->nvec);
1345  fprintf(fp, "exct %-5d\n", StdI->exct);
1346  fprintf(fp, "LanczosEps %-5d\n", StdI->LanczosEps);
1347  fprintf(fp, "LanczosTarget %-5d\n", StdI->LanczosTarget);
1348  fprintf(fp, "LargeValue %-25.15e\n", StdI->LargeValue);
1349  fprintf(fp, "NumAve %-5d\n", StdI->NumAve);
1350  fprintf(fp, "ExpecInterval %-5d\n", StdI->ExpecInterval);
1351  fprintf(fp, "NOmega %-5d\n", StdI->Nomega);
1352  fprintf(fp, "OmegaMax %-25.15e %-25.15e\n", StdI->OmegaMax, StdI->OmegaIm);
1353  fprintf(fp, "OmegaMin %-25.15e %-25.15e\n", StdI->OmegaMin, StdI->OmegaIm);
1354  fprintf(fp, "OmegaOrg 0.0 0.0\n");
1355  if (strcmp(StdI->method, "timeevolution") == 0)
1356  fprintf(fp, "ExpandCoef %-5d\n", StdI->ExpandCoef);
1357 #elif defined(_mVMC)
1358  fprintf(fp, "VMC_Cal_Parameters\n");
1359  fprintf(fp, "--------------------\n");
1360  fprintf(fp, "CDataFileHead %s\n", StdI->CDataFileHead);
1361  fprintf(fp, "CParaFileHead %s\n", StdI->CParaFileHead);
1362  fprintf(fp, "--------------------\n");
1363  fprintf(fp, "NVMCCalMode %d\n", StdI->NVMCCalMode);
1364  /*fprintf(fp, "NLanczosMode %d\n", StdI->NLanczosMode);*/
1365  fprintf(fp, "--------------------\n");
1366  fprintf(fp, "NDataIdxStart %d\n", StdI->NDataIdxStart);
1367  fprintf(fp, "NDataQtySmp %d\n", StdI->NDataQtySmp);
1368  fprintf(fp, "--------------------\n");
1369  fprintf(fp, "Nsite %d\n", StdI->nsite);
1370  fprintf(fp, "Ncond %-5d\n", StdI->nelec);
1371  if (StdI->Sz2 != StdI->NaN_i)
1372  fprintf(fp, "2Sz %d\n", StdI->Sz2);
1373  if (StdI->NSPGaussLeg != StdI->NaN_i)
1374  fprintf(fp, "NSPGaussLeg %d\n", StdI->NSPGaussLeg);
1375  if (StdI->NSPStot != StdI->NaN_i)
1376  fprintf(fp, "NSPStot %d\n", StdI->NSPStot);
1377  fprintf(fp, "NMPTrans %d\n", StdI->NMPTrans);
1378  fprintf(fp, "NSROptItrStep %d\n", StdI->NSROptItrStep);
1379  fprintf(fp, "NSROptItrSmp %d\n", StdI->NSROptItrSmp);
1380  fprintf(fp, "DSROptRedCut %.10f\n", StdI->DSROptRedCut);
1381  fprintf(fp, "DSROptStaDel %.10f\n", StdI->DSROptStaDel);
1382  fprintf(fp, "DSROptStepDt %.10f\n", StdI->DSROptStepDt);
1383  fprintf(fp, "NVMCWarmUp %d\n", StdI->NVMCWarmUp);
1384  fprintf(fp, "NVMCInterval %d\n", StdI->NVMCInterval);
1385  fprintf(fp, "NVMCSample %d\n", StdI->NVMCSample);
1386  fprintf(fp, "NExUpdatePath %d\n", StdI->NExUpdatePath);
1387  fprintf(fp, "RndSeed %d\n", StdI->RndSeed);
1388  fprintf(fp, "NSplitSize %d\n", StdI->NSplitSize);
1389  fprintf(fp, "NStore %d\n", StdI->NStore);
1390  fprintf(fp, "NSRCG %d\n", StdI->NSRCG);
1391 #endif
1392 
1393  fflush(fp);
1394  fclose(fp);
1395  fprintf(stdout, " modpara.def is written.\n");
1396 }/*static void PrintModPara*/
1401 static void Print1Green(struct StdIntList *StdI)
1402 {
1403  FILE *fp;
1404  int ngreen, igreen, store, xkondo;
1405  int isite, jsite, ispin, jspin, SiMax, SjMax;
1406  int **greenindx;
1407  /*
1408  Set Indices of correlation functions
1409  */
1410  ngreen = 0;
1411  if (StdI->ioutputmode != 0) {
1412  for (store = 0; store < 2; store++) {
1413 
1414  if (store == 1) {
1415  greenindx = (int **)malloc(sizeof(int*) * (ngreen + 1));
1416  for (igreen = 0; igreen < ngreen; igreen++) {
1417  greenindx[igreen] = (int *)malloc(sizeof(int) * 4);
1418  }
1419  ngreen = 0;
1420  }/*if (store == 1)*/
1421 
1422  if (strcmp(StdI->model, "kondo") == 0) xkondo = 2;
1423  else xkondo = 1;
1424 
1425  if (StdI->ioutputmode == 1) {
1426  for (isite = 0; isite < StdI->NsiteUC*xkondo; isite++) {
1427 
1428  if (isite >= StdI->NsiteUC) isite += StdI->nsite / 2;
1429 
1430  if (StdI->locspinflag[isite] == 0) SiMax = 1;
1431  else SiMax = StdI->locspinflag[isite];
1432 
1433  for (ispin = 0; ispin <= SiMax; ispin++) {
1434  for (jsite = 0; jsite < StdI->nsite; jsite++) {
1435 
1436  if (StdI->locspinflag[jsite] == 0) SjMax = 1;
1437  else SjMax = StdI->locspinflag[jsite];
1438 
1439  for (jspin = 0; jspin <= SjMax; jspin++) {
1440 
1441  if (isite != jsite &&
1442  (StdI->locspinflag[isite] != 0 && StdI->locspinflag[jsite] != 0)) continue;
1443 
1444  if (ispin == jspin){
1445  if (store == 1) {
1446  greenindx[ngreen][0] = isite;
1447  greenindx[ngreen][1] = ispin;
1448  greenindx[ngreen][2] = jsite;
1449  greenindx[ngreen][3] = jspin;
1450  }
1451  ngreen++;
1452  }
1453 
1454  }/*for (jspin = 0; jspin <= SjMax; jspin++)*/
1455  }/*for (jsite = 0; jsite < StdI->nsite; jsite++)*/
1456  }/*for (ispin = 0; ispin <= SiMax; ispin++)*/
1457  }/*for (isite = 0; isite < StdI->nsite; isite++)*/
1458  }/*if (StdI->ioutputmode == 1)*/
1459  else {
1460  for (isite = 0; isite < StdI->nsite; isite++) {
1461 
1462  if (StdI->locspinflag[isite] == 0) SiMax = 1;
1463  else SiMax = StdI->locspinflag[isite];
1464 
1465  for (ispin = 0; ispin <= SiMax; ispin++) {
1466  for (jsite = 0; jsite < StdI->nsite; jsite++) {
1467 
1468  if (StdI->locspinflag[jsite] == 0) SjMax = 1;
1469  else SjMax = StdI->locspinflag[jsite];
1470 
1471  for (jspin = 0; jspin <= SjMax; jspin++) {
1472 
1473  if (isite != jsite &&
1474  (StdI->locspinflag[isite] != 0 && StdI->locspinflag[jsite] != 0)) continue;
1475 
1476  if (store == 1) {
1477  greenindx[ngreen][0] = isite;
1478  greenindx[ngreen][1] = ispin;
1479  greenindx[ngreen][2] = jsite;
1480  greenindx[ngreen][3] = jspin;
1481  }
1482  ngreen++;
1483 
1484  }/*for (jspin = 0; jspin <= SjMax; jspin++)*/
1485  }/*for (jsite = 0; jsite < StdI->nsite; jsite++)*/
1486  }/*for (ispin = 0; ispin <= SiMax; ispin++)*/
1487  }/*for (isite = 0; isite < StdI->nsite; isite++)*/
1488  }/*if (StdI->ioutputmode == 2)*/
1489  }/*if (StdI->ioutputmode != 0)*/
1490 
1491  fp = fopen("greenone.def", "w");
1492  fprintf(fp, "===============================\n");
1493  fprintf(fp, "NCisAjs %10d\n", ngreen);
1494  fprintf(fp, "===============================\n");
1495  fprintf(fp, "======== Green functions ======\n");
1496  fprintf(fp, "===============================\n");
1497  for (igreen = 0; igreen < ngreen; igreen++) {
1498  fprintf(fp, "%5d %5d %5d %5d\n",
1499  greenindx[igreen][0], greenindx[igreen][1], greenindx[igreen][2], greenindx[igreen][3]);
1500  }
1501  fflush(fp);
1502  fclose(fp);
1503 
1504  fprintf(stdout, " greenone.def is written.\n");
1505 
1506  for (igreen = 0; igreen < ngreen; igreen++) {
1507  free(greenindx[igreen]);
1508  }
1509  free(greenindx);
1510 
1511  }/*if (StdI->ioutputmode != 0) */
1512 }/*static void Print1Green*/
1517 static void Print2Green(struct StdIntList *StdI) {
1518  FILE *fp;
1519  int ngreen, store, igreen, xkondo;
1520  int site1, site2, site3, site4;
1521  int spin1, spin2, spin3, spin4;
1522  int S1Max, S2Max, S3Max, S4Max;
1523  int **greenindx;
1524  /*
1525  Set Indices of correlation functions
1526  */
1527  ngreen = 0;
1528  if (StdI->ioutputmode == 1) {
1529  for (store = 0; store < 2; store++) {
1530 
1531  if (store == 1) {
1532  greenindx = (int **)malloc(sizeof(int*) * (ngreen + 1));
1533  for (igreen = 0; igreen < ngreen; igreen++)
1534  greenindx[igreen] = (int *)malloc(sizeof(int) * 8);
1535  ngreen = 0;
1536  }/*if (store == 1)*/
1537 
1538  if (strcmp(StdI->model, "kondo") == 0) xkondo = 2;
1539  else xkondo = 1;
1540 
1541  for (site1 = 0; site1 < StdI->NsiteUC*xkondo; site1++) {
1542 
1543  if (site1 >= StdI->NsiteUC) site1 += StdI->nsite / 2;
1544 
1545  if (StdI->locspinflag[site1] == 0) S1Max = 1;
1546  else S1Max = StdI->locspinflag[site1];
1547  for (spin1 = 0; spin1 <= S1Max; spin1++) {
1548  for (spin2 = 0; spin2 <= S1Max; spin2++) {
1549 
1550  for (site3 = 0; site3 < StdI->nsite; site3++) {
1551 
1552  if (StdI->locspinflag[site3] == 0) S3Max = 1;
1553  else S3Max = StdI->locspinflag[site3];
1554  for (spin3 = 0; spin3 <= S3Max; spin3++) {
1555  for (spin4 = 0; spin4 <= S3Max; spin4++) {
1556 
1557  if (spin1 - spin2 + spin3 - spin4 == 0) {
1558  if (store == 1) {
1559 #if defined(_mVMC)
1560  if (spin1 != spin2 || spin3 != spin4)
1561  {
1562  greenindx[ngreen][0] = site1;
1563  greenindx[ngreen][1] = spin1;
1564  greenindx[ngreen][2] = site3;
1565  greenindx[ngreen][3] = spin4;
1566  greenindx[ngreen][4] = site3;
1567  greenindx[ngreen][5] = spin3;
1568  greenindx[ngreen][6] = site1;
1569  greenindx[ngreen][7] = spin2;
1570  }
1571  else
1572 #endif
1573  {
1574  greenindx[ngreen][0] = site1;
1575  greenindx[ngreen][1] = spin1;
1576  greenindx[ngreen][2] = site1;
1577  greenindx[ngreen][3] = spin2;
1578  greenindx[ngreen][4] = site3;
1579  greenindx[ngreen][5] = spin3;
1580  greenindx[ngreen][6] = site3;
1581  greenindx[ngreen][7] = spin4;
1582  }
1583  }/*if (store == 1)*/
1584  ngreen++;
1585  }/*if (spin1 - spin2 + spin3 - spin4 == 0)*/
1586 
1587  }/*for (spin4 = 0; spin4 <= S3Max; spin4++)*/
1588  }/*for (spin3 = 0; spin3 <= S3Max; spin3++*/
1589  }/*for (site3 = 0; site3 < StdI->nsite; site3++)*/
1590  }/*for (spin2 = 0; spin2 <= S1Max; spin2++)*/
1591  }/*for (spin1 = 0; spin1 <= S1Max; spin1++)*/
1592  }/*for (site1 = 0; site1 < StdI->nsite; site1++)*/
1593 
1594  }/*for (store = 0; store < 2; store++)*/
1595  }/*if (StdI->ioutputmode == 1)*/
1596  else if (StdI->ioutputmode == 2) {
1597  for (store = 0; store < 2; store++) {
1598 
1599  if (store == 1) {
1600  greenindx = (int **)malloc(sizeof(int*) * (ngreen + 1));
1601  for (igreen = 0; igreen < ngreen; igreen++)
1602  greenindx[igreen] = (int *)malloc(sizeof(int) * 8);
1603  ngreen = 0;
1604  }/*if (store == 1)*/
1605 
1606  for (site1 = 0; site1 < StdI->nsite; site1++) {
1607 
1608  if (StdI->locspinflag[site1] == 0) S1Max = 1;
1609  else S1Max = StdI->locspinflag[site1];
1610  for (spin1 = 0; spin1 <= S1Max; spin1++) {
1611 
1612  for (site2 = 0; site2 < StdI->nsite; site2++) {
1613 
1614  if (StdI->locspinflag[site1] != 0 && StdI->locspinflag[site2] != 0
1615  && site1 != site2) continue;
1616 
1617  if (StdI->locspinflag[site2] == 0) S2Max = 1;
1618  else S2Max = StdI->locspinflag[site2];
1619  for (spin2 = 0; spin2 <= S2Max; spin2++) {
1620 
1621  for (site3 = 0; site3 < StdI->nsite; site3++) {
1622 
1623  if (StdI->locspinflag[site3] == 0) S3Max = 1;
1624  else S3Max = StdI->locspinflag[site3];
1625  for (spin3 = 0; spin3 <= S3Max; spin3++) {
1626 
1627  for (site4 = 0; site4 < StdI->nsite; site4++) {
1628 
1629  if (StdI->locspinflag[site3] != 0 && StdI->locspinflag[site4] != 0
1630  && site3 != site4) continue;
1631 
1632  if (StdI->locspinflag[site4] == 0) S4Max = 1;
1633  else S4Max = StdI->locspinflag[site4];
1634  for (spin4 = 0; spin4 <= S4Max; spin4++) {
1635 
1636  if (store == 1) {
1637  greenindx[ngreen][0] = site1;
1638  greenindx[ngreen][1] = spin1;
1639  greenindx[ngreen][2] = site2;
1640  greenindx[ngreen][3] = spin2;
1641  greenindx[ngreen][4] = site3;
1642  greenindx[ngreen][5] = spin3;
1643  greenindx[ngreen][6] = site4;
1644  greenindx[ngreen][7] = spin4;
1645  }/*if (store == 1)*/
1646  ngreen++;
1647 
1648  }/*for (spin4 = 0; spin4 <= S4Max; spin4++)*/
1649  }/*for (site4 = 0; site4 < StdI->nsite; site4++)*/
1650  }/*for (spin3 = 0; spin3 <= S3Max; spin3++*/
1651  }/*for (site3 = 0; site3 < StdI->nsite; site3++)*/
1652  }/*for (spin2 = 0; spin2 <= S2Max; spin2++)*/
1653  }/*for (site2 = 0; site2 < StdI->nsite; site2++)*/
1654  }/*for (spin1 = 0; spin1 <= S1Max; spin1++)*/
1655  }/*for (site1 = 0; site1 < StdI->nsite; site1++)*/
1656 
1657  }/*for (store = 0; store < 2; store++)*/
1658  }/*if (StdI->ioutputmode == 2)*/
1659  if (StdI->ioutputmode != 0) {
1660  fp = fopen("greentwo.def", "w");
1661  fprintf(fp, "=============================================\n");
1662  fprintf(fp, "NCisAjsCktAltDC %10d\n", ngreen);
1663  fprintf(fp, "=============================================\n");
1664  fprintf(fp, "======== Green functions for Sq AND Nq ======\n");
1665  fprintf(fp, "=============================================\n");
1666  for (igreen = 0; igreen < ngreen; igreen++) {
1667  fprintf(fp, "%5d %5d %5d %5d %5d %5d %5d %5d\n",
1668  greenindx[igreen][0], greenindx[igreen][1], greenindx[igreen][2], greenindx[igreen][3],
1669  greenindx[igreen][4], greenindx[igreen][5], greenindx[igreen][6], greenindx[igreen][7]);
1670  }
1671  fflush(fp);
1672  fclose(fp);
1673 
1674  fprintf(stdout, " greentwo.def is written.\n");
1675 
1676  for (igreen = 0; igreen < ngreen; igreen++) {
1677  free(greenindx[igreen]);
1678  }
1679  free(greenindx);
1680  }/*if (StdI->ioutputmode != 0)*/
1681 }/*static void Print2Green(struct StdIntList *StdI)*/
1686 static void UnsupportedSystem(
1687  char *model,
1688  char *lattice
1689 )
1690 {
1691  fprintf(stdout, "\nSorry, specified combination, \n");
1692  fprintf(stdout, " MODEL : %s \n", model);
1693  fprintf(stdout, " LATTICE : %s, \n", lattice);
1694  fprintf(stdout, "is unsupported in the STANDARD MODE...\n");
1695  fprintf(stdout, "Please use the EXPART MODE, or write a NEW FUNCTION and post us.\n");
1696  StdFace_exit(-1);
1697 }/*static void UnsupportedSystem*/
1702 static void CheckOutputMode(struct StdIntList *StdI)
1703 {
1704  /*
1705  Form for Correlation function
1706  */
1707  if (strcmp(StdI->outputmode, "non") == 0
1708  || strcmp(StdI->outputmode, "none") == 0
1709  || strcmp(StdI->outputmode, "off") == 0) {
1710  StdI->ioutputmode = 0;
1711  fprintf(stdout, " ioutputmode = %-10d\n", StdI->ioutputmode);
1712  }
1713  else if (strcmp(StdI->outputmode, "cor") == 0
1714  || strcmp(StdI->outputmode, "corr") == 0
1715  || strcmp(StdI->outputmode, "correlation") == 0) {
1716  StdI->ioutputmode = 1;
1717  fprintf(stdout, " ioutputmode = %-10d\n", StdI->ioutputmode);
1718  }
1719  else if (strcmp(StdI->outputmode, "****") == 0) {
1720  StdI->ioutputmode = 1;
1721  fprintf(stdout, " ioutputmode = %-10d ###### DEFAULT VALUE IS USED ######\n", StdI->ioutputmode);
1722  }
1723  else if (strcmp(StdI->outputmode, "raw") == 0
1724  || strcmp(StdI->outputmode, "all") == 0
1725  || strcmp(StdI->outputmode, "full") == 0) {
1726  StdI->ioutputmode = 2;
1727  fprintf(stdout, " ioutputmode = %-10d\n", StdI->ioutputmode);
1728  }
1729  else{
1730  fprintf(stdout, "\n ERROR ! Unsupported OutPutMode : %s\n", StdI->outputmode);
1731  StdFace_exit(-1);
1732  }
1733 }/*static void CheckOutputMode*/
1739 static void CheckModPara(struct StdIntList *StdI)
1740 {
1741 
1742 
1743 #if defined(_HPhi)
1744  StdFace_PrintVal_i("Lanczos_max", &StdI->Lanczos_max, 2000);
1745  StdFace_PrintVal_i("initial_iv", &StdI->initial_iv, -1);
1746  /*StdFace_PrintVal_i("nvec", &StdI->nvec, 1);*/
1747  StdFace_PrintVal_i("exct", &StdI->exct, 1);
1748  StdFace_PrintVal_i("LanczosEps", &StdI->LanczosEps, 14);
1749  StdFace_PrintVal_i("LanczosTarget", &StdI->LanczosTarget, 2);
1750  if(StdI->LanczosTarget < StdI->exct) StdI->LanczosTarget = StdI->exct;
1751  StdFace_PrintVal_i("NumAve", &StdI->NumAve, 5);
1752  StdFace_PrintVal_i("ExpecInterval", &StdI->ExpecInterval, 20);
1753  StdFace_PrintVal_i("NOmega", &StdI->Nomega, 200);
1754  StdFace_PrintVal_d("OmegaMax", &StdI->OmegaMax, StdI->LargeValue*StdI->nsite);
1755  StdFace_PrintVal_d("OmegaMin", &StdI->OmegaMin, -StdI->LargeValue*StdI->nsite);
1756  StdFace_PrintVal_d("OmegaIm", &StdI->OmegaIm, 0.01* (int)StdI->LargeValue);
1757 #elif defined(_mVMC)
1758  if (strcmp(StdI->CParaFileHead, "****") == 0) {
1759  strcpy(StdI->CParaFileHead, "zqp\0");
1760  fprintf(stdout, " CParaFileHead = %-12s###### DEFAULT VALUE IS USED ######\n", StdI->CParaFileHead);
1761  }
1762  else fprintf(stdout, " CParaFileHead = %-s\n", StdI->CParaFileHead);
1763 
1764  StdFace_PrintVal_i("NVMCCalMode", &StdI->NVMCCalMode, 0);
1765  StdFace_PrintVal_i("NLanczosMode", &StdI->NLanczosMode, 0);
1766  StdFace_PrintVal_i("NDataIdxStart", &StdI->NDataIdxStart, 1);
1767 
1768  if (StdI->NVMCCalMode == 0) StdFace_NotUsed_i("NDataQtySmp", StdI->NDataQtySmp);
1769  /*else*/StdFace_PrintVal_i("NDataQtySmp", &StdI->NDataQtySmp, 1);
1770 
1771  if (StdI->lGC == 0 && (StdI->Sz2 == 0 || StdI->Sz2 == StdI->NaN_i)) {
1772  StdFace_PrintVal_i("NSPGaussLeg", &StdI->NSPGaussLeg, 8);
1773  StdFace_PrintVal_i("NSPStot", &StdI->NSPStot, 0);
1774  }
1775  else {
1776  StdFace_NotUsed_i("NSPGaussLeg", StdI->NSPGaussLeg);
1777  StdFace_NotUsed_i("NSPStot", StdI->NSPStot);
1778  }
1779 
1780  if (StdI->AntiPeriod[0] == 1 || StdI->AntiPeriod[1] == 1 || StdI->AntiPeriod[2] == 2)
1781  StdFace_PrintVal_i("NMPTrans", &StdI->NMPTrans, -1);
1782  else StdFace_PrintVal_i("NMPTrans", &StdI->NMPTrans, 1);
1783 
1784  StdFace_PrintVal_i("NSROptItrStep", &StdI->NSROptItrStep, 1000);
1785 
1786  if (StdI->NVMCCalMode == 1) StdFace_NotUsed_i("NSROptItrSmp", StdI->NSROptItrSmp);
1787  /*else*/ StdFace_PrintVal_i("NSROptItrSmp", &StdI->NSROptItrSmp, StdI->NSROptItrStep/10);
1788 
1789  StdFace_PrintVal_i("NVMCWarmUp", &StdI->NVMCWarmUp, 10);
1790  StdFace_PrintVal_i("NVMCInterval", &StdI->NVMCInterval, 1);
1791  StdFace_PrintVal_i("NVMCSample", &StdI->NVMCSample, 1000);
1792 
1793  if (strcmp(StdI->model, "hubbard") == 0) StdI->NExUpdatePath = 0;
1794  else if (strcmp(StdI->model, "spin") == 0) StdI->NExUpdatePath = 2;
1795  else if (strcmp(StdI->model, "kondo") == 0) {
1796  if(StdI->lGC==0) StdI->NExUpdatePath = 1;
1797  else StdI->NExUpdatePath = 3;
1798  }
1799  fprintf(stdout, " %15s = %-10d\n", "NExUpdatePath", StdI->NExUpdatePath);
1800 
1801  StdFace_PrintVal_i("RndSeed", &StdI->RndSeed, 123456789);
1802  StdFace_PrintVal_i("NSplitSize", &StdI->NSplitSize, 1);
1803  StdFace_PrintVal_i("NStore", &StdI->NStore, 1);
1804  StdFace_PrintVal_i("NSRCG", &StdI->NSRCG, 0);
1805 
1806  StdFace_PrintVal_d("DSROptRedCut", &StdI->DSROptRedCut, 0.001);
1807  StdFace_PrintVal_d("DSROptStaDel", &StdI->DSROptStaDel, 0.02);
1808  StdFace_PrintVal_d("DSROptStepDt", &StdI->DSROptStepDt, 0.02);
1809 #endif
1810  /*
1811  (Un)Conserved variables (Number of electrons, total Sz)
1812  */
1813  if (strcmp(StdI->model, "hubbard") == 0){
1814 #if defined(_HPhi)
1815  if (StdI->lGC == 0) StdFace_RequiredVal_i("nelec", StdI->nelec);
1816  else {
1817  StdFace_NotUsed_i("nelec", StdI->nelec);
1818  StdFace_NotUsed_i("2Sz", StdI->Sz2);
1819  }
1820 #else
1821  StdFace_RequiredVal_i("nelec", StdI->nelec);
1822  if (StdI->lGC == 0) StdFace_PrintVal_i("2Sz", &StdI->Sz2, 0);
1823  else StdFace_NotUsed_i("2Sz", StdI->Sz2);
1824 #endif
1825  }
1826  else if (strcmp(StdI->model, "spin") == 0) {
1827  StdFace_NotUsed_i("nelec", StdI->nelec);
1828 #if defined(_mVMC)
1829  StdI->nelec = 0;
1830 #endif
1831  if (StdI->lGC == 0) StdFace_RequiredVal_i("2Sz", StdI->Sz2);
1832  else StdFace_NotUsed_i("2Sz", StdI->Sz2);
1833  }/*else if (strcmp(StdI->model, "spin") == 0)*/
1834  else if (strcmp(StdI->model, "kondo") == 0) {
1835 #if defined(_HPhi)
1836  if (StdI->lGC == 0) StdFace_RequiredVal_i("nelec", StdI->nelec);
1837  else {
1838  StdFace_NotUsed_i("nelec", StdI->nelec);
1839  StdFace_NotUsed_i("2Sz", StdI->Sz2);
1840  }
1841 #else
1842  StdFace_RequiredVal_i("nelec", StdI->nelec);
1843  if (StdI->lGC == 0) StdFace_PrintVal_i("2Sz", &StdI->Sz2, 0);
1844  else StdFace_NotUsed_i("2Sz", StdI->Sz2);
1845 #endif
1846  }/*else if (strcmp(StdI->model, "kondo") == 0)*/
1847 }/*static void CheckModPara*/
1852 static void PrintInteractions(struct StdIntList *StdI)
1853 {
1854  FILE *fp;
1855  int nintr0, kintr, jintr;
1856  /*
1857  Coulomb INTRA
1858  */
1859  for (kintr = 0; kintr < StdI->NCintra; kintr++) {
1860  for (jintr = kintr + 1; jintr < StdI->NCintra; jintr++)
1861  if(StdI->CintraIndx[jintr][0] == StdI->CintraIndx[kintr][0])
1862  {
1863  StdI->Cintra[kintr] += StdI->Cintra[jintr];
1864  StdI->Cintra[jintr] = 0.0;
1865  }
1866  }
1867  nintr0 = 0;
1868  for (kintr = 0; kintr < StdI->NCintra; kintr++) {
1869  if (fabs(StdI->Cintra[kintr]) > 0.000001) nintr0 = nintr0 + 1;
1870  }
1871  if (nintr0 == 0 || StdI->lBoost == 1) StdI->LCintra = 0;
1872  else StdI->LCintra = 1;
1873 
1874  if (StdI->LCintra == 1) {
1875  fp = fopen("coulombintra.def", "w");
1876  fprintf(fp, "=============================================\n");
1877  fprintf(fp, "NCoulombIntra %10d\n", nintr0);
1878  fprintf(fp, "=============================================\n");
1879  fprintf(fp, "================== CoulombIntra ================\n");
1880  fprintf(fp, "=============================================\n");
1881  for (kintr = 0; kintr < StdI->NCintra; kintr++) {
1882  if (fabs(StdI->Cintra[kintr]) > 0.000001)
1883  fprintf(fp, "%5d %25.15f\n",
1884  StdI->CintraIndx[kintr][0], StdI->Cintra[kintr]);
1885  }
1886  fflush(fp);
1887  fclose(fp);
1888  fprintf(stdout, " coulombintra.def is written.\n");
1889  }/*if (StdI->LCintra == 1)*/
1890  /*
1891  Coulomb INTER
1892  */
1893  for (kintr = 0; kintr < StdI->NCinter; kintr++) {
1894  for (jintr = kintr + 1; jintr < StdI->NCinter; jintr++)
1895  if (
1896  ( StdI->CinterIndx[jintr][0] == StdI->CinterIndx[kintr][0]
1897  && StdI->CinterIndx[jintr][1] == StdI->CinterIndx[kintr][1])
1898  ||
1899  ( StdI->CinterIndx[jintr][0] == StdI->CinterIndx[kintr][1]
1900  && StdI->CinterIndx[jintr][1] == StdI->CinterIndx[kintr][0])
1901  )
1902  {
1903  StdI->Cinter[kintr] += StdI->Cinter[jintr];
1904  StdI->Cinter[jintr] = 0.0;
1905  }
1906  }/*for (kintr = 0; kintr < StdI->NCinter; kintr++)*/
1907  nintr0 = 0;
1908  for (kintr = 0; kintr < StdI->NCinter; kintr++) {
1909  if (fabs(StdI->Cinter[kintr]) > 0.000001) nintr0 = nintr0 + 1;
1910  }
1911  if (nintr0 == 0 || StdI->lBoost == 1) StdI->LCinter = 0;
1912  else StdI->LCinter = 1;
1913 
1914  if (StdI->LCinter == 1) {
1915  fp = fopen("coulombinter.def", "w");
1916  fprintf(fp, "=============================================\n");
1917  fprintf(fp, "NCoulombInter %10d\n", nintr0);
1918  fprintf(fp, "=============================================\n");
1919  fprintf(fp, "================== CoulombInter ================\n");
1920  fprintf(fp, "=============================================\n");
1921  for (kintr = 0; kintr < StdI->NCinter; kintr++) {
1922  if (fabs(StdI->Cinter[kintr]) > 0.000001)
1923  fprintf(fp, "%5d %5d %25.15f\n",
1924  StdI->CinterIndx[kintr][0], StdI->CinterIndx[kintr][1], StdI->Cinter[kintr]);
1925  }
1926  fflush(fp);
1927  fclose(fp);
1928  fprintf(stdout, " coulombinter.def is written.\n");
1929  }/*if (StdI->LCinter == 1)*/
1930  /*
1931  Hund
1932  */
1933  for (kintr = 0; kintr < StdI->NHund; kintr++) {
1934  for (jintr = kintr + 1; jintr < StdI->NHund; jintr++)
1935  if (
1936  (StdI->HundIndx[jintr][0] == StdI->HundIndx[kintr][0]
1937  && StdI->HundIndx[jintr][1] == StdI->HundIndx[kintr][1])
1938  ||
1939  (StdI->HundIndx[jintr][0] == StdI->HundIndx[kintr][1]
1940  && StdI->HundIndx[jintr][1] == StdI->HundIndx[kintr][0])
1941  )
1942  {
1943  StdI->Hund[kintr] += StdI->Hund[jintr];
1944  StdI->Hund[jintr] = 0.0;
1945  }
1946  }/*for (kintr = 0; kintr < StdI->NHund; kintr++)*/
1947  nintr0 = 0;
1948  for (kintr = 0; kintr < StdI->NHund; kintr++) {
1949  if (fabs(StdI->Hund[kintr]) > 0.000001) nintr0 = nintr0 + 1;
1950  }
1951  if (nintr0 == 0 || StdI->lBoost == 1) StdI->LHund = 0;
1952  else StdI->LHund = 1;
1953 
1954  if (StdI->LHund == 1) {
1955  fp = fopen("hund.def", "w");
1956  fprintf(fp, "=============================================\n");
1957  fprintf(fp, "NHund %10d\n", nintr0);
1958  fprintf(fp, "=============================================\n");
1959  fprintf(fp, "=============== Hund coupling ===============\n");
1960  fprintf(fp, "=============================================\n");
1961  for (kintr = 0; kintr < StdI->NHund; kintr++) {
1962  if (fabs(StdI->Hund[kintr]) > 0.000001)
1963  fprintf(fp, "%5d %5d %25.15f\n",
1964  StdI->HundIndx[kintr][0], StdI->HundIndx[kintr][1], StdI->Hund[kintr]);
1965  }
1966  fflush(fp);
1967  fclose(fp);
1968  fprintf(stdout, " hund.def is written.\n");
1969  }/*if (StdI->LHund == 1)*/
1970  /*
1971  Exchange
1972  */
1973  for (kintr = 0; kintr < StdI->NEx; kintr++) {
1974  for (jintr = kintr + 1; jintr < StdI->NEx; jintr++)
1975  if (
1976  (StdI->ExIndx[jintr][0] == StdI->ExIndx[kintr][0]
1977  && StdI->ExIndx[jintr][1] == StdI->ExIndx[kintr][1])
1978  ||
1979  (StdI->ExIndx[jintr][0] == StdI->ExIndx[kintr][1]
1980  && StdI->ExIndx[jintr][1] == StdI->ExIndx[kintr][0])
1981  )
1982  {
1983  StdI->Ex[kintr] += StdI->Ex[jintr];
1984  StdI->Ex[jintr] = 0.0;
1985  }
1986  }/*for (kintr = 0; kintr < StdI->NEx; kintr++)*/
1987  nintr0 = 0;
1988  for (kintr = 0; kintr < StdI->NEx; kintr++) {
1989  if (fabs(StdI->Ex[kintr]) > 0.000001) nintr0 = nintr0 + 1;
1990  }
1991  if (nintr0 == 0 || StdI->lBoost == 1) StdI->LEx = 0;
1992  else StdI->LEx = 1;
1993 
1994  if (StdI->LEx == 1) {
1995  fp = fopen("exchange.def", "w");
1996  fprintf(fp, "=============================================\n");
1997  fprintf(fp, "NExchange %10d\n", nintr0);
1998  fprintf(fp, "=============================================\n");
1999  fprintf(fp, "====== ExchangeCoupling coupling ============\n");
2000  fprintf(fp, "=============================================\n");
2001  for (kintr = 0; kintr < StdI->NEx; kintr++) {
2002  if (fabs(StdI->Ex[kintr]) > 0.000001)
2003  fprintf(fp, "%5d %5d %25.15f\n",
2004  StdI->ExIndx[kintr][0], StdI->ExIndx[kintr][1], StdI->Ex[kintr]);
2005  }
2006  fflush(fp);
2007  fclose(fp);
2008  fprintf(stdout, " exchange.def is written.\n");
2009  }
2010  /*
2011  PairLift
2012  */
2013  for (kintr = 0; kintr < StdI->NPairLift; kintr++) {
2014  for (jintr = kintr + 1; jintr < StdI->NPairLift; jintr++)
2015  if (
2016  (StdI->PLIndx[jintr][0] == StdI->PLIndx[kintr][0]
2017  && StdI->PLIndx[jintr][1] == StdI->PLIndx[kintr][1])
2018  ||
2019  (StdI->PLIndx[jintr][0] == StdI->PLIndx[kintr][1]
2020  && StdI->PLIndx[jintr][1] == StdI->PLIndx[kintr][0])
2021  )
2022  {
2023  StdI->PairLift[kintr] += StdI->PairLift[jintr];
2024  StdI->PairLift[jintr] = 0.0;
2025  }
2026  }/*for (kintr = 0; kintr < StdI->NPairLift; kintr++)*/
2027  nintr0 = 0;
2028  for (kintr = 0; kintr < StdI->NPairLift; kintr++) {
2029  if (fabs(StdI->PairLift[kintr]) > 0.000001) nintr0 = nintr0 + 1;
2030  }
2031  if (nintr0 == 0 || StdI->lBoost == 1) StdI->LPairLift = 0;
2032  else StdI->LPairLift = 1;
2033 
2034  if (StdI->LPairLift == 1) {
2035  fp = fopen("pairlift.def", "w");
2036  fprintf(fp, "=============================================\n");
2037  fprintf(fp, "NPairLift %10d\n", nintr0);
2038  fprintf(fp, "=============================================\n");
2039  fprintf(fp, "====== Pair-Lift term ============\n");
2040  fprintf(fp, "=============================================\n");
2041  for (kintr = 0; kintr < StdI->NPairLift; kintr++) {
2042  if (fabs(StdI->PairLift[kintr]) > 0.000001)
2043  fprintf(fp, "%5d %5d %25.15f\n",
2044  StdI->PLIndx[kintr][0], StdI->PLIndx[kintr][1], StdI->PairLift[kintr]);
2045  }
2046  fflush(fp);
2047  fclose(fp);
2048  fprintf(stdout, " pairlift.def is written.\n");
2049  }
2050  /*
2051  PairHopp
2052  */
2053  for (kintr = 0; kintr < StdI->NPairHopp; kintr++) {
2054  for (jintr = kintr + 1; jintr < StdI->NPairHopp; jintr++)
2055  if (
2056  (StdI->PHIndx[jintr][0] == StdI->PHIndx[kintr][0]
2057  && StdI->PHIndx[jintr][1] == StdI->PHIndx[kintr][1])
2058  ||
2059  (StdI->PHIndx[jintr][0] == StdI->PHIndx[kintr][1]
2060  && StdI->PHIndx[jintr][1] == StdI->PHIndx[kintr][0])
2061  )
2062  {
2063  StdI->PairHopp[kintr] += StdI->PairHopp[jintr];
2064  StdI->PairHopp[jintr] = 0.0;
2065  }
2066  }/*for (kintr = 0; kintr < StdI->NPairHopp; kintr++)*/
2067  nintr0 = 0;
2068  for (kintr = 0; kintr < StdI->NPairHopp; kintr++) {
2069  if (fabs(StdI->PairHopp[kintr]) > 0.000001) nintr0 = nintr0 + 1;
2070  }
2071  if (nintr0 == 0 || StdI->lBoost == 1) StdI->LPairHopp = 0;
2072  else StdI->LPairHopp = 1;
2073 
2074  if (StdI->LPairHopp == 1) {
2075  fp = fopen("pairhopp.def", "w");
2076  fprintf(fp, "=============================================\n");
2077  fprintf(fp, "NPairHopp %10d\n", nintr0);
2078  fprintf(fp, "=============================================\n");
2079  fprintf(fp, "====== Pair-Hopping term ============\n");
2080  fprintf(fp, "=============================================\n");
2081  for (kintr = 0; kintr < StdI->NPairHopp; kintr++) {
2082  if (fabs(StdI->PairHopp[kintr]) > 0.000001)
2083  fprintf(fp, "%5d %5d %25.15f\n",
2084  StdI->PHIndx[kintr][0], StdI->PHIndx[kintr][1], StdI->PairHopp[kintr]);
2085  }
2086  fflush(fp);
2087  fclose(fp);
2088  fprintf(stdout, " pairhopp.def is written.\n");
2089  }
2090  /*
2091  InterAll
2092  */
2093  for (jintr = 0; jintr < StdI->nintr; jintr++) {
2094  for (kintr = jintr + 1; kintr < StdI->nintr; kintr++) {
2095  if (
2096  (StdI->intrindx[jintr][0] == StdI->intrindx[kintr][0]
2097  && StdI->intrindx[jintr][1] == StdI->intrindx[kintr][1]
2098  && StdI->intrindx[jintr][2] == StdI->intrindx[kintr][2]
2099  && StdI->intrindx[jintr][3] == StdI->intrindx[kintr][3]
2100  && StdI->intrindx[jintr][4] == StdI->intrindx[kintr][4]
2101  && StdI->intrindx[jintr][5] == StdI->intrindx[kintr][5]
2102  && StdI->intrindx[jintr][6] == StdI->intrindx[kintr][6]
2103  && StdI->intrindx[jintr][7] == StdI->intrindx[kintr][7])
2104  ||
2105  (StdI->intrindx[jintr][0] == StdI->intrindx[kintr][4]
2106  && StdI->intrindx[jintr][1] == StdI->intrindx[kintr][5]
2107  && StdI->intrindx[jintr][2] == StdI->intrindx[kintr][6]
2108  && StdI->intrindx[jintr][3] == StdI->intrindx[kintr][7]
2109  && StdI->intrindx[jintr][4] == StdI->intrindx[kintr][0]
2110  && StdI->intrindx[jintr][5] == StdI->intrindx[kintr][1]
2111  && StdI->intrindx[jintr][6] == StdI->intrindx[kintr][2]
2112  && StdI->intrindx[jintr][7] == StdI->intrindx[kintr][3])
2113  ) {
2114  StdI->intr[jintr] = StdI->intr[jintr] + StdI->intr[kintr];
2115  StdI->intr[kintr] = 0.0;
2116  }
2117  else if (
2118  (StdI->intrindx[jintr][0] == StdI->intrindx[kintr][4]
2119  && StdI->intrindx[jintr][1] == StdI->intrindx[kintr][5]
2120  && StdI->intrindx[jintr][2] == StdI->intrindx[kintr][2]
2121  && StdI->intrindx[jintr][3] == StdI->intrindx[kintr][3]
2122  && StdI->intrindx[jintr][4] == StdI->intrindx[kintr][0]
2123  && StdI->intrindx[jintr][5] == StdI->intrindx[kintr][1]
2124  && StdI->intrindx[jintr][6] == StdI->intrindx[kintr][6]
2125  && StdI->intrindx[jintr][7] == StdI->intrindx[kintr][7])
2126  ||
2127  (StdI->intrindx[jintr][0] == StdI->intrindx[kintr][0]
2128  && StdI->intrindx[jintr][1] == StdI->intrindx[kintr][1]
2129  && StdI->intrindx[jintr][2] == StdI->intrindx[kintr][6]
2130  && StdI->intrindx[jintr][3] == StdI->intrindx[kintr][7]
2131  && StdI->intrindx[jintr][4] == StdI->intrindx[kintr][4]
2132  && StdI->intrindx[jintr][5] == StdI->intrindx[kintr][5]
2133  && StdI->intrindx[jintr][6] == StdI->intrindx[kintr][2]
2134  && StdI->intrindx[jintr][7] == StdI->intrindx[kintr][3])
2135  ) {
2136  StdI->intr[jintr] = StdI->intr[jintr] - StdI->intr[kintr];
2137  StdI->intr[kintr] = 0.0;
2138  }
2139  }/*for (kintr = jintr + 1; kintr < StdI->nintr; kintr++)*/
2140  }/*for (jintr = 0; jintr < StdI->nintr; jintr++)*/
2141 
2142  for (jintr = 0; jintr < StdI->nintr; jintr++) {
2143  for (kintr = jintr + 1; kintr < StdI->nintr; kintr++) {
2144  if (StdI->intrindx[jintr][6] == StdI->intrindx[kintr][4]
2145  && StdI->intrindx[jintr][7] == StdI->intrindx[kintr][5]
2146  && StdI->intrindx[jintr][4] == StdI->intrindx[kintr][6]
2147  && StdI->intrindx[jintr][5] == StdI->intrindx[kintr][7]
2148  && StdI->intrindx[jintr][2] == StdI->intrindx[kintr][0]
2149  && StdI->intrindx[jintr][3] == StdI->intrindx[kintr][1]
2150  && StdI->intrindx[jintr][0] == StdI->intrindx[kintr][2]
2151  && StdI->intrindx[jintr][1] == StdI->intrindx[kintr][3]
2152  ) {
2153  StdI->intrindx[kintr][0] = StdI->intrindx[jintr][6];
2154  StdI->intrindx[kintr][1] = StdI->intrindx[jintr][7];
2155  StdI->intrindx[kintr][2] = StdI->intrindx[jintr][4];
2156  StdI->intrindx[kintr][3] = StdI->intrindx[jintr][5];
2157  StdI->intrindx[kintr][4] = StdI->intrindx[jintr][2];
2158  StdI->intrindx[kintr][5] = StdI->intrindx[jintr][3];
2159  StdI->intrindx[kintr][6] = StdI->intrindx[jintr][0];
2160  StdI->intrindx[kintr][7] = StdI->intrindx[jintr][1];
2161  }
2162  else if (
2163  (StdI->intrindx[jintr][6] == StdI->intrindx[kintr][4]
2164  && StdI->intrindx[jintr][7] == StdI->intrindx[kintr][5]
2165  && StdI->intrindx[jintr][4] == StdI->intrindx[kintr][2]
2166  && StdI->intrindx[jintr][5] == StdI->intrindx[kintr][3]
2167  && StdI->intrindx[jintr][2] == StdI->intrindx[kintr][0]
2168  && StdI->intrindx[jintr][3] == StdI->intrindx[kintr][1]
2169  && StdI->intrindx[jintr][0] == StdI->intrindx[kintr][6]
2170  && StdI->intrindx[jintr][1] == StdI->intrindx[kintr][7])
2171  ||
2172  (StdI->intrindx[jintr][6] == StdI->intrindx[kintr][0]
2173  && StdI->intrindx[jintr][7] == StdI->intrindx[kintr][1]
2174  && StdI->intrindx[jintr][4] == StdI->intrindx[kintr][6]
2175  && StdI->intrindx[jintr][5] == StdI->intrindx[kintr][7]
2176  && StdI->intrindx[jintr][2] == StdI->intrindx[kintr][4]
2177  && StdI->intrindx[jintr][3] == StdI->intrindx[kintr][5]
2178  && StdI->intrindx[jintr][0] == StdI->intrindx[kintr][2]
2179  && StdI->intrindx[jintr][1] == StdI->intrindx[kintr][3])
2180  ) {
2181  StdI->intrindx[kintr][0] = StdI->intrindx[jintr][6];
2182  StdI->intrindx[kintr][1] = StdI->intrindx[jintr][7];
2183  StdI->intrindx[kintr][2] = StdI->intrindx[jintr][4];
2184  StdI->intrindx[kintr][3] = StdI->intrindx[jintr][5];
2185  StdI->intrindx[kintr][4] = StdI->intrindx[jintr][2];
2186  StdI->intrindx[kintr][5] = StdI->intrindx[jintr][3];
2187  StdI->intrindx[kintr][6] = StdI->intrindx[jintr][0];
2188  StdI->intrindx[kintr][7] = StdI->intrindx[jintr][1];
2189 
2190  StdI->intr[kintr] = -StdI->intr[kintr];
2191  }
2192  }/*for (kintr = jintr + 1; kintr < StdI->nintr; kintr++)*/
2193  }/*for (jintr = 0; jintr < StdI->nintr; jintr++)*/
2194 
2195  for (jintr = 0; jintr < StdI->nintr; jintr++) {
2196 
2197  if (
2198  (StdI->intrindx[jintr][0] == StdI->intrindx[jintr][4]
2199  && StdI->intrindx[jintr][1] == StdI->intrindx[jintr][5]) ||
2200  (StdI->intrindx[jintr][2] == StdI->intrindx[jintr][6]
2201  && StdI->intrindx[jintr][3] == StdI->intrindx[jintr][7])
2202  ) {
2203 
2204  if (!(
2205  (StdI->intrindx[jintr][0] == StdI->intrindx[jintr][2]
2206  && StdI->intrindx[jintr][1] == StdI->intrindx[jintr][3])
2207  ||
2208  (StdI->intrindx[jintr][0] == StdI->intrindx[jintr][6]
2209  && StdI->intrindx[jintr][1] == StdI->intrindx[jintr][7])
2210  ||
2211  (StdI->intrindx[jintr][4] == StdI->intrindx[jintr][2]
2212  && StdI->intrindx[jintr][5] == StdI->intrindx[jintr][3])
2213  ||
2214  (StdI->intrindx[jintr][4] == StdI->intrindx[jintr][6]
2215  && StdI->intrindx[jintr][5] == StdI->intrindx[jintr][7])
2216  ))
2217  StdI->intr[jintr] = 0.0;
2218  }
2219  }/*for (jintr = 0; jintr < StdI->nintr; jintr++)*/
2220 
2221  nintr0 = 0;
2222  for (kintr = 0; kintr < StdI->nintr; kintr++) {
2223  if (std::abs(StdI->intr[kintr]) > 0.000001) nintr0 = nintr0 + 1;
2224  }
2225  if (nintr0 == 0 || StdI->lBoost == 1) StdI->Lintr = 0;
2226  else StdI->Lintr = 1;
2227 
2228  if (StdI->Lintr == 1) {
2229  fp = fopen("interall.def", "w");
2230  fprintf(fp, "====================== \n");
2231  fprintf(fp, "NInterAll %7d \n", nintr0);
2232  fprintf(fp, "====================== \n");
2233  fprintf(fp, "========zInterAll===== \n");
2234  fprintf(fp, "====================== \n");
2235 
2236  if (StdI->lBoost == 0) {
2237  nintr0 = 0;
2238  for (kintr = 0; kintr < StdI->nintr; kintr++) {
2239  if (std::abs(StdI->intr[kintr]) > 0.000001)
2240  fprintf(fp, "%5d %5d %5d %5d %5d %5d %5d %5d %25.15f %25.15f\n",
2241  StdI->intrindx[kintr][0], StdI->intrindx[kintr][1],
2242  StdI->intrindx[kintr][2], StdI->intrindx[kintr][3],
2243  StdI->intrindx[kintr][4], StdI->intrindx[kintr][5],
2244  StdI->intrindx[kintr][6], StdI->intrindx[kintr][7],
2245  real(StdI->intr[kintr]), imag(StdI->intr[kintr]));
2246  }/*for (kintr = 0; kintr < StdI->nintr; kintr++)*/
2247  }/* if (StdI->lBoost == 0)*/
2248 
2249  fflush(fp);
2250  fclose(fp);
2251  fprintf(stdout, " interall.def is written.\n");
2252  }
2253 }/*static void PrintInteractions*/
2259  char *fname
2260 )
2261 {
2262  struct StdIntList *StdI;
2263  FILE *fp;
2264  int ktrans, kintr;
2265  char ctmpline[256];
2266  char *keyword, *value;
2267 
2268  StdI = (struct StdIntList *)malloc(sizeof(struct StdIntList));
2269 
2270  fprintf(stdout, "\n###### Input Parameter of Standard Intarface ######\n");
2271  if ((fp = fopen(fname, "r")) == NULL) {
2272  fprintf(stdout, "\n ERROR ! Cannot open input file %s !\n\n", fname);
2273  StdFace_exit(-1);
2274  }
2275  else {
2276  fprintf(stdout, "\n Open Standard-Mode Inputfile %s \n\n", fname);
2277  }
2278 
2279  StdFace_ResetVals(StdI);
2280 
2281  while (fgets(ctmpline, 256, fp) != NULL) {
2282 
2283  TrimSpaceQuote(ctmpline);
2284  if (strncmp(ctmpline, "//", 2) == 0) {
2285  fprintf(stdout, " Skipping a line.\n");
2286  continue;
2287  }
2288  else if (ctmpline[0] == '\0') {
2289  fprintf(stdout, " Skipping a line.\n");
2290  continue;
2291  }
2292  keyword = strtok(ctmpline, "=");
2293  value = strtok(NULL, "=");
2294  if (value == NULL) {
2295  fprintf(stdout, "\n ERROR ! \"=\" is NOT found !\n\n");
2296  StdFace_exit(-1);
2297  }
2298  Text2Lower(keyword);
2299  fprintf(stdout, " KEYWORD : %-20s | VALUE : %s \n", keyword, value);
2300 
2301  if (strcmp(keyword, "a") == 0) StoreWithCheckDup_d(keyword, value, &StdI->a);
2302  else if (strcmp(keyword, "a0h") == 0) StoreWithCheckDup_i(keyword, value, &StdI->box[0][2]);
2303  else if (strcmp(keyword, "a0l") == 0) StoreWithCheckDup_i(keyword, value, &StdI->box[0][1]);
2304  else if (strcmp(keyword, "a0w") == 0) StoreWithCheckDup_i(keyword, value, &StdI->box[0][0]);
2305  else if (strcmp(keyword, "a1h") == 0) StoreWithCheckDup_i(keyword, value, &StdI->box[1][2]);
2306  else if (strcmp(keyword, "a1l") == 0) StoreWithCheckDup_i(keyword, value, &StdI->box[1][1]);
2307  else if (strcmp(keyword, "a1w") == 0) StoreWithCheckDup_i(keyword, value, &StdI->box[1][0]);
2308  else if (strcmp(keyword, "a2h") == 0) StoreWithCheckDup_i(keyword, value, &StdI->box[2][2]);
2309  else if (strcmp(keyword, "a2l") == 0) StoreWithCheckDup_i(keyword, value, &StdI->box[2][1]);
2310  else if (strcmp(keyword, "a2w") == 0) StoreWithCheckDup_i(keyword, value, &StdI->box[2][0]);
2311  else if (strcmp(keyword, "cutoff_j") == 0) StoreWithCheckDup_d(keyword, value, &StdI->cutoff_j);
2312  else if (strcmp(keyword, "cutoff_jh") == 0) StoreWithCheckDup_i(keyword, value, &StdI->cutoff_JR[2]);
2313  else if (strcmp(keyword, "cutoff_jl") == 0) StoreWithCheckDup_i(keyword, value, &StdI->cutoff_JR[1]);
2314  else if (strcmp(keyword, "cutoff_jw") == 0) StoreWithCheckDup_i(keyword, value, &StdI->cutoff_JR[0]);
2315  else if (strcmp(keyword, "cutoff_length_j") == 0) StoreWithCheckDup_d(keyword, value, &StdI->cutoff_length_J);
2316  else if (strcmp(keyword, "cutoff_length_u") == 0) StoreWithCheckDup_d(keyword, value, &StdI->cutoff_length_U);
2317  else if (strcmp(keyword, "cutoff_length_t") == 0) StoreWithCheckDup_d(keyword, value, &StdI->cutoff_length_t);
2318  else if (strcmp(keyword, "cutoff_t") == 0) StoreWithCheckDup_d(keyword, value, &StdI->cutoff_t);
2319  else if (strcmp(keyword, "cutoff_th") == 0) StoreWithCheckDup_i(keyword, value, &StdI->cutoff_tR[2]);
2320  else if (strcmp(keyword, "cutoff_tl") == 0) StoreWithCheckDup_i(keyword, value, &StdI->cutoff_tR[1]);
2321  else if (strcmp(keyword, "cutoff_tw") == 0) StoreWithCheckDup_i(keyword, value, &StdI->cutoff_tR[0]);
2322  else if (strcmp(keyword, "cutoff_u") == 0) StoreWithCheckDup_d(keyword, value, &StdI->cutoff_u);
2323  else if (strcmp(keyword, "cutoff_uh") == 0) StoreWithCheckDup_i(keyword, value, &StdI->cutoff_UR[2]);
2324  else if (strcmp(keyword, "cutoff_ul") == 0) StoreWithCheckDup_i(keyword, value, &StdI->cutoff_UR[1]);
2325  else if (strcmp(keyword, "cutoff_uw") == 0) StoreWithCheckDup_i(keyword, value, &StdI->cutoff_UR[0]);
2326  else if (strcmp(keyword, "d") == 0) StoreWithCheckDup_d(keyword, value, &StdI->D[2][2]);
2327  else if (strcmp(keyword, "doublecounting") == 0) StoreWithCheckDup_i(keyword, value, &StdI->double_counting);
2328  else if (strcmp(keyword, "gamma") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Gamma);
2329  else if (strcmp(keyword, "h") == 0) StoreWithCheckDup_d(keyword, value, &StdI->h);
2330  else if (strcmp(keyword, "height") == 0) StoreWithCheckDup_i(keyword, value, &StdI->Height);
2331  else if (strcmp(keyword, "hlength") == 0) StoreWithCheckDup_d(keyword, value, &StdI->length[2]);
2332  else if (strcmp(keyword, "hx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->direct[2][0]);
2333  else if (strcmp(keyword, "hy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->direct[2][1]);
2334  else if (strcmp(keyword, "hz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->direct[2][2]);
2335  else if (strcmp(keyword, "j") == 0) StoreWithCheckDup_d(keyword, value, &StdI->JAll);
2336  else if (strcmp(keyword, "jx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J[0][0]);
2337  else if (strcmp(keyword, "jxy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J[0][1]);
2338  else if (strcmp(keyword, "jxz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J[0][2]);
2339  else if (strcmp(keyword, "jy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J[1][1]);
2340  else if (strcmp(keyword, "jyx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J[1][0]);
2341  else if (strcmp(keyword, "jyz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J[1][2]);
2342  else if (strcmp(keyword, "jz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J[2][2]);
2343  else if (strcmp(keyword, "jzx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J[2][0]);
2344  else if (strcmp(keyword, "jzy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J[2][1]);
2345  else if (strcmp(keyword, "j0") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0All);
2346  else if (strcmp(keyword, "j0x") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0[0][0]);
2347  else if (strcmp(keyword, "j0xy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0[0][1]);
2348  else if (strcmp(keyword, "j0xz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0[0][2]);
2349  else if (strcmp(keyword, "j0y") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0[1][1]);
2350  else if (strcmp(keyword, "j0yx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0[1][0]);
2351  else if (strcmp(keyword, "j0yz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0[1][2]);
2352  else if (strcmp(keyword, "j0z") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0[2][2]);
2353  else if (strcmp(keyword, "j0zx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0[2][0]);
2354  else if (strcmp(keyword, "j0zy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0[2][1]);
2355  else if (strcmp(keyword, "j0'") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0pAll);
2356  else if (strcmp(keyword, "j0'x") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0p[0][0]);
2357  else if (strcmp(keyword, "j0'xy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0p[0][1]);
2358  else if (strcmp(keyword, "j0'xz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0p[0][2]);
2359  else if (strcmp(keyword, "j0'y") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0p[1][1]);
2360  else if (strcmp(keyword, "j0'yx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0p[1][0]);
2361  else if (strcmp(keyword, "j0'yz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0p[1][2]);
2362  else if (strcmp(keyword, "j0'z") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0p[2][2]);
2363  else if (strcmp(keyword, "j0'zx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0p[2][0]);
2364  else if (strcmp(keyword, "j0'zy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0p[2][1]);
2365  else if (strcmp(keyword, "j0''") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0ppAll);
2366  else if (strcmp(keyword, "j0''x") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0pp[0][0]);
2367  else if (strcmp(keyword, "j0''xy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0pp[0][1]);
2368  else if (strcmp(keyword, "j0''xz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0pp[0][2]);
2369  else if (strcmp(keyword, "j0''y") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0pp[1][1]);
2370  else if (strcmp(keyword, "j0''yx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0pp[1][0]);
2371  else if (strcmp(keyword, "j0''yz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0pp[1][2]);
2372  else if (strcmp(keyword, "j0''z") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0pp[2][2]);
2373  else if (strcmp(keyword, "j0''zx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0pp[2][0]);
2374  else if (strcmp(keyword, "j0''zy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J0pp[2][1]);
2375  else if (strcmp(keyword, "j1") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1All);
2376  else if (strcmp(keyword, "j1x") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1[0][0]);
2377  else if (strcmp(keyword, "j1xy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1[0][1]);
2378  else if (strcmp(keyword, "j1xz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1[0][2]);
2379  else if (strcmp(keyword, "j1y") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1[1][1]);
2380  else if (strcmp(keyword, "j1yx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1[1][0]);
2381  else if (strcmp(keyword, "j1yz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1[1][2]);
2382  else if (strcmp(keyword, "j1z") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1[2][2]);
2383  else if (strcmp(keyword, "j1zx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1[2][0]);
2384  else if (strcmp(keyword, "j1zy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1[2][1]);
2385  else if (strcmp(keyword, "j1'") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1pAll);
2386  else if (strcmp(keyword, "j1'x") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1p[0][0]);
2387  else if (strcmp(keyword, "j1'xy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1p[0][1]);
2388  else if (strcmp(keyword, "j1'xz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1p[0][2]);
2389  else if (strcmp(keyword, "j1'y") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1p[1][1]);
2390  else if (strcmp(keyword, "j1'yx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1p[1][0]);
2391  else if (strcmp(keyword, "j1'yz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1p[1][2]);
2392  else if (strcmp(keyword, "j1'z") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1p[2][2]);
2393  else if (strcmp(keyword, "j1'zx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1p[2][0]);
2394  else if (strcmp(keyword, "j1'zy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1p[2][1]);
2395  else if (strcmp(keyword, "j1''") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1ppAll);
2396  else if (strcmp(keyword, "j1''x") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1pp[0][0]);
2397  else if (strcmp(keyword, "j1''xy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1pp[0][1]);
2398  else if (strcmp(keyword, "j1''xz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1pp[0][2]);
2399  else if (strcmp(keyword, "j1''y") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1pp[1][1]);
2400  else if (strcmp(keyword, "j1''yx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1pp[1][0]);
2401  else if (strcmp(keyword, "j1''yz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1pp[1][2]);
2402  else if (strcmp(keyword, "j1''z") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1pp[2][2]);
2403  else if (strcmp(keyword, "j1''zx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1pp[2][0]);
2404  else if (strcmp(keyword, "j1''zy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J1pp[2][1]);
2405  else if (strcmp(keyword, "j2") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2All);
2406  else if (strcmp(keyword, "j2x") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2[0][0]);
2407  else if (strcmp(keyword, "j2xy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2[0][1]);
2408  else if (strcmp(keyword, "j2xz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2[0][2]);
2409  else if (strcmp(keyword, "j2y") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2[1][1]);
2410  else if (strcmp(keyword, "j2yx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2[1][0]);
2411  else if (strcmp(keyword, "j2yz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2[1][2]);
2412  else if (strcmp(keyword, "j2z") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2[2][2]);
2413  else if (strcmp(keyword, "j2zx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2[2][0]);
2414  else if (strcmp(keyword, "j2zy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2[2][1]);
2415  else if (strcmp(keyword, "j2'") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2pAll);
2416  else if (strcmp(keyword, "j2'x") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2p[0][0]);
2417  else if (strcmp(keyword, "j2'xy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2p[0][1]);
2418  else if (strcmp(keyword, "j2'xz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2p[0][2]);
2419  else if (strcmp(keyword, "j2'y") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2p[1][1]);
2420  else if (strcmp(keyword, "j2'yx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2p[1][0]);
2421  else if (strcmp(keyword, "j2'yz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2p[1][2]);
2422  else if (strcmp(keyword, "j2'z") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2p[2][2]);
2423  else if (strcmp(keyword, "j2'zx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2p[2][0]);
2424  else if (strcmp(keyword, "j2'zy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2p[2][1]);
2425  else if (strcmp(keyword, "j2''") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2ppAll);
2426  else if (strcmp(keyword, "j2''x") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2pp[0][0]);
2427  else if (strcmp(keyword, "j2''xy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2pp[0][1]);
2428  else if (strcmp(keyword, "j2''xz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2pp[0][2]);
2429  else if (strcmp(keyword, "j2''y") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2pp[1][1]);
2430  else if (strcmp(keyword, "j2''yx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2pp[1][0]);
2431  else if (strcmp(keyword, "j2''yz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2pp[1][2]);
2432  else if (strcmp(keyword, "j2''z") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2pp[2][2]);
2433  else if (strcmp(keyword, "j2''zx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2pp[2][0]);
2434  else if (strcmp(keyword, "j2''zy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->J2pp[2][1]);
2435  else if (strcmp(keyword, "j'") == 0) StoreWithCheckDup_d(keyword, value, &StdI->JpAll);
2436  else if (strcmp(keyword, "j'x") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jp[0][0]);
2437  else if (strcmp(keyword, "j'xy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jp[0][1]);
2438  else if (strcmp(keyword, "j'xz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jp[0][2]);
2439  else if (strcmp(keyword, "j'y") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jp[1][1]);
2440  else if (strcmp(keyword, "j'yx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jp[1][0]);
2441  else if (strcmp(keyword, "j'yz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jp[1][2]);
2442  else if (strcmp(keyword, "j'z") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jp[2][2]);
2443  else if (strcmp(keyword, "j'zx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jp[2][0]);
2444  else if (strcmp(keyword, "j'zy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jp[2][1]);
2445  else if (strcmp(keyword, "j''") == 0) StoreWithCheckDup_d(keyword, value, &StdI->JppAll);
2446  else if (strcmp(keyword, "j''x") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jpp[0][0]);
2447  else if (strcmp(keyword, "j''xy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jpp[0][1]);
2448  else if (strcmp(keyword, "j''xz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jpp[0][2]);
2449  else if (strcmp(keyword, "j''y") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jpp[1][1]);
2450  else if (strcmp(keyword, "j''yx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jpp[1][0]);
2451  else if (strcmp(keyword, "j''yz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jpp[1][2]);
2452  else if (strcmp(keyword, "j''z") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jpp[2][2]);
2453  else if (strcmp(keyword, "j''zx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jpp[2][0]);
2454  else if (strcmp(keyword, "j''zy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Jpp[2][1]);
2455  else if (strcmp(keyword, "k") == 0) StoreWithCheckDup_d(keyword, value, &StdI->K);
2456  else if (strcmp(keyword, "l") == 0) StoreWithCheckDup_i(keyword, value, &StdI->L);
2457  else if (strcmp(keyword, "lattice") == 0) StoreWithCheckDup_sl(keyword, value, StdI->lattice);
2458  else if (strcmp(keyword, "llength") == 0) StoreWithCheckDup_d(keyword, value, &StdI->length[1]);
2459  else if (strcmp(keyword, "lx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->direct[1][0]);
2460  else if (strcmp(keyword, "ly") == 0) StoreWithCheckDup_d(keyword, value, &StdI->direct[1][1]);
2461  else if (strcmp(keyword, "lz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->direct[1][2]);
2462  else if (strcmp(keyword, "model") == 0) StoreWithCheckDup_sl(keyword, value, StdI->model);
2463  else if (strcmp(keyword, "mu") == 0) StoreWithCheckDup_d(keyword, value, &StdI->mu);
2464  else if (strcmp(keyword, "nelec") == 0) StoreWithCheckDup_i(keyword, value, &StdI->nelec);
2465  else if (strcmp(keyword, "outputmode") == 0) StoreWithCheckDup_sl(keyword, value, StdI->outputmode);
2466  else if (strcmp(keyword, "phase0") == 0) StoreWithCheckDup_d(keyword, value, &StdI->phase[0]);
2467  else if (strcmp(keyword, "phase1") == 0) StoreWithCheckDup_d(keyword, value, &StdI->phase[1]);
2468  else if (strcmp(keyword, "phase2") == 0) StoreWithCheckDup_d(keyword, value, &StdI->phase[2]);
2469  else if (strcmp(keyword, "t") == 0) StoreWithCheckDup_c(keyword, value, &StdI->t);
2470  else if (strcmp(keyword, "t0") == 0) StoreWithCheckDup_c(keyword, value, &StdI->t0);
2471  else if (strcmp(keyword, "t0'") == 0) StoreWithCheckDup_c(keyword, value, &StdI->t0p);
2472  else if (strcmp(keyword, "t0''") == 0) StoreWithCheckDup_c(keyword, value, &StdI->t0pp);
2473  else if (strcmp(keyword, "t1") == 0) StoreWithCheckDup_c(keyword, value, &StdI->t1);
2474  else if (strcmp(keyword, "t1'") == 0) StoreWithCheckDup_c(keyword, value, &StdI->t1p);
2475  else if (strcmp(keyword, "t1''") == 0) StoreWithCheckDup_c(keyword, value, &StdI->t1pp);
2476  else if (strcmp(keyword, "t2") == 0) StoreWithCheckDup_c(keyword, value, &StdI->t2);
2477  else if (strcmp(keyword, "t2'") == 0) StoreWithCheckDup_c(keyword, value, &StdI->t2p);
2478  else if (strcmp(keyword, "t2''") == 0) StoreWithCheckDup_c(keyword, value, &StdI->t2pp);
2479  else if (strcmp(keyword, "t'") == 0) StoreWithCheckDup_c(keyword, value, &StdI->tp);
2480  else if (strcmp(keyword, "t''") == 0) StoreWithCheckDup_c(keyword, value, &StdI->tpp);
2481  else if (strcmp(keyword, "u") == 0) StoreWithCheckDup_d(keyword, value, &StdI->U);
2482  else if (strcmp(keyword, "v") == 0) StoreWithCheckDup_d(keyword, value, &StdI->V);
2483  else if (strcmp(keyword, "v0") == 0) StoreWithCheckDup_d(keyword, value, &StdI->V0);
2484  else if (strcmp(keyword, "v0'") == 0) StoreWithCheckDup_d(keyword, value, &StdI->V0p);
2485  else if (strcmp(keyword, "v0''") == 0) StoreWithCheckDup_d(keyword, value, &StdI->V0pp);
2486  else if (strcmp(keyword, "v1") == 0) StoreWithCheckDup_d(keyword, value, &StdI->V1);
2487  else if (strcmp(keyword, "v1'") == 0) StoreWithCheckDup_d(keyword, value, &StdI->V1p);
2488  else if (strcmp(keyword, "v1''") == 0) StoreWithCheckDup_d(keyword, value, &StdI->V1pp);
2489  else if (strcmp(keyword, "v2") == 0) StoreWithCheckDup_d(keyword, value, &StdI->V2);
2490  else if (strcmp(keyword, "v2'") == 0) StoreWithCheckDup_d(keyword, value, &StdI->V2p);
2491  else if (strcmp(keyword, "v2''") == 0) StoreWithCheckDup_d(keyword, value, &StdI->V2pp);
2492  else if (strcmp(keyword, "v'") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Vp);
2493  else if (strcmp(keyword, "v''") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Vpp);
2494  else if (strcmp(keyword, "w") == 0) StoreWithCheckDup_i(keyword, value, &StdI->W);
2495  else if (strcmp(keyword, "wlength") == 0) StoreWithCheckDup_d(keyword, value, &StdI->length[0]);
2496  else if (strcmp(keyword, "wx") == 0) StoreWithCheckDup_d(keyword, value, &StdI->direct[0][0]);
2497  else if (strcmp(keyword, "wy") == 0) StoreWithCheckDup_d(keyword, value, &StdI->direct[0][1]);
2498  else if (strcmp(keyword, "wz") == 0) StoreWithCheckDup_d(keyword, value, &StdI->direct[0][2]);
2499  else if (strcmp(keyword, "2sz") == 0) StoreWithCheckDup_i(keyword, value, &StdI->Sz2);
2500 
2501 #if defined(_HPhi)
2502  else if (strcmp(keyword, "calcspec") == 0) StoreWithCheckDup_sl(keyword, value, StdI->CalcSpec);
2503  else if (strcmp(keyword, "exct") == 0) StoreWithCheckDup_i(keyword, value, &StdI->exct);
2504  else if (strcmp(keyword, "eigenvecio") == 0) StoreWithCheckDup_sl(keyword, value, StdI->EigenVecIO);
2505  else if (strcmp(keyword, "expandcoef") == 0) StoreWithCheckDup_i(keyword, value, &StdI->ExpandCoef);
2506  else if (strcmp(keyword, "expecinterval") == 0) StoreWithCheckDup_i(keyword, value, &StdI->ExpecInterval);
2507  else if (strcmp(keyword, "cdatafilehead") == 0) StoreWithCheckDup_s(keyword, value, StdI->CDataFileHead);
2508  else if (strcmp(keyword, "dt") == 0) StoreWithCheckDup_d(keyword, value, &StdI->dt);
2509  else if (strcmp(keyword, "flgtemp") == 0) StoreWithCheckDup_i(keyword, value, &StdI->FlgTemp);
2510  else if (strcmp(keyword, "freq") == 0) StoreWithCheckDup_d(keyword, value, &StdI->freq);
2511  else if (strcmp(keyword, "initialvectype") == 0) StoreWithCheckDup_sl(keyword, value, StdI->InitialVecType);
2512  else if (strcmp(keyword, "initial_iv") == 0) StoreWithCheckDup_i(keyword, value, &StdI->initial_iv);
2513  else if (strcmp(keyword, "lanczoseps") == 0) StoreWithCheckDup_i(keyword, value, &StdI->LanczosEps);
2514  else if (strcmp(keyword, "lanczostarget") == 0) StoreWithCheckDup_i(keyword, value, &StdI->LanczosTarget);
2515  else if (strcmp(keyword, "lanczos_max") == 0) StoreWithCheckDup_i(keyword, value, &StdI->Lanczos_max);
2516  else if (strcmp(keyword, "largevalue") == 0) StoreWithCheckDup_d(keyword, value, &StdI->LargeValue);
2517  else if (strcmp(keyword, "method") == 0) StoreWithCheckDup_sl(keyword, value, StdI->method);
2518  else if (strcmp(keyword, "nomega") == 0) StoreWithCheckDup_i(keyword, value, &StdI->Nomega);
2519  else if (strcmp(keyword, "numave") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NumAve);
2520  else if (strcmp(keyword, "nvec") == 0) StoreWithCheckDup_i(keyword, value, &StdI->nvec);
2521  else if (strcmp(keyword, "omegamax") == 0) StoreWithCheckDup_d(keyword, value, &StdI->OmegaMax);
2522  else if (strcmp(keyword, "omegamin") == 0) StoreWithCheckDup_d(keyword, value, &StdI->OmegaMin);
2523  else if (strcmp(keyword, "omegaim") == 0) StoreWithCheckDup_d(keyword, value, &StdI->OmegaIm);
2524  else if (strcmp(keyword, "pumptype") == 0) StoreWithCheckDup_sl(keyword, value, StdI->PumpType);
2525  else if (strcmp(keyword, "restart") == 0) StoreWithCheckDup_sl(keyword, value, StdI->Restart);
2526  else if (strcmp(keyword, "spectrumqh") == 0) StoreWithCheckDup_d(keyword, value, &StdI->SpectrumQ[2]);
2527  else if (strcmp(keyword, "spectrumql") == 0) StoreWithCheckDup_d(keyword, value, &StdI->SpectrumQ[1]);
2528  else if (strcmp(keyword, "spectrumqw") == 0) StoreWithCheckDup_d(keyword, value, &StdI->SpectrumQ[0]);
2529  else if (strcmp(keyword, "spectrumtype") == 0) StoreWithCheckDup_sl(keyword, value, StdI->SpectrumType);
2530  else if (strcmp(keyword, "tdump") == 0) StoreWithCheckDup_d(keyword, value, &StdI->tdump);
2531  else if (strcmp(keyword, "tshift") == 0) StoreWithCheckDup_d(keyword, value, &StdI->tshift);
2532  else if (strcmp(keyword, "uquench") == 0) StoreWithCheckDup_d(keyword, value, &StdI->Uquench);
2533  else if (strcmp(keyword, "vecpoth") == 0) StoreWithCheckDup_d(keyword, value, &StdI->VecPot[2]);
2534  else if (strcmp(keyword, "vecpotl") == 0) StoreWithCheckDup_d(keyword, value, &StdI->VecPot[1]);
2535  else if (strcmp(keyword, "vecpotw") == 0) StoreWithCheckDup_d(keyword, value, &StdI->VecPot[0]);
2536  else if (strcmp(keyword, "2s") == 0) StoreWithCheckDup_i(keyword, value, &StdI->S2);
2537 #elif defined(_mVMC)
2538  else if (strcmp(keyword, "a0hsub") == 0) StoreWithCheckDup_i(keyword, value, &StdI->boxsub[0][2]);
2539  else if (strcmp(keyword, "a0lsub") == 0) StoreWithCheckDup_i(keyword, value, &StdI->boxsub[0][1]);
2540  else if (strcmp(keyword, "a0wsub") == 0) StoreWithCheckDup_i(keyword, value, &StdI->boxsub[0][0]);
2541  else if (strcmp(keyword, "a1hsub") == 0) StoreWithCheckDup_i(keyword, value, &StdI->boxsub[1][2]);
2542  else if (strcmp(keyword, "a1lsub") == 0) StoreWithCheckDup_i(keyword, value, &StdI->boxsub[1][1]);
2543  else if (strcmp(keyword, "a1wsub") == 0) StoreWithCheckDup_i(keyword, value, &StdI->boxsub[1][0]);
2544  else if (strcmp(keyword, "a2hsub") == 0) StoreWithCheckDup_i(keyword, value, &StdI->boxsub[2][2]);
2545  else if (strcmp(keyword, "a2lsub") == 0) StoreWithCheckDup_i(keyword, value, &StdI->boxsub[2][1]);
2546  else if (strcmp(keyword, "a2wsub") == 0) StoreWithCheckDup_i(keyword, value, &StdI->boxsub[2][0]);
2547  else if (strcmp(keyword, "complextype") == 0) StoreWithCheckDup_i(keyword, value, &StdI->ComplexType);
2548  else if (strcmp(keyword, "cparafilehead") == 0) StoreWithCheckDup_s(keyword, value, StdI->CParaFileHead);
2549  else if (strcmp(keyword, "dsroptredcut") == 0) StoreWithCheckDup_d(keyword, value, &StdI->DSROptRedCut);
2550  else if (strcmp(keyword, "dsroptstadel") == 0) StoreWithCheckDup_d(keyword, value, &StdI->DSROptStaDel);
2551  else if (strcmp(keyword, "dsroptstepdt") == 0) StoreWithCheckDup_d(keyword, value, &StdI->DSROptStepDt);
2552  else if (strcmp(keyword, "hsub") == 0) StoreWithCheckDup_i(keyword, value, &StdI->Hsub);
2553  else if (strcmp(keyword, "lsub") == 0) StoreWithCheckDup_i(keyword, value, &StdI->Lsub);
2554  else if (strcmp(keyword, "nvmccalmode") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NVMCCalMode);
2555  else if (strcmp(keyword, "ndataidxstart") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NDataIdxStart);
2556  else if (strcmp(keyword, "ndataqtysmp") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NDataQtySmp);
2557  else if (strcmp(keyword, "nlanczosmode") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NLanczosMode);
2558  else if (strcmp(keyword, "nmptrans") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NMPTrans);
2559  else if (strcmp(keyword, "nspgaussleg") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NSPGaussLeg);
2560  else if (strcmp(keyword, "nsplitsize") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NSplitSize);
2561  else if (strcmp(keyword, "nspstot") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NSPStot);
2562  else if (strcmp(keyword, "nsroptitrsmp") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NSROptItrSmp);
2563  else if (strcmp(keyword, "nsroptitrstep") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NSROptItrStep);
2564  else if (strcmp(keyword, "nstore") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NStore);
2565  else if (strcmp(keyword, "nsrcg") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NSRCG);
2566  else if (strcmp(keyword, "nvmcinterval") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NVMCInterval);
2567  else if (strcmp(keyword, "nvmcsample") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NVMCSample);
2568  else if (strcmp(keyword, "nvmcwarmup") == 0) StoreWithCheckDup_i(keyword, value, &StdI->NVMCWarmUp);
2569  else if (strcmp(keyword, "rndseed") == 0) StoreWithCheckDup_i(keyword, value, &StdI->RndSeed);
2570  else if (strcmp(keyword, "wsub") == 0) StoreWithCheckDup_i(keyword, value, &StdI->Wsub);
2571 #endif
2572  else {
2573  fprintf(stdout, "ERROR ! Unsupported Keyword in Standard mode!\n");
2574  StdFace_exit(-1);
2575  }
2576  }
2577  fflush(fp);
2578  fclose(fp);
2579  fprintf(stdout, "\n");
2580  fprintf(stdout, "####### Construct Model #######\n");
2581  fprintf(stdout, "\n");
2582  /*
2583  Check the model
2584  */
2585  if (strcmp(StdI->CDataFileHead, "****") == 0) {
2586  strcpy(StdI->CDataFileHead, "zvo\0");
2587  fprintf(stdout, " CDataFileHead = %-12s###### DEFAULT VALUE IS USED ######\n", StdI->CDataFileHead);
2588  }
2589  else fprintf(stdout, " CDataFileHead = %-s\n", StdI->CDataFileHead);
2590 
2591  StdI->lGC = 0;
2592  StdI->lBoost = 0;
2593  if (strcmp(StdI->model, "fermionhubbard") == 0
2594  || strcmp(StdI->model, "hubbard") == 0)
2595  strcpy(StdI->model, "hubbard\0");
2596  else if(strcmp(StdI->model, "fermionhubbardgc") == 0
2597  || strcmp(StdI->model, "hubbardgc") == 0) {
2598  strcpy(StdI->model, "hubbard\0");
2599  StdI->lGC = 1;
2600  }
2601  else if (strcmp(StdI->model, "spin") == 0)
2602  strcpy(StdI->model, "spin\0");
2603  else if (strcmp(StdI->model, "spingc") == 0) {
2604  strcpy(StdI->model, "spin\0");
2605  StdI->lGC = 1;
2606  }
2607 #if defined(_HPhi)
2608  else if(strcmp(StdI->model, "spingcboost") == 0 ||
2609  strcmp(StdI->model, "spingccma") == 0) {
2610  strcpy(StdI->model, "spin\0");
2611  StdI->lGC = 1;
2612  StdI->lBoost = 1;
2613  }
2614 #endif
2615  else if (strcmp(StdI->model, "kondolattice") == 0
2616  || strcmp(StdI->model, "kondo") == 0) {
2617  strcpy(StdI->model, "kondo\0");
2618  }
2619  else if(strcmp(StdI->model, "kondolatticegc") == 0
2620  || strcmp(StdI->model, "kondogc") == 0) {
2621  strcpy(StdI->model, "kondo\0");
2622  StdI->lGC = 1;
2623  }
2624  else UnsupportedSystem(StdI->model, StdI->lattice);
2625 #if defined(_HPhi)
2626  /*
2627  Check the method
2628  */
2629  if (strcmp(StdI->method, "direct") == 0
2630  || strcmp(StdI->method, "alldiag") == 0)
2631  strcpy(StdI->method, "fulldiag\0");
2632  else if (strcmp(StdI->method, "te") == 0
2633  || strcmp(StdI->method, "time-evolution") == 0) {
2634  strcpy(StdI->method, "timeevolution\0");
2635  }
2636  /*
2637  Compute vector potential and electrical field
2638  */
2639  if (strcmp(StdI->method, "timeevolution") == 0) VectorPotential(StdI);
2640 #endif
2641  /*>>
2642  Generate Hamiltonian definition files
2643  */
2644  if (strcmp(StdI->lattice, "chain") == 0
2645  || strcmp(StdI->lattice, "chainlattice") == 0) StdFace_Chain(StdI);
2646  else if (strcmp(StdI->lattice, "face-centeredorthorhombic") == 0
2647  || strcmp(StdI->lattice, "fcorthorhombic") == 0
2648  || strcmp(StdI->lattice, "fco") == 0) StdFace_FCOrtho(StdI);
2649  else if (strcmp(StdI->lattice, "honeycomb") == 0
2650  || strcmp(StdI->lattice, "honeycomblattice") == 0) StdFace_Honeycomb(StdI);
2651  else if (strcmp(StdI->lattice, "kagome") == 0
2652  || strcmp(StdI->lattice, "kagomelattice") == 0) StdFace_Kagome(StdI);
2653  else if (strcmp(StdI->lattice, "ladder") == 0
2654  || strcmp(StdI->lattice, "ladderlattice") == 0) StdFace_Ladder(StdI);
2655  else if (strcmp(StdI->lattice, "orthorhombic") == 0
2656  || strcmp(StdI->lattice, "simpleorthorhombic") == 0) StdFace_Orthorhombic(StdI);
2657  else if (strcmp(StdI->lattice, "pyrochlore") == 0) StdFace_Pyrochlore(StdI);
2658  else if (strcmp(StdI->lattice, "tetragonal") == 0
2659  || strcmp(StdI->lattice, "tetragonallattice") == 0
2660  || strcmp(StdI->lattice, "square") == 0
2661  || strcmp(StdI->lattice, "squarelattice") == 0) StdFace_Tetragonal(StdI);
2662  else if (strcmp(StdI->lattice, "triangular") == 0
2663  || strcmp(StdI->lattice, "triangularlattice") == 0) StdFace_Triangular(StdI);
2664  else if (strcmp(StdI->lattice, "wannier90") == 0) StdFace_Wannier90(StdI);
2665  else UnsupportedSystem(StdI->model, StdI->lattice);//<<
2666 
2667 #if defined(_HPhi)
2668  StdFace_LargeValue(StdI);
2669  /*
2670  Generate Hamiltonian for Boost
2671  */
2672  if (StdI->lBoost == 1) {
2673  if (strcmp(StdI->lattice, "chain") == 0
2674  || strcmp(StdI->lattice, "chainlattice") == 0) StdFace_Chain_Boost(StdI);
2675  else if (strcmp(StdI->lattice, "honeycomb") == 0
2676  || strcmp(StdI->lattice, "honeycomblattice") == 0) StdFace_Honeycomb_Boost(StdI);
2677  else if (strcmp(StdI->lattice, "kagome") == 0
2678  || strcmp(StdI->lattice, "kagomelattice") == 0) StdFace_Kagome_Boost(StdI);
2679  else if (strcmp(StdI->lattice, "ladder") == 0
2680  || strcmp(StdI->lattice, "ladderlattice") == 0) StdFace_Ladder_Boost(StdI);
2681  else UnsupportedSystem(StdI->model, StdI->lattice);
2682  }
2683 #endif
2684 
2685  fprintf(stdout, "\n");
2686  fprintf(stdout, "###### Print Expert input files ######\n");
2687  fprintf(stdout, "\n");
2688  PrintLocSpin(StdI);
2689  PrintTrans(StdI);
2690  PrintInteractions(StdI);
2691  CheckModPara(StdI);
2692  PrintModPara(StdI);
2693 #if defined(_HPhi)
2694  PrintExcitation(StdI);
2695  if (strcmp(StdI->method, "timeevolution") == 0) PrintPump(StdI);
2696  PrintCalcMod(StdI);
2697 #elif defined(_mVMC)
2698 
2699  if(StdI->lGC == 0 && (StdI->Sz2 == 0 || StdI->Sz2 == StdI->NaN_i))
2700  StdFace_PrintVal_i("ComplexType", &StdI->ComplexType, 0);
2701  else StdFace_PrintVal_i("ComplexType", &StdI->ComplexType, 1);
2702 
2703  StdFace_generate_orb(StdI);
2704  StdFace_Proj(StdI);
2705  PrintJastrow(StdI);
2706  if(StdI->lGC == 1 || (StdI->Sz2 != 0 && StdI->Sz2 != StdI->NaN_i) )
2707  PrintOrbPara(StdI);
2708  PrintGutzwiller(StdI);
2709  PrintOrb(StdI);
2710 #endif
2711  CheckOutputMode(StdI);
2712  Print1Green(StdI);
2713  Print2Green(StdI);
2714  PrintNamelist(StdI);
2715  /*
2716  Finalize All
2717  */
2718  free(StdI->locspinflag);
2719  for (ktrans = 0; ktrans < StdI->ntrans; ktrans++) {
2720  free(StdI->transindx[ktrans]);
2721  }
2722  free(StdI->transindx);
2723  free(StdI->trans);
2724  for (kintr = 0; kintr < StdI->nintr; kintr++) {
2725  free(StdI->intrindx[kintr]);
2726  }
2727  free(StdI->intrindx);
2728  free(StdI->intr);
2729 
2730  fprintf(stdout, "\n###### Input files are generated. ######\n\n");
2731  free(StdI);
2732 }/*void StdFace_main*/
void StdFace_Ladder(struct StdIntList *StdI)
Setup a Hamiltonian for the generalized Heisenberg model on a square lattice.
Definition: Ladder.cpp:33
void StdFace_Ladder_Boost(struct StdIntList *StdI)
Definition: Ladder.cpp:291
void StdFace_Honeycomb_Boost(struct StdIntList *StdI)
void StdFace_Tetragonal(struct StdIntList *StdI)
Setup a Hamiltonian for the square lattice.
void StdFace_PrintVal_d(const char *valname, double *val, double val0)
Print a valiable (real) read from the input file if it is not specified in the input file (=NaN)...
static void TrimSpaceQuote(char *value)
Remove : space etc. from keyword and value in an iput file.
static void Print2Green(struct StdIntList *StdI)
Print greentwo.def.
void StdFace_RequiredVal_i(const char *valname, int val)
Stop HPhi if a variable (integer) which must be specified is absent in the input file (=2147483647...
void StdFace_Pyrochlore(struct StdIntList *StdI)
Setup a Hamiltonian for the Pyrochlore structure.
Definition: Pyrochlore.cpp:33
static void PrintCalcMod(struct StdIntList *StdI)
Print calcmod.def.
void StdFace_exit(int errorcode)
MPI Abortation wrapper.
void StdFace_Kagome_Boost(struct StdIntList *StdI)
Definition: Kagome.cpp:367
static void PrintNamelist(struct StdIntList *StdI)
Print namelist.def.
void StdFace_FCOrtho(struct StdIntList *StdI)
Setup a Hamiltonian for the Face-Centered Orthorhombic lattice.
Definition: FCOrtho.cpp:33
static void VectorPotential(struct StdIntList *StdI)
static void PrintLocSpin(struct StdIntList *StdI)
Print the locspin file.
void StdFace_Orthorhombic(struct StdIntList *StdI)
Setup a Hamiltonian for the Simple Orthorhombic lattice.
static void UnsupportedSystem(char *model, char *lattice)
Stop HPhi if unsupported model is read.
void StdFace_main(char *fname)
Main routine for the standard mode.
void StdFace_Triangular(struct StdIntList *StdI)
Setup a Hamiltonian for the Triangular lattice.
static void Text2Lower(char *value)
static void PrintModPara(struct StdIntList *StdI)
Print modpara.def.
static void PrintInteractions(struct StdIntList *StdI)
Output .def file for Specific interaction.
static void StoreWithCheckDup_sl(char *keyword, char *valuestring, char *value)
Store an input value into the valiable (string) Force string lower. If duplicated, HPhi will stop.
void StdFace_Chain(struct StdIntList *StdI)
Setup a Hamiltonian for the Hubbard model on a Chain lattice.
void StdFace_Kagome(struct StdIntList *StdI)
Setup a Hamiltonian for the Kagome lattice.
Definition: Kagome.cpp:33
static void StoreWithCheckDup_s(char *keyword, char *valuestring, char *value)
Store an input value into the valiable (string) If duplicated, HPhi will stop.
void StdFace_Wannier90(struct StdIntList *StdI)
Setup a Hamiltonian for the Wannier90 *_hr.dat.
Definition: Wannier90.cpp:424
static void StdFace_ResetVals(struct StdIntList *StdI)
Clear grobal variables in the standard mode All variables refered in this function is modified...
static void StoreWithCheckDup_d(char *keyword, char *valuestring, double *value)
Store an input value into the valiable (double) If duplicated, HPhi will stop.
void StdFace_NotUsed_i(const char *valname, int val)
Stop HPhi if a variable (integer) not used is specified in the input file (!=2147483647, the upper limt of Int).
static void StoreWithCheckDup_i(char *keyword, char *valuestring, int *value)
Store an input value into the valiable (integer) If duplicated, HPhi will stop.
static void Print1Green(struct StdIntList *StdI)
Print greenone.def.
void StdFace_PrintVal_i(const char *valname, int *val, int val0)
Print a valiable (integer) read from the input file if it is not specified in the input file (=214748...
static void PrintPump(struct StdIntList *StdI)
Print single.def or pair.def.
static void CheckModPara(struct StdIntList *StdI)
Summary numerical parameter check the combination of the number of sites, total spin, the number of electrons.
void StdFace_Honeycomb(struct StdIntList *StdI)
Setup a Hamiltonian for the Hubbard model on a Honeycomb lattice.
static void PrintExcitation(struct StdIntList *StdI)
Print single.def or pair.def.
static void PrintTrans(struct StdIntList *StdI)
Print the transfer file.
static void StoreWithCheckDup_c(char *keyword, char *valuestring, std::complex< double > *value)
Store an input value into the valiable (Double complex) If duplicated, HPhi will stop.
void StdFace_Chain_Boost(struct StdIntList *StdI)
Setup a Hamiltonian for the generalized Heisenberg model on a Chain lattice.
static void CheckOutputMode(struct StdIntList *StdI)
Verify outputmode.
static void StdFace_LargeValue(struct StdIntList *StdI)
Set Largevalue (StdIntList::LargeValue) for TPQ. Sum absolute-value of all one- and two- body terms...