HΦ  3.2.0
Ladder.c
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 */
21 #include "StdFace_vals.h"
22 #include "StdFace_ModelUtil.h"
23 #include <stdlib.h>
24 #include <stdio.h>
25 #include <math.h>
26 #include <complex.h>
27 #include <string.h>
28 
34  struct StdIntList *StdI
35 )
36 {
37  FILE *fp;
38  int isite, jsite, ntransMax, nintrMax;
39  int iL, isiteUC;
40  double complex Cphase;
41  double dR[3];
42 
46  fp = fopen("lattice.gp", "w");
47 
48  fprintf(stdout, " @ Lattice Size & Shape\n\n");
49 
50  StdFace_PrintVal_d("a", &StdI->a, 1.0);
51  StdFace_PrintVal_d("Wlength", &StdI->length[0], StdI->a);
52  StdFace_PrintVal_d("Llength", &StdI->length[1], StdI->a);
53  StdFace_PrintVal_d("Wx", &StdI->direct[0][0], StdI->length[0]);
54  StdFace_PrintVal_d("Wy", &StdI->direct[0][1], 0.0);
55  StdFace_PrintVal_d("Lx", &StdI->direct[1][0], 0.0);
56  StdFace_PrintVal_d("Ly", &StdI->direct[1][1], StdI->length[1]);
57 
58  StdFace_RequiredVal_i("L", StdI->L);
59  StdFace_RequiredVal_i("W", StdI->W);
60  StdFace_NotUsed_i("a0W", StdI->box[0][0]);
61  StdFace_NotUsed_i("a0L", StdI->box[0][1]);
62  StdFace_NotUsed_i("a1W", StdI->box[1][0]);
63  StdFace_NotUsed_i("a1L", StdI->box[1][1]);
64 
65  StdFace_PrintVal_d("phase0", &StdI->phase[0], 0.0);
66  StdFace_NotUsed_d("phase1", StdI->phase[1]);
67  StdI->phase[1] = StdI->phase[0];
68  StdI->phase[0] = 0.0;
69 
70  StdI->NsiteUC = StdI->W;
71  StdI->W = 1;
72  StdI->direct[0][0] = (double)StdI->NsiteUC;
73  StdFace_InitSite(StdI, fp, 2);
74  for (isite = 0; isite < StdI->NsiteUC; isite++){
75  StdI->tau[isite][0] = (double)isite / (double)StdI->NsiteUC;
76  StdI->tau[isite][1] = 0.0; StdI->tau[isite][2] = 0.0;
77  }
81  fprintf(stdout, "\n @ Hamiltonian \n\n");
82  StdFace_NotUsed_J("J", StdI->JAll, StdI->J);
83  StdFace_NotUsed_J("J'", StdI->JpAll, StdI->Jp);
84  StdFace_NotUsed_c("t", StdI->t);
85  StdFace_NotUsed_c("t'", StdI->tp);
86  StdFace_NotUsed_d("V", StdI->V);
87  StdFace_NotUsed_d("V'", StdI->Vp);
88  StdFace_NotUsed_d("K", StdI->K);
89  StdFace_PrintVal_d("h", &StdI->h, 0.0);
90  StdFace_PrintVal_d("Gamma", &StdI->Gamma, 0.0);
91 
92  if (strcmp(StdI->model, "spin") == 0 ) {
93  StdFace_PrintVal_i("2S", &StdI->S2, 1);
94  StdFace_PrintVal_d("D", &StdI->D[2][2], 0.0);
95  StdFace_InputSpin(StdI->J0, StdI->J0All, "J0");
96  StdFace_InputSpin(StdI->J1, StdI->J1All, "J1");
97  StdFace_InputSpin(StdI->J2, StdI->J2All, "J2");
98  StdFace_InputSpin(StdI->J1p, StdI->J1pAll, "J1'");
99  StdFace_InputSpin(StdI->J2p, StdI->J2pAll, "J2'");
100 
101  StdFace_NotUsed_d("mu", StdI->mu);
102  StdFace_NotUsed_d("U", StdI->U);
103  StdFace_NotUsed_c("t0", StdI->t0);
104  StdFace_NotUsed_c("t1", StdI->t1);
105  StdFace_NotUsed_c("t2", StdI->t2);
106  StdFace_NotUsed_c("t1'", StdI->t1p);
107  StdFace_NotUsed_c("t2'", StdI->t2p);
108  StdFace_NotUsed_d("V0", StdI->V0);
109  StdFace_NotUsed_d("V1", StdI->V1);
110  StdFace_NotUsed_d("V2", StdI->V2);
111  StdFace_NotUsed_d("V1'", StdI->V1p);
112  StdFace_NotUsed_d("V2'", StdI->V2p);
113  }/*if (strcmp(StdI->model, "spin") == 0 )*/
114  else {
115  StdFace_PrintVal_d("mu", &StdI->mu, 0.0);
116  StdFace_PrintVal_d("U", &StdI->U, 0.0);
117  StdFace_InputHopp(StdI->t, &StdI->t0, "t0");
118  StdFace_InputHopp(StdI->t, &StdI->t1, "t1");
119  StdFace_InputHopp(StdI->t, &StdI->t2, "t2");
120  StdFace_InputHopp(StdI->t, &StdI->t1p, "t1'");
121  StdFace_InputHopp(StdI->t, &StdI->t2p, "t2'");
122  StdFace_InputCoulombV(StdI->V, &StdI->V0, "V0");
123  StdFace_InputCoulombV(StdI->V, &StdI->V1, "V1");
124  StdFace_InputCoulombV(StdI->V, &StdI->V2, "V2");
125  StdFace_InputCoulombV(StdI->V, &StdI->V1p, "V1'");
126  StdFace_InputCoulombV(StdI->V, &StdI->V2p, "V2'");
127 
128  StdFace_NotUsed_J("J0", StdI->J0All, StdI->J0);
129  StdFace_NotUsed_J("J1", StdI->J1All, StdI->J1);
130  StdFace_NotUsed_J("J2", StdI->J2All, StdI->J2);
131  StdFace_NotUsed_J("J1p", StdI->J1pAll, StdI->J1p);
132  StdFace_NotUsed_J("J2p", StdI->J2pAll, StdI->J2p);
133  StdFace_NotUsed_d("D", StdI->D[2][2]);
134 
135  if (strcmp(StdI->model, "hubbard") == 0 ) {
136  StdFace_NotUsed_i("2S", StdI->S2);
137  StdFace_NotUsed_J("J", StdI->JAll, StdI->J);
138  }
139  else {
140  StdFace_PrintVal_i("2S", &StdI->S2, 1);
141  StdFace_InputSpin(StdI->J, StdI->JAll, "J");
142  }
143  }/*if (model != "spin")*/
144  fprintf(stdout, "\n @ Numerical conditions\n\n");
149  StdI->nsite = StdI->L * StdI->NsiteUC;
150  if (strcmp(StdI->model, "kondo") == 0 ) StdI->nsite *= 2;
151  StdI->locspinflag = (int *)malloc(sizeof(int) * StdI->nsite);
152 
153  if (strcmp(StdI->model, "spin") == 0 )
154  for (isite = 0; isite < StdI->nsite; isite++)StdI->locspinflag[isite] = StdI->S2;
155  else if (strcmp(StdI->model, "hubbard") == 0 )
156  for (isite = 0; isite < StdI->nsite; isite++)StdI->locspinflag[isite] = 0;
157  else if (strcmp(StdI->model, "kondo") == 0 )
158  for (isite = 0; isite < StdI->nsite / 2; isite++) {
159  StdI->locspinflag[isite] = StdI->S2;
160  StdI->locspinflag[isite + StdI->nsite / 2] = 0;
161  }
165  if (strcmp(StdI->model, "spin") == 0 ) {
166  ntransMax = StdI->L * StdI->NsiteUC * (StdI->S2 + 1/*h*/ + 2 * StdI->S2/*Gamma*/);
167  nintrMax = StdI->L * StdI->NsiteUC * (1/*D*/ + 1/*J1*/ + 1/*J1'*/)
168  * (3 * StdI->S2 + 1) * (3 * StdI->S2 + 1)
169  + StdI->L * (StdI->NsiteUC - 1) * (1/*J0*/ + 1/*J2*/ + 1/*J2'*/)
170  * (3 * StdI->S2 + 1) * (3 * StdI->S2 + 1);
171  }/*if (strcmp(StdI->model, "spin") == 0 )*/
172  else {
173  ntransMax = StdI->L*StdI->NsiteUC * (2/*mu+h+Gamma*/ + 2/*t1*/ + 2/*t1'*/)
174  + StdI->L*(StdI->NsiteUC - 1) * (2/*t0*/ + 2/*t2*/ + 2/*t2'*/);
175  nintrMax = StdI->L*StdI->NsiteUC * 1/*U*/
176  + StdI->L*StdI->NsiteUC * 4 * (1/*V1*/ + 1/*V1'*/)
177  + StdI->L*(StdI->NsiteUC - 1) * 4 * (1/*V0*/ + 1/*V2*/ + 1/*V2'*/);
178 
179  if (strcmp(StdI->model, "kondo") == 0) {
180  ntransMax += StdI->L * StdI->NsiteUC * (StdI->S2 + 1/*h*/ + 2 * StdI->S2/*Gamma*/);
181  nintrMax += StdI->nsite / 2 * (3 * 1 + 1) * (3 * StdI->S2 + 1);
182  }/*if (strcmp(StdI->model, "kondo") == 0)*/
183  }
184 
185  StdFace_MallocInteractions(StdI, ntransMax, nintrMax);
189  for (iL = 0; iL < StdI->L; iL++) {
190  for (isiteUC = 0; isiteUC < StdI->NsiteUC; isiteUC++) {
191 
192  isite = isiteUC + iL * StdI->NsiteUC;
193  if (strcmp(StdI->model, "kondo") == 0 ) isite += StdI->L * StdI->NsiteUC;
194  /*
195  Local term
196  */
197  if (strcmp(StdI->model, "spin") == 0 ) {
198  StdFace_MagField(StdI, StdI->S2, -StdI->h, -StdI->Gamma, isite);
199  StdFace_GeneralJ(StdI, StdI->D, StdI->S2, StdI->S2, isite, isite);
200  }/*if (strcmp(StdI->model, "spin") == 0 )*/
201  else {
202  StdFace_HubbardLocal(StdI, StdI->mu, -StdI->h, -StdI->Gamma, StdI->U, isite);
203  if (strcmp(StdI->model, "kondo") == 0 ) {
204  jsite = isiteUC + iL * StdI->NsiteUC;
205  StdFace_GeneralJ(StdI, StdI->J, 1, StdI->S2, isite, jsite);
206  StdFace_MagField(StdI, StdI->S2, -StdI->h, -StdI->Gamma, jsite);
207  }/*if (strcmp(StdI->model, "kondo") == 0 )*/
208  }/*if (model != "spin")*/
209  /*
210  Nearest neighbor along the ladder
211  */
212  StdFace_SetLabel(StdI, fp, 0, iL, 0, 1, isiteUC, isiteUC, &isite, &jsite, 1, &Cphase, dR);
213 
214  if (strcmp(StdI->model, "spin") == 0 ) {
215  StdFace_GeneralJ(StdI, StdI->J1, StdI->S2, StdI->S2, isite, jsite);
216  }/*if (strcmp(StdI->model, "spin") == 0 )*/
217  else {
218  StdFace_Hopping(StdI, Cphase * StdI->t1, isite, jsite, dR);
219  StdFace_Coulomb(StdI, StdI->V1, isite, jsite);
220  }/*if (model != "spin")*/
221  /*
222  Second nearest neighbor along the ladder
223  */
224  StdFace_SetLabel(StdI, fp, 0, iL, 0, 2, isiteUC, isiteUC, &isite, &jsite, 2, &Cphase, dR);
225 
226  if (strcmp(StdI->model, "spin") == 0 ) {
227  StdFace_GeneralJ(StdI, StdI->J1p, StdI->S2, StdI->S2, isite, jsite);
228  }/*if (strcmp(StdI->model, "spin") == 0 )*/
229  else {
230  StdFace_Hopping(StdI, Cphase * StdI->t1p, isite, jsite, dR);
231  StdFace_Coulomb(StdI, StdI->V1p, isite, jsite);
232  }/*if (model != "spin")*/
233  /*
234  Across rung
235  */
236  if (isiteUC < StdI->NsiteUC - 1) {
237  /*
238  Vertical
239  */
240  StdFace_SetLabel(StdI, fp, 0, iL, 0, 0, isiteUC, isiteUC + 1, &isite, &jsite, 1, &Cphase, dR);
241 
242  if (strcmp(StdI->model, "spin") == 0 ) {
243  StdFace_GeneralJ(StdI, StdI->J0, StdI->S2, StdI->S2, isite, jsite);
244  }/*if (strcmp(StdI->model, "spin") == 0 )*/
245  else {
246  StdFace_Hopping(StdI, Cphase * StdI->t0, isite, jsite, dR);
247  StdFace_Coulomb(StdI, StdI->V0, isite, jsite);
248  }/*if (model != "spin")*/
249  /*
250  Diagonal 1
251  */
252  StdFace_SetLabel(StdI, fp, 0, iL, 0, 1, isiteUC, isiteUC + 1, &isite, &jsite, 1, &Cphase, dR);
253 
254  if (strcmp(StdI->model, "spin") == 0 ) {
255  StdFace_GeneralJ(StdI, StdI->J2, StdI->S2, StdI->S2, isite, jsite);
256  }/*if (strcmp(StdI->model, "spin") == 0 )*/
257  else {
258  StdFace_Hopping(StdI, Cphase * StdI->t2, isite, jsite, dR);
259  StdFace_Coulomb(StdI, StdI->V2, isite, jsite);
260  }/*if (model != "spin")*/
261  /*
262  Diagonal 2
263  */
264  StdFace_SetLabel(StdI, fp, 0, iL, 0, -1, isiteUC, isiteUC + 1, &isite, &jsite, 1, &Cphase, dR);
265 
266  if (strcmp(StdI->model, "spin") == 0 ) {
267  StdFace_GeneralJ(StdI, StdI->J2p, StdI->S2, StdI->S2, isite, jsite);
268  }/*if (strcmp(StdI->model, "spin") == 0 )*/
269  else {
270  StdFace_Hopping(StdI, Cphase * StdI->t2p, isite, jsite, dR);
271  StdFace_Coulomb(StdI, StdI->V2p, isite, jsite);
272  }/*if (model != "spin")*/
273 
274  }/*if (isiteUC < StdI->NsiteUC - 1)*/
275 
276  }/*for (isiteUC = 0; isiteUC < StdI->NsiteUC; isiteUC++)*/
277  }/*for (iL = 0; iL < StdI->L; iL++)*/
278 
279  fprintf(fp, "plot \'-\' w d lc 7\n0.0 0.0\nend\npause -1\n");
280  fclose(fp);
281  StdFace_PrintGeometry(StdI);
282 }/*void StdFace_Ladder*/
283 
284 #if defined(_HPhi)
285 
292 {
293  int isite, ipivot;
294  int kintr;
295  FILE *fp;
296 
297  StdI->W = StdI->NsiteUC;
298  StdI->NsiteUC = 1;
299  /*
300  Magnetic field
301  */
302  fp = fopen("boost.def", "w");
303  fprintf(fp, "# Magnetic field\n");
304  fprintf(fp, "%25.15e %25.15e %25.15e\n",
305  -0.5 * StdI->Gamma, 0.0, -0.5 * StdI->h);
306  /*
307  Interaction
308  */
309  fprintf(fp, "%d # Number of type of J\n", 5);
310  fprintf(fp, "# J 1 (inter chain, vertical)\n");
311  fprintf(fp, "%25.15e %25.15e %25.15e\n",
312  0.25 * StdI->J0[0][0], 0.25 * StdI->J0[0][1], 0.25 * StdI->J0[0][2]);
313  fprintf(fp, "%25.15e %25.15e %25.15e\n",
314  0.25 * StdI->J0[0][1], 0.25 * StdI->J0[1][1], 0.25 * StdI->J0[1][2]);
315  fprintf(fp, "%25.15e %25.15e %25.15e\n",
316  0.25 * StdI->J0[0][2], 0.25 * StdI->J0[1][2], 0.25 * StdI->J0[2][2]);
317  fprintf(fp, "# J 2 (Nearest neighbor, along chain)\n");
318  fprintf(fp, "%25.15e %25.15e %25.15e\n",
319  0.25 * StdI->J1[0][0], 0.25 * StdI->J1[0][1], 0.25 * StdI->J1[0][2]);
320  fprintf(fp, "%25.15e %25.15e %25.15e\n",
321  0.25 * StdI->J1[0][1], 0.25 * StdI->J1[1][1], 0.25 * StdI->J1[1][2]);
322  fprintf(fp, "%25.15e %25.15e %25.15e\n",
323  0.25 * StdI->J1[0][2], 0.25 * StdI->J1[1][2], 0.25 * StdI->J1[2][2]);
324  fprintf(fp, "# J 3 (Second nearest neighbor, along chain)\n");
325  fprintf(fp, "%25.15e %25.15e %25.15e\n",
326  0.25 * StdI->J1p[0][0], 0.25 * StdI->J1p[0][1], 0.25 * StdI->J1p[0][2]);
327  fprintf(fp, "%25.15e %25.15e %25.15e\n",
328  0.25 * StdI->J1p[0][1], 0.25 * StdI->J1p[1][1], 0.25 * StdI->J1p[1][2]);
329  fprintf(fp, "%25.15e %25.15e %25.15e\n",
330  0.25 * StdI->J1p[0][2], 0.25 * StdI->J1p[1][2], 0.25 * StdI->J1p[2][2]);
331  fprintf(fp, "# J 4 (inter chain, diagonal1)\n");
332  fprintf(fp, "%25.15e %25.15e %25.15e\n",
333  0.25 * StdI->J2[0][0], 0.25 * StdI->J2[0][1], 0.25 * StdI->J2[0][2]);
334  fprintf(fp, "%25.15e %25.15e %25.15e\n",
335  0.25 * StdI->J2[0][1], 0.25 * StdI->J2[1][1], 0.25 * StdI->J2[1][2]);
336  fprintf(fp, "%25.15e %25.15e %25.15e\n",
337  0.25 * StdI->J2[0][2], 0.25 * StdI->J2[1][2], 0.25 * StdI->J2[2][2]);
338  fprintf(fp, "# J 5 (inter chain, diagonal2)\n");
339  fprintf(fp, "%25.15e %25.15e %25.15e\n",
340  0.25 * StdI->J2p[0][0], 0.25 * StdI->J2p[0][1], 0.25 * StdI->J2p[0][2]);
341  fprintf(fp, "%25.15e %25.15e %25.15e\n",
342  0.25 * StdI->J2p[0][1], 0.25 * StdI->J2p[1][1], 0.25 * StdI->J2p[1][2]);
343  fprintf(fp, "%25.15e %25.15e %25.15e\n",
344  0.25 * StdI->J2p[0][2], 0.25 * StdI->J2p[1][2], 0.25 * StdI->J2p[2][2]);
345  /*
346  Topology
347  */
348  if (StdI->S2 != 1) {
349  fprintf(stdout, "\n ERROR! S2 must be 1 in Boost. \n\n");
350  StdFace_exit(-1);
351  }
352  StdI->ishift_nspin = 2;
353  if (StdI->W != 2) {
354  fprintf(stdout, "\n ERROR! W != 2 \n\n");
355  StdFace_exit(-1);
356  }
357  if (StdI->L % 2 != 0) {
358  fprintf(stdout, "\n ERROR! L %% 2 != 0 \n\n");
359  StdFace_exit(-1);
360  }
361  if (StdI->L < 4) {
362  fprintf(stdout, "\n ERROR! L < 4 \n\n");
363  StdFace_exit(-1);
364  }
365  StdI->W = StdI->L;
366  StdI->L = 2;
367  StdI->num_pivot = StdI->W / 2;
368 
369  fprintf(fp, "# W0 R0 StdI->num_pivot StdI->ishift_nspin\n");
370  fprintf(fp, "%d %d %d %d\n", StdI->W, StdI->L, StdI->num_pivot, StdI->ishift_nspin);
371 
372  StdI->list_6spin_star = (int **)malloc(sizeof(int*) * StdI->num_pivot);
373  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
374  StdI->list_6spin_star[ipivot] = (int *)malloc(sizeof(int) * 7);
375  }
376 
377  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
378  StdI->list_6spin_star[ipivot][0] = 7; // num of J
379  StdI->list_6spin_star[ipivot][1] = 1;
380  StdI->list_6spin_star[ipivot][2] = 1;
381  StdI->list_6spin_star[ipivot][3] = 1;
382  StdI->list_6spin_star[ipivot][4] = 1;
383  StdI->list_6spin_star[ipivot][5] = 1;
384  StdI->list_6spin_star[ipivot][6] = 1; // flag
385  }
386 
387  fprintf(fp, "# StdI->list_6spin_star\n");
388  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
389  fprintf(fp, "# pivot %d\n", ipivot);
390  for (isite = 0; isite < 7; isite++) {
391  fprintf(fp, "%d ", StdI->list_6spin_star[ipivot][isite]);
392  }
393  fprintf(fp, "\n");
394  }
395 
396  StdI->list_6spin_pair = (int ***)malloc(sizeof(int**) * StdI->num_pivot);
397  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
398  StdI->list_6spin_pair[ipivot] = (int **)malloc(sizeof(int*) * 7);
399  for (isite = 0; isite < 7; isite++) {
400  StdI->list_6spin_pair[ipivot][isite] = (int *)malloc(sizeof(int) * StdI->list_6spin_star[ipivot][0]);
401  }
402  }
403 
404  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
405  StdI->list_6spin_pair[ipivot][0][0] = 0;
406  StdI->list_6spin_pair[ipivot][1][0] = 1;
407  StdI->list_6spin_pair[ipivot][2][0] = 2;
408  StdI->list_6spin_pair[ipivot][3][0] = 3;
409  StdI->list_6spin_pair[ipivot][4][0] = 4;
410  StdI->list_6spin_pair[ipivot][5][0] = 5;
411  StdI->list_6spin_pair[ipivot][6][0] = 1; // type of J
412  StdI->list_6spin_pair[ipivot][0][1] = 0;
413  StdI->list_6spin_pair[ipivot][1][1] = 2;
414  StdI->list_6spin_pair[ipivot][2][1] = 1;
415  StdI->list_6spin_pair[ipivot][3][1] = 3;
416  StdI->list_6spin_pair[ipivot][4][1] = 4;
417  StdI->list_6spin_pair[ipivot][5][1] = 5;
418  StdI->list_6spin_pair[ipivot][6][1] = 2; // type of J
419  StdI->list_6spin_pair[ipivot][0][2] = 1;
420  StdI->list_6spin_pair[ipivot][1][2] = 3;
421  StdI->list_6spin_pair[ipivot][2][2] = 0;
422  StdI->list_6spin_pair[ipivot][3][2] = 2;
423  StdI->list_6spin_pair[ipivot][4][2] = 4;
424  StdI->list_6spin_pair[ipivot][5][2] = 5;
425  StdI->list_6spin_pair[ipivot][6][2] = 2; // type of J
426  StdI->list_6spin_pair[ipivot][0][3] = 0;
427  StdI->list_6spin_pair[ipivot][1][3] = 4;
428  StdI->list_6spin_pair[ipivot][2][3] = 1;
429  StdI->list_6spin_pair[ipivot][3][3] = 2;
430  StdI->list_6spin_pair[ipivot][4][3] = 3;
431  StdI->list_6spin_pair[ipivot][5][3] = 5;
432  StdI->list_6spin_pair[ipivot][6][3] = 3; // type of J
433  StdI->list_6spin_pair[ipivot][0][4] = 1;
434  StdI->list_6spin_pair[ipivot][1][4] = 5;
435  StdI->list_6spin_pair[ipivot][2][4] = 0;
436  StdI->list_6spin_pair[ipivot][3][4] = 2;
437  StdI->list_6spin_pair[ipivot][4][4] = 3;
438  StdI->list_6spin_pair[ipivot][5][4] = 4;
439  StdI->list_6spin_pair[ipivot][6][4] = 3; // type of J
440  StdI->list_6spin_pair[ipivot][0][5] = 0;
441  StdI->list_6spin_pair[ipivot][1][5] = 3;
442  StdI->list_6spin_pair[ipivot][2][5] = 1;
443  StdI->list_6spin_pair[ipivot][3][5] = 2;
444  StdI->list_6spin_pair[ipivot][4][5] = 4;
445  StdI->list_6spin_pair[ipivot][5][5] = 5;
446  StdI->list_6spin_pair[ipivot][6][5] = 4; // type of J
447  StdI->list_6spin_pair[ipivot][0][6] = 1;
448  StdI->list_6spin_pair[ipivot][1][6] = 2;
449  StdI->list_6spin_pair[ipivot][2][6] = 0;
450  StdI->list_6spin_pair[ipivot][3][6] = 3;
451  StdI->list_6spin_pair[ipivot][4][6] = 4;
452  StdI->list_6spin_pair[ipivot][5][6] = 5;
453  StdI->list_6spin_pair[ipivot][6][6] = 5; // type of J
454  }
455 
456  fprintf(fp, "# StdI->list_6spin_pair\n");
457  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
458  fprintf(fp, "# pivot %d\n", ipivot);
459  for (kintr = 0; kintr < StdI->list_6spin_star[ipivot][0]; kintr++) {
460  for (isite = 0; isite < 7; isite++) {
461  fprintf(fp, "%d ", StdI->list_6spin_pair[ipivot][isite][kintr]);
462  }
463  fprintf(fp, "\n");
464  }
465  }
466  fclose(fp);
467 
468  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
469  free(StdI->list_6spin_star[ipivot]);
470  }
471  free(StdI->list_6spin_star);
472 
473  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
474  for (isite = 0; isite < 7; isite++) {
475  free(StdI->list_6spin_pair[ipivot][isite]);
476  }
477  free(StdI->list_6spin_pair[ipivot]);
478  }
479  free(StdI->list_6spin_pair);
480 }
481 #endif
void StdFace_PrintVal_i(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...
double V2
Anisotropic Coulomb potential (1st), input parameter.
Definition: StdFace_vals.h:83
double Jp[3][3]
Isotropic, diagonal/off-diagonal spin coupling (2nd Near.), input parameter J&#39;x, J&#39;y, J&#39;z, J&#39;xy, etc.
Definition: StdFace_vals.h:114
int box[3][3]
The shape of the super-cell. Input parameter a0W, a0L, a0H, etc. or defined from StdIntList::W, etc. in StdFace_InitSite().
Definition: StdFace_vals.h:44
double complex t2p
Anisotropic hopping (2nd), input parameter.
Definition: StdFace_vals.h:71
double J[3][3]
Isotropic, diagonal/off-diagonal spin coupling (1st Near.), input parameter Jx, Jy, Jz, Jxy, etc.
Definition: StdFace_vals.h:112
void StdFace_GeneralJ(struct StdIntList *StdI, double J[3][3], int Si2, int Sj2, int isite, int jsite)
Treat J as a 3*3 matrix [(6S + 1)*(6S&#39; + 1) interactions].
int L
Number of sites along the 2nd axis, input parameter.
Definition: StdFace_vals.h:40
void StdFace_HubbardLocal(struct StdIntList *StdI, double mu0, double h0, double Gamma0, double U0, int isite)
Add intra-Coulomb, magnetic field, chemical potential for the itenerant electron. ...
double D[3][3]
Coefficient for input parameter D. Only D[2][2] is used.
Definition: StdFace_vals.h:145
double J1p[3][3]
Isotropic, diagonal/off-diagonal spin coupling (2nd Near.), input parameter J1&#39;x, J1&#39;y...
Definition: StdFace_vals.h:128
void StdFace_Ladder_Boost(struct StdIntList *StdI)
Definition: Ladder.c:291
void StdFace_PrintGeometry(struct StdIntList *StdI)
Print geometry of sites for the pos-process of correlation function.
double JpAll
Isotropic, diagonal spin coupling (2nd Near), input parameter Jp.
Definition: StdFace_vals.h:90
void StdFace_InputHopp(double complex t, double complex *t0, char *t0name)
Input hopping integral from the input file, if it is not specified, use the default value(0 or the is...
double J1[3][3]
Isotropic, diagonal/off-diagonal spin coupling (1st Near.), input parameter J1x, J1y, J1z, J1xy, etc. or set in StdFace_InputSpinNN().
Definition: StdFace_vals.h:125
double J2p[3][3]
Isotropic, diagonal/off-diagonal spin coupling (2nd Near.), input parameter J2&#39;x, J2&#39;y...
Definition: StdFace_vals.h:137
int ** list_6spin_star
Definition: StdFace_vals.h:279
void StdFace_Hopping(struct StdIntList *StdI, double complex trans0, int isite, int jsite, double *dR)
Add Hopping for the both spin.
double complex t
Nearest-neighbor hopping, input parameter.
Definition: StdFace_vals.h:62
void StdFace_MallocInteractions(struct StdIntList *StdI, int ntransMax, int nintrMax)
Malloc Arrays for interactions.
double JAll
Isotropic, diagonal spin coupling (1st Near.), input parameter J.
Definition: StdFace_vals.h:88
int S2
Total spin |S| of a local spin, input from file.
Definition: StdFace_vals.h:236
int NsiteUC
Number of sites in the unit cell. Defined in the beginning of each lattice function.
Definition: StdFace_vals.h:53
void StdFace_InputSpin(double Jp[3][3], double JpAll, char *Jpname)
Input spin-spin interaction other than nearest-neighbor.
double J1All
Anisotropic, diagonal spin coupling (1st Near), input parameter J1.
Definition: StdFace_vals.h:98
void StdFace_InitSite(struct StdIntList *StdI, FILE *fp, int dim)
Initialize the super-cell where simulation is performed.
int W
Number of sites along the 1st axis, input parameter.
Definition: StdFace_vals.h:39
double V2p
Anisotropic Coulomb potential (2nd), input parameter.
Definition: StdFace_vals.h:84
char model[256]
Name of model, input parameter.
Definition: StdFace_vals.h:60
void StdFace_NotUsed_J(char *valname, double JAll, double J[3][3])
Stop HPhi if variables (real) not used is specified in the input file (!=NaN).
double Gamma
Transvars magnetic field, input parameter.
Definition: StdFace_vals.h:148
double V1
Anisotropic Coulomb potential (1st), input parameter.
Definition: StdFace_vals.h:80
double J0All
Anisotropic, diagonal spin coupling (1st Near), input parameter J0.
Definition: StdFace_vals.h:92
double J0[3][3]
Isotropic, diagonal/off-diagonal spin coupling (1st Near.), input parameter J0x, J0y, J0z, J0xy, etc. or set in StdFace_InputSpinNN().
Definition: StdFace_vals.h:116
double U
On-site Coulomb potential, input parameter.
Definition: StdFace_vals.h:74
double phase[3]
Boundary phase, input parameter phase0, etc.
Definition: StdFace_vals.h:154
int ishift_nspin
Definition: StdFace_vals.h:281
double length[3]
Anisotropic lattice constant, input parameter wlength, llength, hlength.
Definition: StdFace_vals.h:37
int * locspinflag
[StdIntList::nsite] LocSpin in Expert mode, malloc and set in each lattice file.
Definition: StdFace_vals.h:162
double complex tp
2nd-nearest hopping, input parameter
Definition: StdFace_vals.h:63
double complex t1
Anisotropic hopping (1st), input parameter.
Definition: StdFace_vals.h:67
double V
Off-site Coulomb potential (1st), input parameter.
Definition: StdFace_vals.h:75
double complex t0
Anisotropic hopping (1st), input parameter.
Definition: StdFace_vals.h:64
double direct[3][3]
The unit direct lattice vector. Set in StdFace_InitSite().
Definition: StdFace_vals.h:42
void StdFace_MagField(struct StdIntList *StdI, int S2, double h, double Gamma, int isite)
Add longitudinal and transvars magnetic field to the list.
void StdFace_NotUsed_d(char *valname, double val)
Stop HPhi if a variable (real) not used is specified in the input file (!=NaN).
double V0
Anisotropic Coulomb potential (1st), input parameter.
Definition: StdFace_vals.h:77
void StdFace_Coulomb(struct StdIntList *StdI, double V, int isite, int jsite)
Add onsite/offsite Coulomb term to the list StdIntList::Cinter and StdIntList::CinterIndx, and increase the number of them (StdIntList::NCinter).
void StdFace_InputCoulombV(double V, double *V0, char *V0name)
Input off-site Coulomb interaction from the input file, if it is not specified, use the default value...
void StdFace_SetLabel(struct StdIntList *StdI, FILE *fp, int iW, int iL, int diW, int diL, int isiteUC, int jsiteUC, int *isite, int *jsite, int connect, double complex *Cphase, double *dR)
Set Label in the gnuplot display (Only used in 2D system)
double complex t1p
Anisotropic hopping (2nd), input parameter.
Definition: StdFace_vals.h:68
double complex t2
Anisotropic hopping (1st), input parameter.
Definition: StdFace_vals.h:70
double mu
Chemical potential, input parameter.
Definition: StdFace_vals.h:61
void StdFace_PrintVal_d(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)...
double Vp
Off-site Coulomb potential (2nd), input parameter.
Definition: StdFace_vals.h:76
double J1pAll
Anisotropic, diagonal spin coupling (2nd Near), input parameter J1&#39;.
Definition: StdFace_vals.h:100
void StdFace_RequiredVal_i(char *valname, int val)
Stop HPhi if a variable (integer) which must be specified is absent in the input file (=2147483647...
int *** list_6spin_pair
Definition: StdFace_vals.h:278
void StdFace_NotUsed_c(char *valname, double complex val)
Stop HPhi if a variable (complex) not used is specified in the input file (!=NaN).
int nsite
Number of sites, set in the each lattice file.
Definition: StdFace_vals.h:161
double J2All
Anisotropic, diagonal spin coupling (1st Near), input parameter J2.
Definition: StdFace_vals.h:104
Variables used in the Standard mode. These variables are passed as a pointer of the structure(StdIntL...
void StdFace_NotUsed_i(char *valname, int val)
Stop HPhi if a variable (integer) not used is specified in the input file (!=2147483647, the upper limt of Int).
double h
Longitudinal magnetic field, input parameter.
Definition: StdFace_vals.h:147
double a
The lattice constant. Input parameter.
Definition: StdFace_vals.h:36
double J2[3][3]
Isotropic, diagonal/off-diagonal spin coupling (1st Near.), input parameter J2x, J2y, J2z, J2xy, etc. or set in StdFace_InputSpinNN().
Definition: StdFace_vals.h:134
void StdFace_Ladder(struct StdIntList *StdI)
Setup a Hamiltonian for the generalized Heisenberg model on a square lattice.
Definition: Ladder.c:33
double J2pAll
Anisotropic, diagonal spin coupling (2nd Near), input parameter J2&#39;.
Definition: StdFace_vals.h:106
double ** tau
Cell-internal site position in the fractional coordinate. Defined in the beginning of each lattice fu...
Definition: StdFace_vals.h:55
void StdFace_exit(int errorcode)
MPI Abortation wrapper.
double V1p
Anisotropic Coulomb potential (2nd), input parameter.
Definition: StdFace_vals.h:81
double K
4-spin term. Not used.
Definition: StdFace_vals.h:149