HΦ  3.2.0
SquareLattice.c File Reference

Standard mode for the tetragonal lattice. More...

#include "StdFace_vals.h"
#include "StdFace_ModelUtil.h"
#include <stdlib.h>
#include <stdio.h>
#include <math.h>
#include <complex.h>
#include <string.h>
+ Include dependency graph for SquareLattice.c:

Go to the source code of this file.

Functions

void StdFace_Tetragonal (struct StdIntList *StdI)
 Setup a Hamiltonian for the square lattice. More...
 

Detailed Description

Standard mode for the tetragonal lattice.

Definition in file SquareLattice.c.

Function Documentation

◆ StdFace_Tetragonal()

void StdFace_Tetragonal ( struct StdIntList StdI)

Setup a Hamiltonian for the square lattice.

Author
Mitsuaki Kawamura (The University of Tokyo)

(1) Compute the shape of the super-cell and sites in the super-cell

(2) check & store parameters of Hamiltonian

(3) Set local spin flag (StdIntList::locspinflag) and the number of sites (StdIntList::nsite)

(4) Compute the upper limit of the number of Transfer & Interaction and malloc them.

(5) Set Transfer & Interaction

Definition at line 33 of file SquareLattice.c.

References StdIntList::a, StdIntList::Cell, StdIntList::D, StdIntList::direct, StdIntList::Gamma, StdIntList::h, StdIntList::J, StdIntList::J0, StdIntList::J0All, StdIntList::J0p, StdIntList::J0pAll, StdIntList::J0pp, StdIntList::J0ppAll, StdIntList::J1, StdIntList::J1All, StdIntList::J1p, StdIntList::J1pAll, StdIntList::J1pp, StdIntList::J1ppAll, StdIntList::J2, StdIntList::J2All, StdIntList::J2p, StdIntList::J2pAll, StdIntList::JAll, StdIntList::Jp, StdIntList::JpAll, StdIntList::Jpp, StdIntList::JppAll, StdIntList::K, StdIntList::length, StdIntList::locspinflag, StdIntList::model, StdIntList::mu, StdIntList::NCell, StdIntList::nsite, StdIntList::NsiteUC, StdIntList::phase, StdIntList::S2, StdFace_Coulomb(), StdFace_GeneralJ(), StdFace_Hopping(), StdFace_HubbardLocal(), StdFace_InitSite(), StdFace_InputCoulombV(), StdFace_InputHopp(), StdFace_InputSpin(), StdFace_InputSpinNN(), StdFace_MagField(), StdFace_MallocInteractions(), StdFace_NotUsed_c(), StdFace_NotUsed_d(), StdFace_NotUsed_i(), StdFace_NotUsed_J(), StdFace_PrintGeometry(), StdFace_PrintVal_d(), StdFace_PrintVal_i(), StdFace_SetLabel(), StdIntList::t, StdIntList::t0, StdIntList::t0p, StdIntList::t0pp, StdIntList::t1, StdIntList::t1p, StdIntList::t1pp, StdIntList::t2, StdIntList::t2p, StdIntList::tau, StdIntList::tp, StdIntList::tpp, StdIntList::U, StdIntList::V, StdIntList::V0, StdIntList::V0p, StdIntList::V0pp, StdIntList::V1, StdIntList::V1p, StdIntList::V1pp, StdIntList::V2, StdIntList::V2p, StdIntList::Vp, and StdIntList::Vpp.

Referenced by StdFace_main().

34 {
35  int isite, jsite, ntransMax, nintrMax;
36  int iL, iW, kCell;
37  FILE *fp;
38  double complex Cphase;
39  double dR[3];
40 
44  fp = fopen("lattice.gp", "w");
45 
46  StdI->NsiteUC = 1;
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_PrintVal_d("phase0", &StdI->phase[0], 0.0);
59  StdFace_PrintVal_d("phase1", &StdI->phase[1], 0.0);
60 
61  StdFace_InitSite(StdI, fp, 2);
62  StdI->tau[0][0] = 0.0; StdI->tau[0][1] = 0.0; StdI->tau[0][2] = 0.0;
66  fprintf(stdout, "\n @ Hamiltonian \n\n");
67  StdFace_NotUsed_J("J2", StdI->J2All, StdI->J2);
68  StdFace_NotUsed_J("J2'", StdI->J2pAll, StdI->J2p);
69  StdFace_NotUsed_c("t2", StdI->t2);
70  StdFace_NotUsed_d("t2'", StdI->t2p);
71  StdFace_NotUsed_d("V2", StdI->V2);
72  StdFace_NotUsed_d("V2'", StdI->V2p);
73  StdFace_NotUsed_d("K", StdI->K);
74  StdFace_PrintVal_d("h", &StdI->h, 0.0);
75  StdFace_PrintVal_d("Gamma", &StdI->Gamma, 0.0);
76 
77  if (strcmp(StdI->model, "spin") == 0 ) {
78  StdFace_PrintVal_i("2S", &StdI->S2, 1);
79  StdFace_PrintVal_d("D", &StdI->D[2][2], 0.0);
80  StdFace_InputSpinNN(StdI->J, StdI->JAll, StdI->J0, StdI->J0All, "J0");
81  StdFace_InputSpinNN(StdI->J, StdI->JAll, StdI->J1, StdI->J1All, "J1");
82  StdFace_InputSpinNN(StdI->Jp, StdI->JpAll, StdI->J0p, StdI->J0pAll, "J0'");
83  StdFace_InputSpinNN(StdI->Jp, StdI->JpAll, StdI->J1p, StdI->J1pAll, "J1'");
84  StdFace_InputSpinNN(StdI->Jpp, StdI->JppAll, StdI->J0pp, StdI->J0ppAll, "J0''");
85  StdFace_InputSpinNN(StdI->Jpp, StdI->JppAll, StdI->J1pp, StdI->J1ppAll, "J1''");
86 
87  StdFace_NotUsed_d("mu", StdI->mu);
88  StdFace_NotUsed_d("U", StdI->U);
89  StdFace_NotUsed_c("t", StdI->t);
90  StdFace_NotUsed_c("t0", StdI->t0);
91  StdFace_NotUsed_c("t1", StdI->t1);
92  StdFace_NotUsed_c("t'", StdI->tp);
93  StdFace_NotUsed_c("t0'", StdI->t0p);
94  StdFace_NotUsed_c("t1'", StdI->t1p);
95  StdFace_NotUsed_c("t''", StdI->tpp);
96  StdFace_NotUsed_c("t0''", StdI->t0pp);
97  StdFace_NotUsed_c("t1''", StdI->t1pp);
98  StdFace_NotUsed_d("V", StdI->V);
99  StdFace_NotUsed_d("V0", StdI->V0);
100  StdFace_NotUsed_d("V1", StdI->V1);
101  StdFace_NotUsed_d("V'", StdI->Vp);
102  StdFace_NotUsed_d("V0'", StdI->V0p);
103  StdFace_NotUsed_d("V1'", StdI->V1p);
104  StdFace_NotUsed_d("V''", StdI->Vpp);
105  StdFace_NotUsed_d("V0''", StdI->V0pp);
106  StdFace_NotUsed_d("V1''", StdI->V1pp);
107  }/*if (strcmp(StdI->model, "spin") == 0 )*/
108  else {
109  StdFace_PrintVal_d("mu", &StdI->mu, 0.0);
110  StdFace_PrintVal_d("U", &StdI->U, 0.0);
111  StdFace_InputHopp(StdI->t, &StdI->t0, "t0");
112  StdFace_InputHopp(StdI->t, &StdI->t1, "t1");
113  StdFace_InputHopp(StdI->tp, &StdI->t0p, "t0'");
114  StdFace_InputHopp(StdI->tp, &StdI->t1p, "t1'");
115  StdFace_InputHopp(StdI->tpp, &StdI->t0pp, "t0''");
116  StdFace_InputHopp(StdI->tpp, &StdI->t1pp, "t1''");
117  StdFace_InputCoulombV(StdI->V, &StdI->V0, "V0");
118  StdFace_InputCoulombV(StdI->V, &StdI->V1, "V1");
119  StdFace_InputCoulombV(StdI->Vp, &StdI->V0p, "V0'");
120  StdFace_InputCoulombV(StdI->Vp, &StdI->V1p, "V1'");
121  StdFace_InputCoulombV(StdI->Vpp, &StdI->V0pp, "V0''");
122  StdFace_InputCoulombV(StdI->Vpp, &StdI->V1pp, "V1''");
123  StdFace_PrintVal_d("V'", &StdI->Vp, 0.0);
124 
125  StdFace_NotUsed_J("J0", StdI->J0All, StdI->J0);
126  StdFace_NotUsed_J("J1", StdI->J1All, StdI->J1);
127  StdFace_NotUsed_J("J'", StdI->JpAll, StdI->Jp);
128  StdFace_NotUsed_d("D", StdI->D[2][2]);
129 
130  if (strcmp(StdI->model, "hubbard") == 0 ) {
131  StdFace_NotUsed_i("2S", StdI->S2);
132  StdFace_NotUsed_J("J", StdI->JAll, StdI->J);
133  }/*if (strcmp(StdI->model, "hubbard") == 0 )*/
134  else {
135  StdFace_PrintVal_i("2S", &StdI->S2, 1);
136  StdFace_InputSpin(StdI->J, StdI->JAll, "J");
137  }/*if (model != "hubbard")*/
138 
139  }/*if (model != "spin")*/
140  fprintf(stdout, "\n @ Numerical conditions\n\n");
145  StdI->nsite = StdI->NsiteUC * StdI->NCell;
146  if (strcmp(StdI->model, "kondo") == 0 ) StdI->nsite *= 2;
147  StdI->locspinflag = (int *)malloc(sizeof(int) * StdI->nsite);
148 
149  if(strcmp(StdI->model, "spin") == 0 )
150  for (isite = 0; isite < StdI->nsite; isite++) StdI->locspinflag[isite] = StdI->S2;
151  else if(strcmp(StdI->model, "hubbard") == 0 )
152  for (isite = 0; isite < StdI->nsite; isite++) StdI->locspinflag[isite] = 0;
153  else
154  for (iL = 0; iL < StdI->nsite / 2; iL++) {
155  StdI->locspinflag[iL] = StdI->S2;
156  StdI->locspinflag[iL + StdI->nsite / 2] = 0;
157  }
161  if (strcmp(StdI->model, "spin") == 0 ) {
162  ntransMax = StdI->nsite * (StdI->S2 + 1/*h*/ + 2 * StdI->S2/*Gamma*/);
163  nintrMax = StdI->NCell * (StdI->NsiteUC/*D*/ + 2/*J*/ + 2/*J'*/ + 2/*J''*/)
164  * (3 * StdI->S2 + 1) * (3 * StdI->S2 + 1);
165  }
166  else {
167  ntransMax = StdI->NCell * 2/*spin*/ * (
168  2 * StdI->NsiteUC/*mu+h+Gamma*/ + 4/*t*/ + 4/*t'*/ + 4/*t''*/);
169  nintrMax = StdI->NCell * (StdI->NsiteUC/*U*/ + 4 * (2/*V*/ + 2/*V'*/ + 2/*V''*/));
170 
171  if (strcmp(StdI->model, "kondo") == 0) {
172  ntransMax += StdI->nsite / 2 * (StdI->S2 + 1/*h*/ + 2 * StdI->S2/*Gamma*/);
173  nintrMax += StdI->nsite / 2 * (3 * StdI->S2 + 1) * (3 * StdI->S2 + 1);
174  }/*if (strcmp(StdI->model, "kondo") == 0)*/
175  }
176 
177  StdFace_MallocInteractions(StdI, ntransMax, nintrMax);
181  for (kCell = 0; kCell < StdI->NCell; kCell++){
182 
183  iW = StdI->Cell[kCell][0];
184  iL = StdI->Cell[kCell][1];
185  /*
186  Local term
187  */
188  isite = kCell;
189  if (strcmp(StdI->model, "kondo") == 0 ) isite += StdI->NCell;
190 
191  if (strcmp(StdI->model, "spin") == 0 ) {
192  StdFace_MagField(StdI, StdI->S2, -StdI->h, -StdI->Gamma, isite);
193  StdFace_GeneralJ(StdI, StdI->D, StdI->S2, StdI->S2, isite, isite);
194  }/*if (strcmp(StdI->model, "spin") == 0 )*/
195  else {
196  StdFace_HubbardLocal(StdI, StdI->mu, -StdI->h, -StdI->Gamma, StdI->U, isite);
197  if (strcmp(StdI->model, "kondo") == 0 ) {
198  jsite = kCell;
199  StdFace_GeneralJ(StdI, StdI->J, 1, StdI->S2, isite, jsite);
200  StdFace_MagField(StdI, StdI->S2, -StdI->h, -StdI->Gamma, jsite);
201  }/*if (strcmp(StdI->model, "kondo") == 0 )*/
202  }
203  /*
204  Nearest neighbor along W
205  */
206  StdFace_SetLabel(StdI, fp, iW, iL, 1, 0, 0, 0, &isite, &jsite, 1, &Cphase, dR);
207 
208  if (strcmp(StdI->model, "spin") == 0 ) {
209  StdFace_GeneralJ(StdI, StdI->J0, StdI->S2, StdI->S2, isite, jsite);
210  }/*if (strcmp(StdI->model, "spin") == 0 )*/
211  else {
212  StdFace_Hopping(StdI, Cphase * StdI->t0, isite, jsite, dR);
213  StdFace_Coulomb(StdI, StdI->V0, isite, jsite);
214  }
215  /*
216  Nearest neighbor along L
217  */
218  StdFace_SetLabel(StdI, fp, iW, iL, 0, 1, 0, 0, &isite, &jsite, 1, &Cphase, dR);
219 
220  if (strcmp(StdI->model, "spin") == 0 ) {
221  StdFace_GeneralJ(StdI, StdI->J1, StdI->S2, StdI->S2, isite, jsite);
222  }
223  else {
224  StdFace_Hopping(StdI, Cphase * StdI->t1, isite, jsite, dR);
225  StdFace_Coulomb(StdI, StdI->V1, isite, jsite);
226  }
227  /*
228  Second nearest neighbor W+L
229  */
230  StdFace_SetLabel(StdI, fp, iW, iL, 1, 1, 0, 0, &isite, &jsite, 2, &Cphase, dR);
231 
232  if (strcmp(StdI->model, "spin") == 0 ) {
233  StdFace_GeneralJ(StdI, StdI->J0p, StdI->S2, StdI->S2, isite, jsite);
234  }/*if (strcmp(StdI->model, "spin") == 0 )*/
235  else {
236  StdFace_Hopping(StdI, Cphase * StdI->t0p, isite, jsite, dR);
237  StdFace_Coulomb(StdI, StdI->V0p, isite, jsite);
238  }
239  /*
240  Second nearest neighbor W-L
241  */
242  StdFace_SetLabel(StdI, fp, iW, iL, 1, -1, 0, 0, &isite, &jsite, 2, &Cphase, dR);
243 
244  if (strcmp(StdI->model, "spin") == 0 ) {
245  StdFace_GeneralJ(StdI, StdI->J1p, StdI->S2, StdI->S2, isite, jsite);
246  }/*if (strcmp(StdI->model, "spin") == 0 )*/
247  else {
248  StdFace_Hopping(StdI, Cphase * StdI->t1p, isite, jsite, dR);
249  StdFace_Coulomb(StdI, StdI->V1p, isite, jsite);
250  }/*if (model != "spin")*/
251  /*
252  Nearest neighbor along 2W
253  */
254  StdFace_SetLabel(StdI, fp, iW, iL, 2, 0, 0, 0, &isite, &jsite, 3, &Cphase, dR);
255 
256  if (strcmp(StdI->model, "spin") == 0) {
257  StdFace_GeneralJ(StdI, StdI->J0pp, StdI->S2, StdI->S2, isite, jsite);
258  }/*if (strcmp(StdI->model, "spin") == 0 )*/
259  else {
260  StdFace_Hopping(StdI, Cphase * StdI->t0pp, isite, jsite, dR);
261  StdFace_Coulomb(StdI, StdI->V0pp, isite, jsite);
262  }
263  /*
264  Nearest neighbor along L
265  */
266  StdFace_SetLabel(StdI, fp, iW, iL, 0, 1, 0, 0, &isite, &jsite, 3, &Cphase, dR);
267 
268  if (strcmp(StdI->model, "spin") == 0) {
269  StdFace_GeneralJ(StdI, StdI->J1pp, StdI->S2, StdI->S2, isite, jsite);
270  }
271  else {
272  StdFace_Hopping(StdI, Cphase * StdI->t1pp, isite, jsite, dR);
273  StdFace_Coulomb(StdI, StdI->V1pp, isite, jsite);
274  }
275  }/*for (kCell = 0; kCell < StdI->NCell; kCell++)*/
276 
277  fprintf(fp, "plot \'-\' w d lc 7\n0.0 0.0\nend\npause -1\n");
278  fclose(fp);
279  StdFace_PrintGeometry(StdI);
280 }
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 J0pp[3][3]
Isotropic, diagonal/off-diagonal spin coupling (3rd Near.), input parameter J0&#39;&#39;x, J0&#39;&#39;y, J0&#39;&#39;z, J0&#39;&#39;xy, etc. or set in StdFace_InputSpin().
Definition: StdFace_vals.h:122
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
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].
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_PrintGeometry(struct StdIntList *StdI)
Print geometry of sites for the pos-process of correlation function.
double J1ppAll
Anisotropic, diagonal spin coupling (3rd Near), input parameter J1&#39;&#39;.
Definition: StdFace_vals.h:102
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
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
double V0p
Anisotropic Coulomb potential (2nd), input parameter.
Definition: StdFace_vals.h:78
double J1pp[3][3]
Isotropic, diagonal/off-diagonal spin coupling (3rd Near.), input parameter J1&#39;&#39;x, J1&#39;&#39;y, J1&#39;&#39;z, J1&#39;&#39;xy, etc. or set in StdFace_InputSpin().
Definition: StdFace_vals.h:131
void StdFace_InitSite(struct StdIntList *StdI, FILE *fp, int dim)
Initialize the super-cell where simulation is performed.
double V2p
Anisotropic Coulomb potential (2nd), input parameter.
Definition: StdFace_vals.h:84
double complex tpp
3rd-nearest hopping, input parameter
Definition: StdFace_vals.h:73
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 JppAll
Isotropic, diagonal spin coupling (3rd Near), input parameter J&#39;&#39;.
Definition: StdFace_vals.h:110
double complex t1pp
Anisotropic hopping (3rd), input parameter.
Definition: StdFace_vals.h:69
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 J0ppAll
Anisotropic, diagonal spin coupling (3rd Near), input parameter J0&#39;&#39;.
Definition: StdFace_vals.h:96
double U
On-site Coulomb potential, input parameter.
Definition: StdFace_vals.h:74
double complex t0p
Anisotropic hopping (2nd), input parameter.
Definition: StdFace_vals.h:65
int ** Cell
[StdIntList][3] The cell position in the fractional coordinate. Malloc and Set in StdFace_InitSite()...
Definition: StdFace_vals.h:51
double phase[3]
Boundary phase, input parameter phase0, etc.
Definition: StdFace_vals.h:154
double V1pp
Anisotropic Coulomb potential (3rd), input parameter.
Definition: StdFace_vals.h:82
double J0pAll
Anisotropic, diagonal spin coupling (2nd Near), input parameter J0&#39;.
Definition: StdFace_vals.h:94
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 Jpp[3][3]
Isotropic, diagonal/off-diagonal spin coupling (3rd Near.), input parameter J&#39;&#39;x, J&#39;&#39;y...
Definition: StdFace_vals.h:143
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 complex t0pp
Anisotropic hopping (3rd), input parameter.
Definition: StdFace_vals.h:66
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 Vpp
Off-site Coulomb potential (3rd), input parameter.
Definition: StdFace_vals.h:86
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)...
int NCell
The number of the unit cell in the super-cell (determinant of StdIntList::box). Set in StdFace_InitSi...
Definition: StdFace_vals.h:49
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_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
double V0pp
Anisotropic Coulomb potential (3rd), input parameter.
Definition: StdFace_vals.h:79
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
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
double J0p[3][3]
Isotropic, diagonal/off-diagonal spin coupling (2nd Near.), input parameter J0&#39;x, J0&#39;y...
Definition: StdFace_vals.h:119
void StdFace_InputSpinNN(double J[3][3], double JAll, double J0[3][3], double J0All, char *J0name)
Input nearest-neighbor spin-spin interaction.
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
+ Here is the call graph for this function:
+ Here is the caller graph for this function: