HΦ  3.2.0
ChainLattice.c File Reference

Standard mode for the chain lattice. More...

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

Go to the source code of this file.

Functions

void StdFace_Chain (struct StdIntList *StdI)
 Setup a Hamiltonian for the Hubbard model on a Chain lattice. More...
 
void StdFace_Chain_Boost (struct StdIntList *StdI)
 Setup a Hamiltonian for the generalized Heisenberg model on a Chain lattice. More...
 

Detailed Description

Standard mode for the chain lattice.

Definition in file ChainLattice.c.

Function Documentation

◆ StdFace_Chain()

void StdFace_Chain ( struct StdIntList StdI)

Setup a Hamiltonian for the Hubbard model on a Chain 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

Parameters
[in,out]StdI

Definition at line 33 of file ChainLattice.c.

References StdIntList::a, 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::J2, StdIntList::J2All, StdIntList::J2p, StdIntList::J2pAll, StdIntList::JAll, StdIntList::Jp, StdIntList::JpAll, StdIntList::Jpp, StdIntList::JppAll, StdIntList::K, StdIntList::L, StdIntList::length, StdIntList::locspinflag, StdIntList::model, StdIntList::mu, 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_RequiredVal_i(), StdFace_SetLabel(), StdIntList::t, StdIntList::t0, StdIntList::t0p, StdIntList::t0pp, StdIntList::t1, StdIntList::t1p, 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::V2, StdIntList::V2p, StdIntList::Vp, StdIntList::Vpp, and StdIntList::W.

Referenced by StdFace_main().

36 {
37  FILE *fp;
38  int isite, jsite, ntransMax, nintrMax;
39  int iL;
40  double complex Cphase;
41  double dR[3];
42 
46  fp = fopen("lattice.gp", "w");
47 
48  StdI->NsiteUC = 1;
49 
50  fprintf(stdout, " @ Lattice Size & Shape\n\n");
51 
52  StdFace_PrintVal_d("a", &StdI->a, 1.0);
53  StdFace_PrintVal_d("Wlength", &StdI->length[0], StdI->a);
54  StdFace_PrintVal_d("Llength", &StdI->length[1], StdI->a);
55  StdFace_PrintVal_d("Wx", &StdI->direct[0][0], StdI->length[0]);
56  StdFace_PrintVal_d("Wy", &StdI->direct[0][1], 0.0);
57  StdFace_PrintVal_d("Lx", &StdI->direct[1][0], 0.0);
58  StdFace_PrintVal_d("Ly", &StdI->direct[1][1], StdI->length[1]);
59 
60  StdFace_PrintVal_d("phase0", &StdI->phase[0], 0.0);
61  StdFace_NotUsed_d("phase1", StdI->phase[1]);
62  StdI->phase[1] = StdI->phase[0];
63  StdI->phase[0] = 0.0;
64 
65  StdFace_RequiredVal_i("L", StdI->L);
66  StdFace_NotUsed_i("W", StdI->W);
67  StdI->W = 1;
68  StdFace_InitSite(StdI, fp, 2);
69  StdI->tau[0][0] = 0.0; StdI->tau[0][1] = 0.0; StdI->tau[0][2] = 0.0;
73  fprintf(stdout, "\n @ Hamiltonian \n\n");
74  StdFace_NotUsed_J("J1", StdI->J1All, StdI->J1);
75  StdFace_NotUsed_J("J2", StdI->J2All, StdI->J2);
76  StdFace_NotUsed_J("J1'", StdI->J1pAll, StdI->J1p);
77  StdFace_NotUsed_J("J2'", StdI->J2pAll, StdI->J2p);
78  StdFace_NotUsed_c("t1", StdI->t1);
79  StdFace_NotUsed_c("t2", StdI->t2);
80  StdFace_NotUsed_d("t1'", StdI->t1p);
81  StdFace_NotUsed_d("t2'", StdI->t2p);
82  StdFace_NotUsed_d("V1", StdI->V1);
83  StdFace_NotUsed_d("V2", StdI->V2);
84  StdFace_NotUsed_d("V1'", StdI->V1p);
85  StdFace_NotUsed_d("V2'", StdI->V2p);
86  StdFace_NotUsed_d("K", StdI->K);
87  StdFace_PrintVal_d("h", &StdI->h, 0.0);
88  StdFace_PrintVal_d("Gamma", &StdI->Gamma, 0.0);
89 
90  if (strcmp(StdI->model, "spin") == 0 ) {
91  StdFace_PrintVal_i("2S", &StdI->S2, 1);
92  StdFace_PrintVal_d("D", &StdI->D[2][2], 0.0);
93  StdFace_InputSpinNN(StdI->J, StdI->JAll, StdI->J0, StdI->J0All, "J0");
94  StdFace_InputSpinNN(StdI->Jp, StdI->JpAll, StdI->J0p, StdI->J0pAll, "J0'");
95  StdFace_InputSpinNN(StdI->Jpp, StdI->JppAll, StdI->J0pp, StdI->J0ppAll, "J0'");
96 
97  StdFace_NotUsed_d("mu", StdI->mu);
98  StdFace_NotUsed_d("U", StdI->U);
99  StdFace_NotUsed_c("t", StdI->t);
100  StdFace_NotUsed_c("t0", StdI->t0);
101  StdFace_NotUsed_c("t'", StdI->tp);
102  StdFace_NotUsed_d("V", StdI->V);
103  StdFace_NotUsed_d("V0", StdI->V0);
104  StdFace_NotUsed_d("V'", StdI->Vp);
105  }/*if (strcmp(StdI->model, "spin") == 0 )*/
106  else {
107  StdFace_PrintVal_d("mu", &StdI->mu, 0.0);
108  StdFace_PrintVal_d("U", &StdI->U, 0.0);
109  StdFace_InputHopp(StdI->t, &StdI->t0, "t0");
110  StdFace_InputHopp(StdI->tp, &StdI->t0p, "t0'");
111  StdFace_InputHopp(StdI->tpp, &StdI->t0pp, "t0''");
112  StdFace_InputCoulombV(StdI->V, &StdI->V0, "V0");
113  StdFace_InputCoulombV(StdI->Vp, &StdI->V0p, "V0'");
114  StdFace_InputCoulombV(StdI->Vpp, &StdI->V0pp, "V0''");
115 
116  StdFace_NotUsed_J("J0", StdI->J0All, StdI->J0);
117  StdFace_NotUsed_J("J0'", StdI->J0pAll, StdI->J0p);
118  StdFace_NotUsed_J("J0''", StdI->J0ppAll, StdI->J0pp);
119  StdFace_NotUsed_d("D", StdI->D[2][2]);
120 
121  if (strcmp(StdI->model, "hubbard") == 0 ) {
122  StdFace_NotUsed_i("2S", StdI->S2);
123  StdFace_NotUsed_J("J", StdI->JAll, StdI->J);
124  }
125  else if (strcmp(StdI->model, "kondo") == 0 ) {
126  StdFace_PrintVal_i("2S", &StdI->S2, 1);
127  StdFace_InputSpin(StdI->J, StdI->JAll, "J");
128  }
129  }/*if (strcmp(StdI->model, "spin") != 0 )*/
130  fprintf(stdout, "\n @ Numerical conditions\n\n");
135  StdI->nsite = StdI->L;
136  if (strcmp(StdI->model, "kondo") == 0 ) StdI->nsite *= 2;
137  StdI->locspinflag = (int *)malloc(sizeof(int) * StdI->nsite);
138 
139  if (strcmp(StdI->model, "spin") == 0 )
140  for (isite = 0; isite < StdI->nsite; isite++)StdI->locspinflag[isite] = StdI->S2;
141  else if (strcmp(StdI->model, "hubbard") == 0 )
142  for (isite = 0; isite < StdI->nsite; isite++)StdI->locspinflag[isite] = 0;
143  else if (strcmp(StdI->model, "kondo") == 0 )
144  for (isite = 0; isite < StdI->nsite / 2; isite++) {
145  StdI->locspinflag[isite] = StdI->S2;
146  StdI->locspinflag[isite + StdI->nsite / 2] = 0;
147  }
151  if (strcmp(StdI->model, "spin") == 0 ) {
152  ntransMax = StdI->L * (StdI->S2 + 1/*h*/ + 2 * StdI->S2/*Gamma*/);
153  nintrMax = StdI->L * (StdI->NsiteUC/*D*/ + 1/*J*/ + 1/*J'*/ + 1/*J''*/)
154  * (3 * StdI->S2 + 1) * (3 * StdI->S2 + 1);
155  }
156  else {
157  ntransMax = StdI->L * 2/*spin*/ * (2 * StdI->NsiteUC/*mu+h+Gamma*/ + 2/*t*/ + 2/*t'*/ + 2/*t''*/);
158  nintrMax = StdI->L * (StdI->NsiteUC/*U*/ + 4 * (1/*V*/ + 1/*V'*/ + 1/*V''*/));
159 
160  if (strcmp(StdI->model, "kondo") == 0) {
161  ntransMax += StdI->L * (StdI->S2 + 1/*h*/ + 2 * StdI->S2/*Gamma*/);
162  nintrMax += StdI->nsite / 2 * (3 * 1 + 1) * (3 * StdI->S2 + 1);
163  }/*if (strcmp(StdI->model, "kondo") == 0)*/
164  }
165 
166  StdFace_MallocInteractions(StdI, ntransMax, nintrMax);
170  for (iL = 0; iL < StdI->L; iL++){
171 
172  isite = iL;
173  if (strcmp(StdI->model, "kondo") == 0 ) isite += StdI->L;
174  /*
175  Local term
176  */
177  if (strcmp(StdI->model, "spin") == 0 ) {
178  StdFace_MagField(StdI, StdI->S2, -StdI->h, -StdI->Gamma, isite);
179  StdFace_GeneralJ(StdI, StdI->D, StdI->S2, StdI->S2, isite, isite);
180  }/*if (strcmp(StdI->model, "spin") == 0 )*/
181  else {
182  StdFace_HubbardLocal(StdI, StdI->mu, -StdI->h, -StdI->Gamma, StdI->U, isite);
183  if (strcmp(StdI->model, "kondo") == 0 ) {
184  jsite = iL;
185  StdFace_GeneralJ(StdI, StdI->J, 1, StdI->S2, isite, jsite);
186  StdFace_MagField(StdI, StdI->S2, -StdI->h, -StdI->Gamma, jsite);
187  }/*if (strcmp(StdI->model, "kondo") == 0 )*/
188  }/*if (model != "spin")*/
189  /*
190  Nearest neighbor
191  */
192  StdFace_SetLabel(StdI, fp, 0, iL, 0, 1, 0, 0, &isite, &jsite, 1, &Cphase, dR);
193 
194  if (strcmp(StdI->model, "spin") == 0 ) {
195  StdFace_GeneralJ(StdI, StdI->J0, StdI->S2, StdI->S2, isite, jsite);
196  }
197  else {
198  StdFace_Hopping(StdI, Cphase * StdI->t0, isite, jsite, dR);
199  StdFace_Coulomb(StdI, StdI->V0, isite, jsite);
200  }
201  /*
202  Second nearest neighbor
203  */
204  StdFace_SetLabel(StdI, fp, 0, iL, 0, 2, 0, 0, &isite, &jsite, 2, &Cphase, dR);
205 
206  if (strcmp(StdI->model, "spin") == 0 ) {
207  StdFace_GeneralJ(StdI, StdI->J0p, StdI->S2, StdI->S2, isite, jsite);
208  }
209  else {
210  StdFace_Hopping(StdI, Cphase * StdI->t0p, isite, jsite, dR);
211  StdFace_Coulomb(StdI, StdI->V0p, isite, jsite);
212  }
213  /*
214  Third nearest neighbor
215  */
216  StdFace_SetLabel(StdI, fp, 0, iL, 0, 3, 0, 0, &isite, &jsite, 3, &Cphase, dR);
217 
218  if (strcmp(StdI->model, "spin") == 0) {
219  StdFace_GeneralJ(StdI, StdI->J0pp, StdI->S2, StdI->S2, isite, jsite);
220  }
221  else {
222  StdFace_Hopping(StdI, Cphase * StdI->t0pp, isite, jsite, dR);
223  StdFace_Coulomb(StdI, StdI->V0pp, isite, jsite);
224  }
225  }/*for (iL = 0; iL < StdI->L; iL++)*/
226 
227  fprintf(fp, "plot \'-\' w d lc 7\n0.0 0.0\nend\npause -1\n");
228  fclose(fp);
229  StdFace_PrintGeometry(StdI);
230 }/*void StdFace_Chain*/
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].
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_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
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
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
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 J0All
Anisotropic, diagonal spin coupling (1st Near), input parameter J0.
Definition: StdFace_vals.h:92
double V1
Anisotropic Coulomb potential (1st), input parameter.
Definition: StdFace_vals.h:80
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
double phase[3]
Boundary phase, input parameter phase0, etc.
Definition: StdFace_vals.h:154
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)...
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...
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:

◆ StdFace_Chain_Boost()

void StdFace_Chain_Boost ( struct StdIntList StdI)

Setup a Hamiltonian for the generalized Heisenberg model on a Chain lattice.

Author
Mitsuaki Kawamura (The University of Tokyo)

Definition at line 237 of file ChainLattice.c.

References StdIntList::Gamma, StdIntList::h, StdIntList::ishift_nspin, StdIntList::J0, StdIntList::Jp, StdIntList::L, StdIntList::list_6spin_pair, StdIntList::list_6spin_star, StdIntList::NsiteUC, StdIntList::num_pivot, StdIntList::S2, StdFace_exit(), and StdIntList::W.

Referenced by StdFace_main().

238 {
239  int isite, ipivot;
240  int kintr;
241  FILE *fp;
242 
243  StdI->NsiteUC = 1;
244  /*
245  Magnetic field
246  */
247  fp = fopen("boost.def", "w");
248  fprintf(fp, "# Magnetic field\n");
249  fprintf(fp, "%25.15e %25.15e %25.15e\n",
250  -0.5 * StdI->Gamma, 0.0, -0.5 *StdI->h);
251  /*
252  Interaction
253  */
254  fprintf(fp, "%d # Number of type of J\n", 2);
255  fprintf(fp, "# J 1\n");
256  fprintf(fp, "%25.15e %25.15e %25.15e\n",
257  0.25 * StdI->J0[0][0], 0.25 * StdI->J0[0][1], 0.25 * StdI->J0[0][2]);
258  fprintf(fp, "%25.15e %25.15e %25.15e\n",
259  0.25 * StdI->J0[1][0], 0.25 * StdI->J0[1][1], 0.25 * StdI->J0[1][2]);
260  fprintf(fp, "%25.15e %25.15e %25.15e\n",
261  0.25 * StdI->J0[2][0], 0.25 * StdI->J0[2][1], 0.25 * StdI->J0[2][2]);
262  fprintf(fp, "# J 2\n");
263  fprintf(fp, "%25.15e %25.15e %25.15e\n",
264  0.25 * StdI->Jp[0][0], 0.25 * StdI->Jp[0][1], 0.25 * StdI->Jp[0][2]);
265  fprintf(fp, "%25.15e %25.15e %25.15e\n",
266  0.25 * StdI->Jp[1][0], 0.25 * StdI->Jp[1][1], 0.25 * StdI->Jp[1][2]);
267  fprintf(fp, "%25.15e %25.15e %25.15e\n",
268  0.25 * StdI->Jp[2][0], 0.25 * StdI->Jp[2][1], 0.25 * StdI->Jp[2][2]);
269  /*
270  Topology
271  */
272  if (StdI->S2 != 1) {
273  fprintf(stdout, "\n ERROR! S2 must be 1 in Boost. \n\n");
274  StdFace_exit(-1);
275  }
276  StdI->ishift_nspin = 4;
277  if(StdI->L % 8 != 0){
278  fprintf(stdout, "\n ERROR! L %% 8 != 0 \n\n");
279  StdFace_exit(-1);
280  }
281  StdI->W = StdI->L / 2;
282  StdI->L = 2;
283  StdI->num_pivot = StdI->W / 4;
284 
285  fprintf(fp, "# W0 R0 StdI->num_pivot StdI->ishift_nspin\n");
286  fprintf(fp, "%d %d %d %d\n", StdI->W, StdI->L, StdI->num_pivot, StdI->ishift_nspin);
287 
288  StdI->list_6spin_star = (int **)malloc(sizeof(int*) * StdI->num_pivot);
289  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
290  StdI->list_6spin_star[ipivot] = (int *)malloc(sizeof(int) * 7);
291  }
292 
293  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
294  StdI->list_6spin_star[ipivot][0] = 8; // num of J
295  StdI->list_6spin_star[ipivot][1] = 1;
296  StdI->list_6spin_star[ipivot][2] = 1;
297  StdI->list_6spin_star[ipivot][3] = 1;
298  StdI->list_6spin_star[ipivot][4] = 1;
299  StdI->list_6spin_star[ipivot][5] = 1;
300  StdI->list_6spin_star[ipivot][6] = 1; // flag
301  }
302 
303  fprintf(fp, "# StdI->list_6spin_star\n");
304  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
305  fprintf(fp, "# pivot %d\n", ipivot);
306  for (isite = 0; isite < 7; isite++) {
307  fprintf(fp, "%d ", StdI->list_6spin_star[ipivot][isite]);
308  }
309  fprintf(fp, "\n");
310  }
311 
312  StdI->list_6spin_pair = (int ***)malloc(sizeof(int**) * StdI->num_pivot);
313  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
314  StdI->list_6spin_pair[ipivot] = (int **)malloc(sizeof(int*) * 7);
315  for (isite = 0; isite < 7; isite++) {
316  StdI->list_6spin_pair[ipivot][isite] = (int *)malloc(sizeof(int) * StdI->list_6spin_star[ipivot][0]);
317  }
318  }
319 
320  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
321  StdI->list_6spin_pair[ipivot][0][0] = 0;
322  StdI->list_6spin_pair[ipivot][1][0] = 1;
323  StdI->list_6spin_pair[ipivot][2][0] = 2;
324  StdI->list_6spin_pair[ipivot][3][0] = 3;
325  StdI->list_6spin_pair[ipivot][4][0] = 4;
326  StdI->list_6spin_pair[ipivot][5][0] = 5;
327  StdI->list_6spin_pair[ipivot][6][0] = 1; // type of J
328  StdI->list_6spin_pair[ipivot][0][1] = 1;
329  StdI->list_6spin_pair[ipivot][1][1] = 2;
330  StdI->list_6spin_pair[ipivot][2][1] = 0;
331  StdI->list_6spin_pair[ipivot][3][1] = 3;
332  StdI->list_6spin_pair[ipivot][4][1] = 4;
333  StdI->list_6spin_pair[ipivot][5][1] = 5;
334  StdI->list_6spin_pair[ipivot][6][1] = 1; // type of J
335  StdI->list_6spin_pair[ipivot][0][2] = 2;
336  StdI->list_6spin_pair[ipivot][1][2] = 3;
337  StdI->list_6spin_pair[ipivot][2][2] = 0;
338  StdI->list_6spin_pair[ipivot][3][2] = 1;
339  StdI->list_6spin_pair[ipivot][4][2] = 4;
340  StdI->list_6spin_pair[ipivot][5][2] = 5;
341  StdI->list_6spin_pair[ipivot][6][2] = 1; // type of J
342  StdI->list_6spin_pair[ipivot][0][3] = 3;
343  StdI->list_6spin_pair[ipivot][1][3] = 4;
344  StdI->list_6spin_pair[ipivot][2][3] = 0;
345  StdI->list_6spin_pair[ipivot][3][3] = 1;
346  StdI->list_6spin_pair[ipivot][4][3] = 2;
347  StdI->list_6spin_pair[ipivot][5][3] = 5;
348  StdI->list_6spin_pair[ipivot][6][3] = 1; // type of J
349  StdI->list_6spin_pair[ipivot][0][4] = 0;
350  StdI->list_6spin_pair[ipivot][1][4] = 2;
351  StdI->list_6spin_pair[ipivot][2][4] = 1;
352  StdI->list_6spin_pair[ipivot][3][4] = 3;
353  StdI->list_6spin_pair[ipivot][4][4] = 4;
354  StdI->list_6spin_pair[ipivot][5][4] = 5;
355  StdI->list_6spin_pair[ipivot][6][4] = 2; // type of J
356  StdI->list_6spin_pair[ipivot][0][5] = 1;
357  StdI->list_6spin_pair[ipivot][1][5] = 3;
358  StdI->list_6spin_pair[ipivot][2][5] = 0;
359  StdI->list_6spin_pair[ipivot][3][5] = 2;
360  StdI->list_6spin_pair[ipivot][4][5] = 4;
361  StdI->list_6spin_pair[ipivot][5][5] = 5;
362  StdI->list_6spin_pair[ipivot][6][5] = 2; // type of J
363  StdI->list_6spin_pair[ipivot][0][6] = 2;
364  StdI->list_6spin_pair[ipivot][1][6] = 4;
365  StdI->list_6spin_pair[ipivot][2][6] = 0;
366  StdI->list_6spin_pair[ipivot][3][6] = 1;
367  StdI->list_6spin_pair[ipivot][4][6] = 3;
368  StdI->list_6spin_pair[ipivot][5][6] = 5;
369  StdI->list_6spin_pair[ipivot][6][6] = 2; // type of J
370  StdI->list_6spin_pair[ipivot][0][7] = 3;
371  StdI->list_6spin_pair[ipivot][1][7] = 5;
372  StdI->list_6spin_pair[ipivot][2][7] = 0;
373  StdI->list_6spin_pair[ipivot][3][7] = 1;
374  StdI->list_6spin_pair[ipivot][4][7] = 2;
375  StdI->list_6spin_pair[ipivot][5][7] = 4;
376  StdI->list_6spin_pair[ipivot][6][7] = 2; // type of J
377  }
378 
379  fprintf(fp, "# StdI->list_6spin_pair\n");
380  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
381  fprintf(fp, "# pivot %d\n", ipivot);
382  for (kintr = 0; kintr < StdI->list_6spin_star[ipivot][0]; kintr++) {
383  for (isite = 0; isite < 7; isite++) {
384  fprintf(fp, "%d ", StdI->list_6spin_pair[ipivot][isite][kintr]);
385  }
386  fprintf(fp, "\n");
387  }
388  }
389  fclose(fp);
390 
391  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
392  free(StdI->list_6spin_star[ipivot]);
393  }
394  free(StdI->list_6spin_star);
395 
396  for (ipivot = 0; ipivot < StdI->num_pivot; ipivot++) {
397  for (isite = 0; isite < 7; isite++) {
398  free(StdI->list_6spin_pair[ipivot][isite]);
399  }
400  free(StdI->list_6spin_pair[ipivot]);
401  }
402  free(StdI->list_6spin_pair);
403 
404 }
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 L
Number of sites along the 2nd axis, input parameter.
Definition: StdFace_vals.h:40
int ** list_6spin_star
Definition: StdFace_vals.h:279
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
int W
Number of sites along the 1st axis, input parameter.
Definition: StdFace_vals.h:39
double Gamma
Transvars magnetic field, input parameter.
Definition: StdFace_vals.h:148
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
int ishift_nspin
Definition: StdFace_vals.h:281
int *** list_6spin_pair
Definition: StdFace_vals.h:278
double h
Longitudinal magnetic field, input parameter.
Definition: StdFace_vals.h:147
void StdFace_exit(int errorcode)
MPI Abortation wrapper.
+ Here is the call graph for this function:
+ Here is the caller graph for this function: