HΦ  3.1.0
StdFace_ModelUtil.c File Reference

Various utility for constructing models. More...

#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <complex.h>
#include "StdFace_vals.h"
#include <string.h>
#include <mpi.h>

Go to the source code of this file.

Functions

void StdFace_exit (int errorcode)
 MPI Abortation wrapper. More...
 
void StdFace_trans (struct StdIntList *StdI, double complex trans0, int isite, int ispin, int jsite, int jspin)
 Add transfer to the list set StdIntList::trans and StdIntList::transindx and increment StdIntList::ntrans. More...
 
void StdFace_Hopping (struct StdIntList *StdI, double complex trans0, int isite, int jsite, double *dR)
 Add Hopping for the both spin. More...
 
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. More...
 
void StdFace_MagField (struct StdIntList *StdI, int S2, double h, double Gamma, int isite)
 Add longitudinal and transvars magnetic field to the list. More...
 
void StdFace_intr (struct StdIntList *StdI, double complex intr0, int site1, int spin1, int site2, int spin2, int site3, int spin3, int site4, int spin4)
 Add interaction (InterAll) to the list Set StdIntList::intr and StdIntList::intrindx and increase the number of that (StdIntList::nintr). More...
 
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' + 1) interactions]. More...
 
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). More...
 
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), set the default value. More...
 
void StdFace_PrintVal_dd (char *valname, double *val, double val0, double val1)
 Print a valiable (real) read from the input file if it is not specified in the input file (=NaN), set the default value. More...
 
void StdFace_PrintVal_c (char *valname, double complex *val, double complex val0)
 Print a valiable (complex) read from the input file if it is not specified in the input file (=NaN), set the default value. More...
 
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 (=2147483647, the upper limt of Int) set the default value. More...
 
void StdFace_NotUsed_d (char *valname, double val)
 Stop HPhi if a variable (real) not used is specified in the input file (!=NaN). More...
 
void StdFace_NotUsed_c (char *valname, double complex val)
 Stop HPhi if a variable (complex) not used is specified in the input file (!=NaN). More...
 
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). More...
 
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). More...
 
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, the upper limt of Int). More...
 
static void StdFace_FoldSite (struct StdIntList *StdI, int iCellV[3], int nBox[3], int iCellV_fold[3])
 Move a site into the original supercell if it is outside the original supercell. More...
 
void StdFace_InitSite (struct StdIntList *StdI, FILE *fp, int dim)
 Initialize the super-cell where simulation is performed. More...
 
void StdFace_FindSite (struct StdIntList *StdI, int iW, int iL, int iH, int diW, int diL, int diH, int isiteUC, int jsiteUC, int *isite, int *jsite, double complex *Cphase, double *dR)
 Find the index of transfer and interaction. More...
 
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) More...
 
void StdFace_PrintXSF (struct StdIntList *StdI)
 Print lattice.xsf (XCrysDen format) More...
 
void StdFace_InputSpinNN (struct StdIntList *StdI, double J0[3][3], double J0All, char *J0name)
 Input nearest-neighbor spin-spin interaction. More...
 
void StdFace_InputSpin (struct StdIntList *StdI, double Jp[3][3], double JpAll, char *Jpname)
 Input spin-spin interaction other than nearest-neighbor. More...
 
void StdFace_InputCoulombV (struct StdIntList *StdI, double *V0, char *V0name)
 Input off-site Coulomb interaction from the input file, if it is not specified, use the default value (0 or the isotropic Coulomb interaction StdIntList::V). More...
 
void StdFace_InputHopp (struct StdIntList *StdI, double complex *t0, char *t0name)
 Input hopping integral from the input file, if it is not specified, use the default value(0 or the isotropic hopping StdIntList::V). More...
 
void StdFace_PrintGeometry (struct StdIntList *StdI)
 Print geometry of sites for the pos-process of correlation function. More...
 
void StdFace_MallocInteractions (struct StdIntList *StdI, int ntransMax, int nintrMax)
 Malloc Arrays for interactions. More...
 

Detailed Description

Various utility for constructing models.

Definition in file StdFace_ModelUtil.c.

Function Documentation

◆ StdFace_Coulomb()

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).

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in,out]StdI
[in]VCoulomb integral U, V, etc.
[in]isitei of n_i
[in]jsitej of n_j

Definition at line 396 of file StdFace_ModelUtil.c.

References StdIntList::Cinter, StdIntList::CinterIndx, StdIntList::NCinter, and StdIntList::V.

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), StdFace_Triangular(), and StdFace_Wannier90().

402 {
403  StdI->Cinter[StdI->NCinter] = V;
404  StdI->CinterIndx[StdI->NCinter][0] = isite;
405  StdI->CinterIndx[StdI->NCinter][1] = jsite;
406  StdI->NCinter += 1;
407 }/*void StdFace_Coulomb*/
int NCinter
Number of inter-site Coulomb interaction, counted in each lattice file.
Definition: StdFace_vals.h:167
double V
Off-site Coulomb potential (1st), input parameter.
Definition: StdFace_vals.h:72
double * Cinter
[StdIntList::NCinter] Coefficient of inter-site Coulomb term, malloc in StdFace_MallocInteractions() ...
Definition: StdFace_vals.h:173
int ** CinterIndx
[StdIntList::NCinter][2] Site indices of inter-site Coulomb term, malloc in StdFace_MallocInteraction...
Definition: StdFace_vals.h:170

◆ StdFace_exit()

void StdFace_exit ( int  errorcode)

MPI Abortation wrapper.

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]errorcode

Definition at line 35 of file StdFace_ModelUtil.c.

Referenced by CheckOutputMode(), PrintCalcMod(), PrintExcitation(), StdFace_Chain_Boost(), StdFace_Honeycomb_Boost(), StdFace_InitSite(), StdFace_InputCoulombV(), StdFace_InputHopp(), StdFace_InputSpin(), StdFace_InputSpinNN(), StdFace_Kagome_Boost(), StdFace_Ladder_Boost(), StdFace_main(), StdFace_NotUsed_c(), StdFace_NotUsed_d(), StdFace_NotUsed_i(), StdFace_RequiredVal_i(), StdFace_Wannier90(), StoreWithCheckDup_c(), StoreWithCheckDup_d(), StoreWithCheckDup_i(), StoreWithCheckDup_s(), StoreWithCheckDup_sl(), UnsupportedSystem(), and VectorPotential().

37 {
38  int ierr = 0;
39  fflush(stdout);
40  fflush(stderr);
41 #ifdef MPI
42  fprintf(stdout, "\n\n ####### You DO NOT have to WORRY about the following MPI-ERROR MESSAGE. #######\n\n");
43  ierr = MPI_Abort(MPI_COMM_WORLD, errorcode);
44  ierr = MPI_Finalize();
45 #endif
46  if (ierr != 0) fprintf(stderr, "\n MPI_Finalize() = %d\n\n", ierr);
47  exit(errorcode);
48 }

◆ StdFace_FindSite()

void StdFace_FindSite ( struct StdIntList StdI,
int  iW,
int  iL,
int  iH,
int  diW,
int  diL,
int  diH,
int  isiteUC,
int  jsiteUC,
int *  isite,
int *  jsite,
double complex *  Cphase,
double *  dR 
)

Find the index of transfer and interaction.

Parameters
[in,out]StdI
[in]iWposition of initial site
[in]iLposition of initial site
[in]iHposition of initial site
[in]diWTranslation from the initial site
[in]diLTranslation from the initial site
[in]diHTranslation from the initial site
[in]isiteUCIntrinsic site index of initial site
[in]jsiteUCIntrinsic site index of final site
[out]isiteinitial site
[out]jsitefinal site
[out]CphaseBoundary phase, if it across boundary
[out]dRR_i - R_j in the fractional coordinate

Definition at line 823 of file StdFace_ModelUtil.c.

References StdIntList::Cell, StdIntList::ExpPhase, StdIntList::model, StdIntList::NCell, StdIntList::NsiteUC, StdFace_FoldSite(), and StdIntList::tau.

Referenced by StdFace_FCOrtho(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_SetLabel(), and StdFace_Wannier90().

838 {
839  int iCell, jCell, kCell, ii;
840  int nBox[3], jCellV[3];
841 
842  dR[0] = - (double)diW + StdI->tau[isiteUC][0] - StdI->tau[jsiteUC][0];
843  dR[1] = - (double)diL + StdI->tau[isiteUC][1] - StdI->tau[jsiteUC][1];
844  dR[2] = - (double)diH + StdI->tau[isiteUC][2] - StdI->tau[jsiteUC][2];
845 
846  jCellV[0] = iW + diW;
847  jCellV[1] = iL + diL;
848  jCellV[2] = iH + diH;
849  StdFace_FoldSite(StdI, jCellV, nBox, jCellV);
850  *Cphase = 1.0;
851  for (ii = 0; ii < 3; ii++) *Cphase *= cpow(StdI->ExpPhase[ii], (double)nBox[ii]);
852 
853  for (kCell = 0; kCell < StdI->NCell; kCell++) {
854  if (jCellV[0] == StdI->Cell[kCell][0] &&
855  jCellV[1] == StdI->Cell[kCell][1] &&
856  jCellV[2] == StdI->Cell[kCell][2])
857  {
858  jCell = kCell;
859  }
860  if (iW == StdI->Cell[kCell][0] &&
861  iL == StdI->Cell[kCell][1] &&
862  iH == StdI->Cell[kCell][2])
863  {
864  iCell = kCell;
865  }
866  }/*for (iCell = 0; iCell < StdI->NCell; iCell++)*/
867  *isite = iCell * StdI->NsiteUC + isiteUC;
868  *jsite = jCell * StdI->NsiteUC + jsiteUC;
869  if (strcmp(StdI->model, "kondo") == 0) {
870  *isite += StdI->NCell * StdI->NsiteUC;
871  *jsite += StdI->NCell * StdI->NsiteUC;
872  }
873 }/*void StdFace_FindSite*/
double complex ExpPhase[3]
.
Definition: StdFace_vals.h:134
int NsiteUC
Number of sites in the unit cell. Defined in the beginning of each lattice function.
Definition: StdFace_vals.h:53
char model[256]
Name of model, input parameter.
Definition: StdFace_vals.h:60
int ** Cell
[StdIntList][3] The cell position in the fractional coordinate. Malloc and Set in StdFace_InitSite()...
Definition: StdFace_vals.h:51
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
static void StdFace_FoldSite(struct StdIntList *StdI, int iCellV[3], int nBox[3], int iCellV_fold[3])
Move a site into the original supercell if it is outside the original supercell.
double ** tau
Cell-internal site position in the fractional coordinate. Defined in the beginning of each lattice fu...
Definition: StdFace_vals.h:55

◆ StdFace_FoldSite()

static void StdFace_FoldSite ( struct StdIntList StdI,
int  iCellV[3],
int  nBox[3],
int  iCellV_fold[3] 
)
static

Move a site into the original supercell if it is outside the original supercell.

Author
Mitsuaki Kawamura (The University of Tokyo)

(1) Transform to fractional coordinate (times NCell).

(2) Search which supercell contains this cell.

(3) Fractional coordinate (times NCell) in the original supercell

Parameters
[in,out]StdI
[in]iCellVThe fractional coordinate of a site
[out]nBoxthe index of supercell
[out]iCellV_foldThe fractional coordinate of a site which is moved into the original cell

Definition at line 595 of file StdFace_ModelUtil.c.

References StdIntList::box, StdIntList::NCell, and StdIntList::rbox.

Referenced by StdFace_FindSite(), and StdFace_InitSite().

602 {
603  int ii, jj, iCellV_frac[3];
607  for (ii = 0; ii < 3; ii++) {
608  iCellV_frac[ii] = 0;
609  for (jj = 0; jj < 3; jj++)iCellV_frac[ii] += StdI->rbox[ii][jj] * iCellV[jj];
610  }
614  for (ii = 0; ii < 3; ii++)
615  nBox[ii] = (iCellV_frac[ii] + StdI->NCell * 1000) / StdI->NCell - 1000;
619  for (ii = 0; ii < 3; ii++)
620  iCellV_frac[ii] -= StdI->NCell*(nBox[ii]);
621 
622  for (ii = 0; ii < 3; ii++) {
623  iCellV_fold[ii] = 0;
624  for (jj = 0; jj < 3; jj++) iCellV_fold[ii] += StdI->box[jj][ii] * iCellV_frac[jj];
625  iCellV_fold[ii] = (iCellV_fold[ii] + StdI->NCell * 1000) / StdI->NCell - 1000;
626  }/*for (ii = 0; ii < 3; ii++)*/
627 }/*static void StdFace_FoldSite*/
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
int rbox[3][3]
The inversion of StdIntList::box. Set in StdFace_InitSite().
Definition: StdFace_vals.h:47
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

◆ StdFace_GeneralJ()

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' + 1) interactions].

Author
Mitsuaki Kawamura (The University of Tokyo)

If both spins are S=1/2 ...

Set StdIntList::Hund, StdIntList::HundIndx, StdIntList::Cinter, StdIntList::CinterIndx for \(S_{iz} S_{jz}\) term, and increase the number of them (StdIntList::NHund, StdIntList::NCinter). And ...

If there is no off-diagonal term, set StdIntList::Ex, StdIntList::ExIndx, StdIntList::PairLift, StdIntList::PLIndx and increase the number of them (StdIntList::NEx, StdIntList::NPairLift).

If S != 1/2 or there is an off-diagonal interaction, use the InterAll format.

(1)

\[ J_z S_{i z} * S_{j z} = J_z \sum_{\sigma, \sigma' = -S}^S \sigma \sigma' c_{i\sigma}^\dagger c_{i\sigma} c_{j\sigma'}^\dagger c_{j\sigma'} \]

(2)

\[ I S_i^+ S_j^- + I^* S_j^+ S_i^- = \sum_{\sigma, \sigma' = -S}^{S-1} \sqrt{S_i(S_i+1) - \sigma(\sigma+1)}\sqrt{S_j(S_j+1) - \sigma'(\sigma'+1)} (I c_{i\sigma+1}^\dagger c_{i\sigma} c_{j\sigma'}^\dagger c_{j\sigma'+1} + I^* c_{i\sigma}^\dagger c_{i\sigma+1} c_{j\sigma'+1}^\dagger c_{j\sigma'}) \\ I \equiv \frac{J_x + J_y + i(J_{xy} - J_{yx})}{4} \]

(3)

\[ I S_i^+ S_j^+ + I^* S_j^- S_i^- = \sum_{\sigma, \sigma' = -S}^{S-1} \sqrt{S_i(S_i+1) - \sigma(\sigma+1)}\sqrt{S_j(S_j+1) - \sigma'(\sigma'+1)} (I c_{i\sigma+1}^\dagger c_{i\sigma} c_{j\sigma'+1}^\dagger c_{j\sigma'} + I^* c_{i\sigma}^\dagger c_{i\sigma+1} c_{j\sigma'}^\dagger c_{j\sigma'+1}) \\ I \equiv \frac{J_x - J_y - i(J_{xy} + J_{yx})}{4} \]

(4)

\[ I S_i^+ S_{j z} + I^* S_{j z} S_i^-= \sum_{\sigma=-S}^{S-1} \sum_{\sigma' = -S}^{S} \sqrt{S_i(S_i+1) - \sigma(\sigma+1)} (I c_{i\sigma+1}^\dagger c_{i\sigma} c_{j\sigma'}^\dagger c_{j\sigma'} + I^* c_{j\sigma'}^\dagger c_{j\sigma'} c_{i\sigma}^\dagger c_{i\sigma+1}) \\ I \equiv \frac{J_{xz} - i J_{yz}}{2} \]

(5)

\[ I S_{i z} S_j^+ + I^* S_j^- S_{i z} = \sum_{\sigma=-S}^{S} \sum_{\sigma' = -S}^{S-1} \sqrt{S_j(S_j+1) - \sigma'(\sigma'+1)} (I c_{i\sigma}^\dagger c_{i\sigma} c_{j\sigma'+1}^\dagger c_{j\sigma'} + I^* c_{j\sigma'}^\dagger c_{j\sigma'+1} c_{i\sigma}^\dagger c_{i\sigma}) \\ I \equiv \frac{J_{zx} - i J_{zy}}{2} \]

Parameters
[in,out]StdI
[in]JThe Spin interaction \(J_x, J_{xy}\), ...
[in]Si2Spin moment in \(i\) site
[in]Sj2Spin moment in \(j\) site
[in]isite\(i\) of \(S_i\)
[in]jsite\(j\) of \(S_j\)

Definition at line 223 of file StdFace_ModelUtil.c.

References StdIntList::Cinter, StdIntList::CinterIndx, StdIntList::Ex, StdIntList::ExIndx, StdIntList::Hund, StdIntList::HundIndx, StdIntList::J, StdIntList::model, StdIntList::NCinter, StdIntList::NEx, StdIntList::NHund, StdIntList::NPairLift, StdIntList::PairLift, StdIntList::PLIndx, and StdFace_intr().

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), StdFace_Triangular(), and StdFace_Wannier90().

231 {
232  int ispin, jspin, ZGeneral, ExGeneral;
233  double Si, Sj, Siz, Sjz;
234  double complex intr0;
238  ZGeneral = 1;
239  ExGeneral = 1;
240  if (Si2 == 1 || Sj2 == 1) {
247  ZGeneral = 0;
248 
249  StdI->Hund[StdI->NHund] = -0.5 * J[2][2];
250  StdI->HundIndx[StdI->NHund][0] = isite;
251  StdI->HundIndx[StdI->NHund][1] = jsite;
252  StdI->NHund += 1;
253 
254  StdI->Cinter[StdI->NCinter] = -0.25 * J[2][2];
255  StdI->CinterIndx[StdI->NCinter][0] = isite;
256  StdI->CinterIndx[StdI->NCinter][1] = jsite;
257  StdI->NCinter += 1;
258 
259  if (fabs(J[0][1]) < 0.000001 && fabs(J[1][0]) < 0.000001
260 #if defined(_mVMC)
261  && abs(J[0][0] - J[1][1]) < 0.000001
262 #endif
263  ) {
270  ExGeneral = 0;
271 
272 #if defined(_mVMC)
273  StdI->Ex[StdI->NEx] = - 0.25 * (J[0][0] + J[1][1]);
274 #else
275  if (strcmp(StdI->model, "kondo") == 0)
276  StdI->Ex[StdI->NEx] = -0.25 * (J[0][0] + J[1][1]);
277  else
278  StdI->Ex[StdI->NEx] = 0.25 * (J[0][0] + J[1][1]);
279 #endif
280  StdI->ExIndx[StdI->NEx][0] = isite;
281  StdI->ExIndx[StdI->NEx][1] = jsite;
282  StdI->NEx += 1;
283 
284  StdI->PairLift[StdI->NPairLift] = 0.25 * (J[0][0] - J[1][1]);
285  StdI->PLIndx[StdI->NPairLift][0] = isite;
286  StdI->PLIndx[StdI->NPairLift][1] = jsite;
287  StdI->NPairLift += 1;
288  }/*if (fabs(J[0][1]) < 0.000001 && fabs(J[1][0]) < 0.000001)*/
289  }
294  Si = 0.5 * (double)Si2;
295  Sj = 0.5 * (double)Sj2;
296 
297  for (ispin = 0; ispin <= Si2; ispin++) {
298  Siz = (double)ispin - Si;
299  for (jspin = 0; jspin <= Sj2; jspin++) {
300  Sjz = (double)jspin - Sj;
308  if (ZGeneral == 1) {
309  intr0 = J[2][2] * Siz * Sjz;
310  StdFace_intr(StdI, intr0,
311  isite, ispin, isite, ispin, jsite, jspin, jsite, jspin);
312  }
324  if ((ispin < Si2 && jspin < Sj2) && ExGeneral == 1) {
325  intr0 = 0.25 * (J[0][0] + J[1][1] + I*(J[0][1] - J[1][0]))
326  * sqrt(Si * (Si + 1.0) - Siz * (Siz + 1.0))
327  * sqrt(Sj * (Sj + 1.0) - Sjz * (Sjz + 1.0));
328  StdFace_intr(StdI, intr0,
329  isite, ispin + 1, isite, ispin, jsite, jspin, jsite, jspin + 1);
330  StdFace_intr(StdI, conj(intr0),
331  isite, ispin, isite, ispin + 1, jsite, jspin + 1, jsite, jspin);
332  }
344  if ((ispin < Si2 && jspin < Sj2) && ExGeneral == 1) {
345  intr0 = 0.5 * 0.5 * (J[0][0] - J[1][1] - I*(J[0][1] + J[1][0]))
346  * sqrt(Si * (Si + 1.0) - Siz * (Siz + 1.0))
347  * sqrt(Sj * (Sj + 1.0) - Sjz * (Sjz + 1.0));
348  StdFace_intr(StdI, intr0,
349  isite, ispin + 1, isite, ispin, jsite, jspin + 1, jsite, jspin);
350  StdFace_intr(StdI, conj(intr0),
351  isite, ispin, isite, ispin + 1, jsite, jspin, jsite, jspin + 1);
352  }
363  if (ispin < Si2) {
364  intr0 = 0.5 * (J[0][2] - I * J[1][2]) * sqrt(Si * (Si + 1.0) - Siz * (Siz + 1.0)) * Sjz;
365  StdFace_intr(StdI, intr0,
366  isite, ispin + 1, isite, ispin, jsite, jspin, jsite, jspin);
367  StdFace_intr(StdI, conj(intr0),
368  jsite, jspin, jsite, jspin, isite, ispin, isite, ispin + 1);
369  }/*if (ispin < Si2)*/
380  if (jspin < Sj2) {
381  intr0 = 0.5 * (J[2][0] - I * J[2][1]) * Siz * sqrt(Sj * (Sj + 1.0) - Sjz * (Sjz + 1.0));
382  StdFace_intr(StdI, intr0,
383  isite, ispin, isite, ispin, jsite, jspin + 1, jsite, jspin);
384  StdFace_intr(StdI, conj(intr0),
385  jsite, jspin, jsite, jspin + 1, isite, ispin, isite, ispin);
386  }/*if (jspin < Sj2)*/
387  }/*for (jspin = 0; jspin <= Sj2; jspin++)*/
388  }/*for (ispin = 0; ispin <= Si2; ispin++)*/
389 }/*StdFace_GeneralJ*/
int NHund
Number of Hund term, counted in each lattice file.
Definition: StdFace_vals.h:176
double J[3][3]
Isotropic, diagonal/off-diagonal spin coupling (1st Near.), input parameter Jx, Jy, Jz, Jxy, etc.
Definition: StdFace_vals.h:100
double * Ex
[StdIntList::NEx] Coefficient of exchange term, malloc in StdFace_MallocInteractions() and set in Std...
Definition: StdFace_vals.h:189
int ** PLIndx
[StdIntList::NPairLift][2] Site indices of pair-lift term, malloc in StdFace_MallocInteractions() and...
Definition: StdFace_vals.h:194
int ** ExIndx
[StdIntList::NEx][2] Site indices of exchange term, malloc in StdFace_MallocInteractions() and set in...
Definition: StdFace_vals.h:186
char model[256]
Name of model, input parameter.
Definition: StdFace_vals.h:60
int NCinter
Number of inter-site Coulomb interaction, counted in each lattice file.
Definition: StdFace_vals.h:167
double * Hund
[StdIntList::NHund] Coefficient of Hund term, malloc in StdFace_MallocInteractions() and set in StdFa...
Definition: StdFace_vals.h:181
double * PairLift
[StdIntList::NPairLift] Coefficient of pair-lift term, malloc in StdFace_MallocInteractions() and set...
Definition: StdFace_vals.h:197
int NPairLift
Number of pair-lift term, counted in each lattice file.
Definition: StdFace_vals.h:192
int ** HundIndx
[StdIntList::NHund][2] Site indices of Hund term, malloc in StdFace_MallocInteractions() and set in S...
Definition: StdFace_vals.h:178
int NEx
Number of exchange term, counted in each lattice file.
Definition: StdFace_vals.h:184
void StdFace_intr(struct StdIntList *StdI, double complex intr0, int site1, int spin1, int site2, int spin2, int site3, int spin3, int site4, int spin4)
Add interaction (InterAll) to the list Set StdIntList::intr and StdIntList::intrindx and increase the...
double * Cinter
[StdIntList::NCinter] Coefficient of inter-site Coulomb term, malloc in StdFace_MallocInteractions() ...
Definition: StdFace_vals.h:173
int ** CinterIndx
[StdIntList::NCinter][2] Site indices of inter-site Coulomb term, malloc in StdFace_MallocInteraction...
Definition: StdFace_vals.h:170

◆ StdFace_Hopping()

void StdFace_Hopping ( struct StdIntList StdI,
double complex  trans0,
int  isite,
int  jsite,
double *  dR 
)

Add Hopping for the both spin.

Author
Mitsuaki Kawamura (The University of Tokyo)

Both \(c_{i \sigma}^\dagger c_{j \sigma}\) and \(c_{j \sigma}^\dagger c_{i \sigma}\) for every spin channel ( \(\sigma\)) is specified

Parameters
[in,out]StdI
[in]trans0Hopping integral \(t\)
[in]isite\(i\) for \(c_{i \sigma}^\dagger\)
[in]jsite\(j\) for \(c_{j \sigma}\)
[in]dRR_i - R_j

Definition at line 75 of file StdFace_ModelUtil.c.

References StdIntList::At, StdIntList::Lanczos_max, StdIntList::method, StdIntList::npump, StdIntList::pump, StdIntList::PumpBody, StdIntList::pumpindx, and StdFace_trans().

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), StdFace_Triangular(), and StdFace_Wannier90().

82 {
83  int ispin, it, ii;
84  double complex Cphase, coef;
90 #if defined(_HPhi)
91  if (strcmp(StdI->method, "timeevolution") == 0 && StdI->PumpBody == 1) {
92  for (it = 0; it < StdI->Lanczos_max; it++) {
93  Cphase = 0.0f;
94  for (ii = 0; ii < 3; ii++) Cphase += /*2.0*StdI->pi */ StdI->At[it][ii] * dR[ii];
95  coef = cos(Cphase) + I * sin(-Cphase);
96  for (ispin = 0; ispin < 2; ispin++) {
97  StdI->pump[it][StdI->npump[it]] = coef * trans0;
98  StdI->pumpindx[it][StdI->npump[it]][0] = isite;
99  StdI->pumpindx[it][StdI->npump[it]][1] = ispin;
100  StdI->pumpindx[it][StdI->npump[it]][2] = jsite;
101  StdI->pumpindx[it][StdI->npump[it]][3] = ispin;
102  StdI->npump[it] = StdI->npump[it] + 1;
103 
104  StdI->pump[it][StdI->npump[it]] = conj(coef * trans0);
105  StdI->pumpindx[it][StdI->npump[it]][0] = jsite;
106  StdI->pumpindx[it][StdI->npump[it]][1] = ispin;
107  StdI->pumpindx[it][StdI->npump[it]][2] = isite;
108  StdI->pumpindx[it][StdI->npump[it]][3] = ispin;
109  StdI->npump[it] = StdI->npump[it] + 1;
110  }/*for (ispin = 0; ispin < 2; ispin++)*/
111  }/*for (it = 0; it < StdI->Lanczos_max; it++)*/
112  }/*if (strcmp(StdI->method, "timeevolution") == 0)*/
113  else
114 #endif
115  for (ispin = 0; ispin < 2; ispin++) {
116  StdFace_trans(StdI, trans0, jsite, ispin, isite, ispin);
117  StdFace_trans(StdI, conj(trans0), isite, ispin, jsite, ispin);
118  }/*for (ispin = 0; ispin < 2; ispin++)*/
119 }/*void StdFace_Hopping*/
void StdFace_trans(struct StdIntList *StdI, double complex trans0, int isite, int ispin, int jsite, int jspin)
Add transfer to the list set StdIntList::trans and StdIntList::transindx and increment StdIntList::nt...
double ** At
[StdIntList::nt][3] Vector potential.
Definition: StdFace_vals.h:285
int *** pumpindx
[StdIntList::nt][StdIntList::npump][4] Site/spin indices of one-body term, malloc in StdFace_MallocIn...
Definition: StdFace_vals.h:279
int PumpBody
one- or two-body pumping, defined from StdIntList::PumpType
Definition: StdFace_vals.h:276
int * npump
[StdIntList::nt] Number of transfer, counted in each lattice file.
Definition: StdFace_vals.h:278
int Lanczos_max
The maxixmum number of iterations, input from file.
Definition: StdFace_vals.h:237
double complex ** pump
[StdIntList::nt][StdIntList::npump] Coefficient of one-body term, malloc in StdFace_MallocInteraction...
Definition: StdFace_vals.h:282
char method[256]
The name of method, input from file.
Definition: StdFace_vals.h:232

◆ StdFace_HubbardLocal()

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.

Author
Mitsuaki Kawamura (The University of Tokyo)

Set StdIntList::Cintra and StdIntList::CintraIndx with the input argument and increase the number of them (StdIntList::NCintra).

Parameters
[in,out]StdI
[in]mu0Chemical potential
[in]h0Longitudinal magnetic feild
[in]Gamma0Transvers magnetic feild
[in]U0Intra-site Coulomb potential
[in]isitei for \(c_{i \sigma}^\dagger\)

Definition at line 125 of file StdFace_ModelUtil.c.

References StdIntList::Cintra, StdIntList::CintraIndx, StdIntList::NCintra, and StdFace_trans().

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), StdFace_Triangular(), and StdFace_Wannier90().

133 {
134  StdFace_trans(StdI, mu0 + 0.5 * h0, isite, 0, isite, 0);
135  StdFace_trans(StdI, mu0 - 0.5 * h0, isite, 1, isite, 1);
136  StdFace_trans(StdI, -0.5 * Gamma0, isite, 1, isite, 0);
137  StdFace_trans(StdI, -0.5 * Gamma0, isite, 0, isite, 1);
143  StdI->Cintra[StdI->NCintra] = U0;
144  StdI->CintraIndx[StdI->NCintra][0] = isite;
145  StdI->NCintra += 1;
146 }/*void StdFace_HubbardLocal*/
void StdFace_trans(struct StdIntList *StdI, double complex trans0, int isite, int ispin, int jsite, int jspin)
Add transfer to the list set StdIntList::trans and StdIntList::transindx and increment StdIntList::nt...
int ** CintraIndx
[StdIntList::NCintra][1] Site indices of intra-site Coulomb term, malloc in StdFace_MallocInteraction...
Definition: StdFace_vals.h:161
int NCintra
Number of intra-site Coulomb interaction, counted in each lattice file.
Definition: StdFace_vals.h:158
double * Cintra
[StdIntList::NCintra] Coefficient of intra-site Coulomb term, malloc in StdFace_MallocInteractions() ...
Definition: StdFace_vals.h:164

◆ StdFace_InitSite()

void StdFace_InitSite ( struct StdIntList StdI,
FILE *  fp,
int  dim 
)

Initialize the super-cell where simulation is performed.

Author
Mitsuaki Kawamura (The University of Tokyo)

(1) Check Input parameters about the shape of super-cell

(2) Define the phase factor at each boundary. Set anti-period flag (StdIntList::AntiPeriod).

(3) Malloc StdIntList::tau, intrinsic structure of unit-cell

(4) Calculate reciprocal lattice vectors (times StdIntList::NCell, the number of unit-cells) for folding sites. store in StdIntList::rbox.

(5) Find Cells in the super-cell (5-1) Find the lower- and the upper bound for surrowndinf the super-cell

(5-2) Find cells in the super-cell (StdIntList::Cell, the fractional coordinate)

(6) For 2D system, print header of lattice.gp

Parameters
[in,out]StdI
[in]fpFile pointer to lattice.gp
[in]dimdimension of system, if = 2, print lattice.gp

Definition at line 632 of file StdFace_ModelUtil.c.

References StdIntList::AntiPeriod, StdIntList::box, StdIntList::Cell, StdIntList::direct, StdIntList::ExpPhase, StdIntList::Height, StdIntList::L, StdIntList::NaN_i, StdIntList::NCell, StdIntList::NsiteUC, StdIntList::phase, StdIntList::pi180, StdIntList::rbox, StdFace_exit(), StdFace_FoldSite(), StdFace_PrintVal_i(), StdIntList::tau, and StdIntList::W.

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), StdFace_Triangular(), and StdFace_Wannier90().

637 {
638  int bound[3][2], edge, ii, jj;
639  int ipos;
640  int nBox[3], iCellV_fold[3], iCellV[3];
641  double pos[4][2], xmin, xmax/*, offset[2], scale*/;
642 
643  fprintf(stdout, "\n @ Super-Lattice setting\n\n");
647  if (
648  (StdI->L != StdI->NaN_i || StdI->W != StdI->NaN_i || StdI->Height != StdI->NaN_i)
649  &&
650  (StdI->box[0][0] != StdI->NaN_i || StdI->box[0][1] != StdI->NaN_i || StdI->box[0][2] != StdI->NaN_i ||
651  StdI->box[1][0] != StdI->NaN_i || StdI->box[1][1] != StdI->NaN_i || StdI->box[1][2] != StdI->NaN_i ||
652  StdI->box[2][0] != StdI->NaN_i || StdI->box[2][1] != StdI->NaN_i || StdI->box[2][2] != StdI->NaN_i)
653  )
654  {
655  fprintf(stdout, "\nERROR ! (L, W, Height) and (a0W, ..., a2H) conflict !\n\n");
656  StdFace_exit(-1);
657  }
658  else if (StdI->L != StdI->NaN_i || StdI->W != StdI->NaN_i || StdI->Height != StdI->NaN_i)
659  {
660  StdFace_PrintVal_i("L", &StdI->L, 1);
661  StdFace_PrintVal_i("W", &StdI->W, 1);
662  StdFace_PrintVal_i("Height", &StdI->Height, 1);
663  for (ii = 0; ii < 3; ii++) for (jj = 0; jj < 3; jj++)
664  StdI->box[ii][jj] = 0;
665  StdI->box[0][0] = StdI->W;
666  StdI->box[1][1] = StdI->L;
667  StdI->box[2][2] = StdI->Height;
668  }
669  else
670  {
671  StdFace_PrintVal_i("a0W", &StdI->box[0][0], 1);
672  StdFace_PrintVal_i("a0L", &StdI->box[0][1], 0);
673  StdFace_PrintVal_i("a0H", &StdI->box[0][2], 0);
674  StdFace_PrintVal_i("a1W", &StdI->box[1][0], 0);
675  StdFace_PrintVal_i("a1L", &StdI->box[1][1], 1);
676  StdFace_PrintVal_i("a1H", &StdI->box[1][2], 0);
677  StdFace_PrintVal_i("a2W", &StdI->box[2][0], 0);
678  StdFace_PrintVal_i("a2L", &StdI->box[2][1], 0);
679  StdFace_PrintVal_i("a2H", &StdI->box[2][2], 1);
680  }
681  /*
682  Parameters for the 3D system will not used.
683  */
684  if (dim == 2) {
685  StdI->direct[0][2] = 0.0;
686  StdI->direct[1][2] = 0.0;
687  StdI->direct[2][0] = 0.0;
688  StdI->direct[2][1] = 0.0;
689  StdI->direct[2][2] = 1.0;
690  }
695  if (dim == 2) StdI->phase[2] = 0.0;
696  for (ii = 0; ii < 3; ii++) {
697  StdI->ExpPhase[ii] = cos(StdI->pi180 * StdI->phase[ii]) + I*sin(StdI->pi180 * StdI->phase[ii]);
698  if (cabs(StdI->ExpPhase[ii] + 1.0) < 0.000001) StdI->AntiPeriod[ii] = 1;
699  else StdI->AntiPeriod[ii] = 0;
700  }
704  StdI->tau = (double **)malloc(sizeof(double*) * StdI->NsiteUC);
705  for (ii = 0; ii < StdI->NsiteUC; ii++) {
706  StdI->tau[ii] = (double *)malloc(sizeof(double) * 3);
707  }
712  StdI->NCell = 0;
713  for (ii = 0; ii < 3; ii++) {
714  StdI->NCell += StdI->box[0][ii]
715  * StdI->box[1][(ii + 1) % 3]
716  * StdI->box[2][(ii + 2) % 3]
717  - StdI->box[0][ii]
718  * StdI->box[1][(ii + 2) % 3]
719  * StdI->box[2][(ii + 1) % 3];
720  }
721  printf(" Number of Cell = %d\n", abs(StdI->NCell));
722  if (StdI->NCell == 0) {
723  StdFace_exit(-1);
724  }
725 
726  for (ii = 0; ii < 3; ii++) {
727  for (jj = 0; jj < 3; jj++) {
728  StdI->rbox[ii][jj] = StdI->box[(ii + 1) % 3][(jj + 1) % 3] * StdI->box[(ii + 2) % 3][(jj + 2) % 3]
729  - StdI->box[(ii + 1) % 3][(jj + 2) % 3] * StdI->box[(ii + 2) % 3][(jj + 1) % 3];
730  }
731  }
732  if (StdI->NCell < 0) {
733  for (ii = 0; ii < 3; ii++)
734  for (jj = 0; jj < 3; jj++)
735  StdI->rbox[ii][jj] *= -1;
736  StdI->NCell *= -1;
737  }/*if (StdI->NCell < 0)*/
742  for (ii = 0; ii < 3; ii++) {
743  bound[ii][0] = 0;
744  bound[ii][1] = 0;
745  for (nBox[2] = 0; nBox[2] < 2; nBox[2]++) {
746  for (nBox[1] = 0; nBox[1] < 2; nBox[1]++) {
747  for (nBox[0] = 0; nBox[0] < 2; nBox[0]++) {
748  edge = 0;
749  for (jj = 0; jj < 3; jj++) edge += nBox[jj] * StdI->box[jj][ii];
750  if (edge < bound[ii][0]) bound[ii][0] = edge;
751  if (edge > bound[ii][1]) bound[ii][1] = edge;
752  }
753  }
754  }
755  }
759  StdI->Cell = (int **)malloc(sizeof(int*) * StdI->NCell);
760  for (ii = 0; ii < StdI->NCell; ii++) {
761  StdI->Cell[ii] = (int *)malloc(sizeof(int) * 3);
762  }/*for (ii = 0; ii < StdI->NCell; ii++)*/
763  jj = 0;
764  for (iCellV[2] = bound[2][0]; iCellV[2] <= bound[2][1]; iCellV[2]++) {
765  for (iCellV[1] = bound[1][0]; iCellV[1] <= bound[1][1]; iCellV[1]++) {
766  for (iCellV[0] = bound[0][0]; iCellV[0] <= bound[0][1]; iCellV[0]++) {
767  StdFace_FoldSite(StdI, iCellV, nBox, iCellV_fold);
768  if (nBox[0] == 0 && nBox[1] == 0 && nBox[2] == 0) {
769  for (ii = 0; ii < 3; ii++)
770  StdI->Cell[jj][ii] = iCellV[ii];
771  jj += 1;
772  }/*if (lUC == 1)*/
773  }/*for (iCellV[0] = bound[0][0]; iCellV[0] <= bound[0][1]; iCellV[0]++*/
774  }/*for (iCellV[1] = bound[1][0]; iCellV[1] <= bound[1][1]; iCellV[1]++)*/
775  }/*for (iCellV[2] = bound[2][0]; iCellV[2] <= bound[2][1]; iCellV[2]++)*/
779  if (dim == 2) {
780  pos[0][0] = 0.0;
781  pos[0][1] = 0.0;
782  pos[1][0] = StdI->direct[0][0] * (double)StdI->box[0][0] + StdI->direct[1][0] * (double)StdI->box[0][1];
783  pos[1][1] = StdI->direct[0][1] * (double)StdI->box[0][0] + StdI->direct[1][1] * (double)StdI->box[0][1];
784  pos[2][0] = StdI->direct[0][0] * (double)StdI->box[1][0] + StdI->direct[1][0] * (double)StdI->box[1][1];
785  pos[2][1] = StdI->direct[0][1] * (double)StdI->box[1][0] + StdI->direct[1][1] * (double)StdI->box[1][1];
786  pos[3][0] = pos[1][0] + pos[2][0];
787  pos[3][1] = pos[1][1] + pos[2][1];
788 
789  xmin = 0.0;
790  xmax = 0.0;
791  for (ipos = 0; ipos < 4; ipos++) {
792  if (pos[ipos][0] < xmin) xmin = pos[ipos][0];
793  if (pos[ipos][0] > xmax) xmax = pos[ipos][0];
794  if (pos[ipos][1] < xmin) xmin = pos[ipos][1];
795  if (pos[ipos][1] > xmax) xmax = pos[ipos][1];
796  }
797  xmin -= 2.0;
798  xmax += 2.0;
799 
800  fprintf(fp, "#set terminal pdf color enhanced \\\n");
801  fprintf(fp, "#dashed dl 1.0 size 20.0cm, 20.0cm \n");
802  fprintf(fp, "#set output \"lattice.pdf\"\n");
803  fprintf(fp, "set xrange [%f: %f]\n", xmin, xmax);
804  fprintf(fp, "set yrange [%f: %f]\n", xmin, xmax);
805  fprintf(fp, "set size square\n");
806  fprintf(fp, "unset key\n");
807  fprintf(fp, "unset tics\n");
808  fprintf(fp, "unset border\n");
809 
810  fprintf(fp, "set style line 1 lc 1 lt 1\n");
811  fprintf(fp, "set style line 2 lc 5 lt 1\n");
812  fprintf(fp, "set style line 3 lc 0 lt 1\n");
813 
814  fprintf(fp, "set arrow from %f, %f to %f, %f nohead front ls 3\n", pos[0][0], pos[0][1], pos[1][0], pos[1][1]);
815  fprintf(fp, "set arrow from %f, %f to %f, %f nohead front ls 3\n", pos[1][0], pos[1][1], pos[3][0], pos[3][1]);
816  fprintf(fp, "set arrow from %f, %f to %f, %f nohead front ls 3\n", pos[3][0], pos[3][1], pos[2][0], pos[2][1]);
817  fprintf(fp, "set arrow from %f, %f to %f, %f nohead front ls 3\n", pos[2][0], pos[2][1], pos[0][0], pos[0][1]);
818  }/*if (dim == 2)*/
819 }/*void StdFace_InitSite2D*/
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...
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 ExpPhase[3]
.
Definition: StdFace_vals.h:134
int L
Number of sites along the 2nd axis, input parameter.
Definition: StdFace_vals.h:40
double pi180
, set in StdFace_ResetVals().
Definition: StdFace_vals.h:132
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
int AntiPeriod[3]
If corresponding StdIntList::phase = 180, it becomes 1.
Definition: StdFace_vals.h:135
int rbox[3][3]
The inversion of StdIntList::box. Set in StdFace_InitSite().
Definition: StdFace_vals.h:47
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:133
double direct[3][3]
The unit direct lattice vector. Set in StdFace_InitSite().
Definition: StdFace_vals.h:42
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
static void StdFace_FoldSite(struct StdIntList *StdI, int iCellV[3], int nBox[3], int iCellV_fold[3])
Move a site into the original supercell if it is outside the original supercell.
int Height
Number of sites along the 3rd axis, input parameter.
Definition: StdFace_vals.h:41
int NaN_i
It is used for initializing input parameter. This means that a parameter wich is not specified in inp...
Definition: StdFace_vals.h:28
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.

◆ StdFace_InputCoulombV()

void StdFace_InputCoulombV ( struct StdIntList StdI,
double *  V0,
char *  V0name 
)

Input off-site Coulomb interaction from the input file, if it is not specified, use the default value (0 or the isotropic Coulomb interaction StdIntList::V).

Parameters
[in,out]StdI
[in]V0
[in]V0nameE.g. V1

Definition at line 1109 of file StdFace_ModelUtil.c.

References StdFace_exit(), StdIntList::V, and StdIntList::V0.

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), and StdFace_Triangular().

1114 {
1115  if (isnan(StdI->V) == 0 && isnan(*V0) == 0) {
1116  fprintf(stdout, "\n ERROR! V and %s conflict !\n\n", V0name);
1117  StdFace_exit(-1);
1118  }
1119  else if (isnan(*V0) == 0)
1120  fprintf(stdout, " %15s = %-10.5f\n", V0name, *V0);
1121  else if (isnan(StdI->V) == 0) {
1122  *V0 = StdI->V;
1123  fprintf(stdout, " %15s = %-10.5f\n", V0name, *V0);
1124  }
1125  else {
1126  *V0 = 0.0;
1127  }
1128 }/*void StdFace_InputCoulombV*/
double V
Off-site Coulomb potential (1st), input parameter.
Definition: StdFace_vals.h:72
double V0
Anisotropic Coulomb potential (1st), input parameter.
Definition: StdFace_vals.h:74
void StdFace_exit(int errorcode)
MPI Abortation wrapper.

◆ StdFace_InputHopp()

void StdFace_InputHopp ( struct StdIntList StdI,
double complex *  t0,
char *  t0name 
)

Input hopping integral from the input file, if it is not specified, use the default value(0 or the isotropic hopping StdIntList::V).

Parameters
[in,out]StdI
[in]t0
[in]t0nameE.g. t1

Definition at line 1134 of file StdFace_ModelUtil.c.

References StdFace_exit(), StdIntList::t, and StdIntList::t0.

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), and StdFace_Triangular().

1139 {
1140  if (isnan(creal(StdI->t)) == 0 && isnan(creal(*t0)) == 0) {
1141  fprintf(stdout, "\n ERROR! t and %s conflict !\n\n", t0name);
1142  StdFace_exit(-1);
1143  }
1144  else if (isnan(creal(*t0)) == 0)
1145  fprintf(stdout, " %15s = %-10.5f %-10.5f\n", t0name, creal(*t0), cimag(*t0));
1146  else if (isnan(creal(StdI->t)) == 0) {
1147  *t0 = StdI->t;
1148  fprintf(stdout, " %15s = %-10.5f %-10.5f\n", t0name, creal(*t0), cimag(*t0));
1149  }
1150  else {
1151  *t0 = 0.0;
1152  }
1153 }/*void StdFace_InputHopp*/
double complex t
Nearest-neighbor hopping, input parameter.
Definition: StdFace_vals.h:62
double complex t0
Anisotropic hopping (1st), input parameter.
Definition: StdFace_vals.h:64
void StdFace_exit(int errorcode)
MPI Abortation wrapper.

◆ StdFace_InputSpin()

void StdFace_InputSpin ( struct StdIntList StdI,
double  Jp[3][3],
double  JpAll,
char *  Jpname 
)

Input spin-spin interaction other than nearest-neighbor.

Parameters
[in,out]StdI
[in]JpFully anisotropic spin interaction
[in]JpAllThe isotropic interaction
JpnameThe name of this spin interaction(e.g.J')

Definition at line 1060 of file StdFace_ModelUtil.c.

References StdIntList::Jp, StdIntList::JpAll, and StdFace_exit().

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), and StdFace_Triangular().

1066 {
1067  int i1, i2;
1068  char Jname[3][3][10];
1069 
1070  strcpy(Jname[0][0], "x\0");
1071  strcpy(Jname[0][1], "xy\0");
1072  strcpy(Jname[0][2], "xz\0");
1073  strcpy(Jname[1][0], "yx\0");
1074  strcpy(Jname[1][1], "y\0");
1075  strcpy(Jname[1][2], "yz\0");
1076  strcpy(Jname[2][0], "zx\0");
1077  strcpy(Jname[2][1], "zy\0");
1078  strcpy(Jname[2][2], "z\0");
1079 
1080  for (i1 = 0; i1 < 3; i1++) {
1081  for (i2 = 0; i2 < 3; i2++) {
1082  if (isnan(JpAll) == 0 && isnan(Jp[i1][i2]) == 0) {
1083  fprintf(stdout, "\n ERROR! %s and %s%s conflict !\n\n", Jpname,
1084  Jpname, Jname[i1][i2]);
1085  StdFace_exit(-1);
1086  }
1087  }/*for (j = 0; j < 3; j++)*/
1088  }/*for (i = 0; i < 3; i++)*/
1089 
1090  for (i1 = 0; i1 < 3; i1++) {
1091  for (i2 = 0; i2 < 3; i2++) {
1092  if (isnan(Jp[i1][i2]) == 0)
1093  fprintf(stdout, " %14s%s = %-10.5f\n", Jpname, Jname[i1][i2], Jp[i1][i2]);
1094  else if (i1 == i2 && isnan(JpAll) == 0) {
1095  Jp[i1][i2] = JpAll;
1096  fprintf(stdout, " %14s%s = %-10.5f\n", Jpname, Jname[i1][i2], Jp[i1][i2]);
1097  }
1098  else {
1099  Jp[i1][i2] = 0.0;
1100  }
1101  }/*for (i2 = 0; i2 < 3; i2++)*/
1102  }/*for (i = 0; i < 3; i++)*/
1103 }/*void StdFace_InputSpin*/
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:102
double JpAll
Isotropic, diagonal spin coupling (2nd Near), input parameter Jp.
Definition: StdFace_vals.h:84
void StdFace_exit(int errorcode)
MPI Abortation wrapper.

◆ StdFace_InputSpinNN()

void StdFace_InputSpinNN ( struct StdIntList StdI,
double  J0[3][3],
double  J0All,
char *  J0name 
)

Input nearest-neighbor spin-spin interaction.

Parameters
[in,out]StdI
[in]J0The anisotropic spin interaction
[in]J0AllThe isotropic interaction
[in]J0nameThe name of this spin interaction (e.g. J1)

Definition at line 973 of file StdFace_ModelUtil.c.

References StdIntList::J, StdIntList::J0, StdIntList::J0All, StdIntList::JAll, and StdFace_exit().

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), and StdFace_Triangular().

979 {
980  int i1, i2, i3, i4;
981  char Jname[3][3][10];
982 
983  strcpy(Jname[0][0], "x\0");
984  strcpy(Jname[0][1], "xy\0");
985  strcpy(Jname[0][2], "xz\0");
986  strcpy(Jname[1][0], "yx\0");
987  strcpy(Jname[1][1], "y\0");
988  strcpy(Jname[1][2], "yz\0");
989  strcpy(Jname[2][0], "zx\0");
990  strcpy(Jname[2][1], "zy\0");
991  strcpy(Jname[2][2], "z\0");
992 
993  if (isnan(StdI->JAll) == 0 && isnan(J0All) == 0) {
994  fprintf(stdout, "\n ERROR! J and %s conflict !\n\n", J0name);
995  StdFace_exit(-1);
996  }
997  for (i1 = 0; i1 < 3; i1++) {
998  for (i2 = 0; i2 < 3; i2++) {
999  if (isnan(StdI->JAll) == 0 && isnan(StdI->J[i1][i2]) == 0) {
1000  fprintf(stdout, "\n ERROR! J and J%s conflict !\n\n", Jname[i1][i2]);
1001  StdFace_exit(-1);
1002  }
1003  else if (isnan(J0All) == 0 && isnan(StdI->J[i1][i2]) == 0) {
1004  fprintf(stdout, "\n ERROR! %s and J%s conflict !\n\n",
1005  J0name, Jname[i1][i2]);
1006  StdFace_exit(-1);
1007  }
1008  else if (isnan(J0All) == 0 && isnan(J0[i1][i2]) == 0) {
1009  fprintf(stdout, "\n ERROR! %s and %s%s conflict !\n\n", J0name,
1010  J0name, Jname[i1][i2]);
1011  StdFace_exit(-1);
1012  }
1013  else if (isnan(J0[i1][i2]) == 0 && isnan(StdI->JAll) == 0) {
1014  fprintf(stdout, "\n ERROR! %s%s and J conflict !\n\n",
1015  J0name, Jname[i1][i2]);
1016  StdFace_exit(-1);
1017  }
1018  }/*for (j = 0; j < 3; j++)*/
1019  }/*for (i = 0; i < 3; i++)*/
1020 
1021  for (i1 = 0; i1 < 3; i1++) {
1022  for (i2 = 0; i2 < 3; i2++) {
1023  for (i3 = 0; i3 < 3; i3++) {
1024  for (i4 = 0; i4 < 3; i4++) {
1025  if (isnan(J0[i1][i2]) == 0 && isnan(StdI->J[i3][i4]) == 0) {
1026  fprintf(stdout, "\n ERROR! %s%s and J%s conflict !\n\n",
1027  J0name, Jname[i1][i2], Jname[i3][i4]);
1028  StdFace_exit(-1);
1029  }
1030  }/*for (i4 = 0; i4 < 3; i4++)*/
1031  }/*for (i3 = 0; i3 < 3; i3++)*/
1032  }/*for (j = 0; j < 3; j++)*/
1033  }/*for (i = 0; i < 3; i++)*/
1034 
1035  for (i1 = 0; i1 < 3; i1++) {
1036  for (i2 = 0; i2 < 3; i2++) {
1037  if (isnan(J0[i1][i2]) == 0)
1038  fprintf(stdout, " %14s%s = %-10.5f\n", J0name, Jname[i1][i2], J0[i1][i2]);
1039  else if (isnan(StdI->J[i1][i2]) == 0) {
1040  J0[i1][i2] = StdI->J[i1][i2];
1041  fprintf(stdout, " %14s%s = %-10.5f\n", J0name, Jname[i1][i2], J0[i1][i2]);
1042  }
1043  else if (i1 == i2 && isnan(J0All) == 0) {
1044  J0[i1][i2] = J0All;
1045  fprintf(stdout, " %14s%s = %-10.5f\n", J0name, Jname[i1][i2], J0[i1][i2]);
1046  }
1047  else if (i1 == i2 && isnan(StdI->JAll) == 0) {
1048  J0[i1][i2] = StdI->JAll;
1049  fprintf(stdout, " %14s%s = %-10.5f\n", J0name, Jname[i1][i2], J0[i1][i2]);
1050  }
1051  else {
1052  J0[i1][i2] = 0.0;
1053  }
1054  }/*for (i2 = 0; i2 < 3; i2++)*/
1055  }/*for (i = 0; i < 3; i++)*/
1056 }/*void StdFace_InputSpinNN*/
double J[3][3]
Isotropic, diagonal/off-diagonal spin coupling (1st Near.), input parameter Jx, Jy, Jz, Jxy, etc.
Definition: StdFace_vals.h:100
double JAll
Isotropic, diagonal spin coupling (1st Near.), input parameter J.
Definition: StdFace_vals.h:82
double J0All
Anisotropic, diagonal spin coupling (1nd Near), input parameter J0.
Definition: StdFace_vals.h:86
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:104
void StdFace_exit(int errorcode)
MPI Abortation wrapper.

◆ StdFace_intr()

void StdFace_intr ( struct StdIntList StdI,
double complex  intr0,
int  site1,
int  spin1,
int  site2,
int  spin2,
int  site3,
int  spin3,
int  site4,
int  spin4 
)

Add interaction (InterAll) to the list Set StdIntList::intr and StdIntList::intrindx and increase the number of that (StdIntList::nintr).

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in,out]StdI
[in]intr0Interaction \(U, V, J\), etc.
[in]site1\(i_1\) for \(c_{i_1 \sigma_1}^\dagger\)
[in]spin1\(sigma1_1\) for \(c_{i_1 \sigma_1}^\dagger\)
[in]site2\(i_2\) for \(c_{i_2 \sigma_2}\)
[in]spin2\(sigma1_2\) for \(c_{i_2 \sigma_2}\)
[in]site3\(i_3\) for \(c_{i_3 \sigma_3}^\dagger\)
[in]spin3\(sigma1_3\) for \(c_{i_3 \sigma_3}^\dagger\)
[in]site4\(i_2\) for \(c_{i_2 \sigma_2}\)
[in]spin4\(sigma1_2\) for \(c_{i_2 \sigma_2}\)

Definition at line 199 of file StdFace_ModelUtil.c.

References StdIntList::intr, StdIntList::intrindx, and StdIntList::nintr.

Referenced by StdFace_GeneralJ().

211 {
212  StdI->intr[StdI->nintr] = intr0;
213  StdI->intrindx[StdI->nintr][0] = site1; StdI->intrindx[StdI->nintr][1] = spin1;
214  StdI->intrindx[StdI->nintr][2] = site2; StdI->intrindx[StdI->nintr][3] = spin2;
215  StdI->intrindx[StdI->nintr][4] = site3; StdI->intrindx[StdI->nintr][5] = spin3;
216  StdI->intrindx[StdI->nintr][6] = site4; StdI->intrindx[StdI->nintr][7] = spin4;
217  StdI->nintr = StdI->nintr + 1;
218 }/*void StdFace_intr*/
double complex * intr
[StdIntList::nintr] Coefficient of general two-body term, malloc in StdFace_MallocInteractions() and ...
Definition: StdFace_vals.h:155
int ** intrindx
[StdIntList::nintr][8] Site/spin indices of two-body term, malloc in StdFace_MallocInteractions() and...
Definition: StdFace_vals.h:152
int nintr
Number of InterAll, counted in each lattice file.
Definition: StdFace_vals.h:150

◆ StdFace_MagField()

void StdFace_MagField ( struct StdIntList StdI,
int  S2,
double  h,
double  Gamma,
int  isite 
)

Add longitudinal and transvars magnetic field to the list.

Author
Mitsuaki Kawamura (The University of Tokyo)

Use Bogoliubov representation.

Londitudinal part

\[ \sum_{\sigma = -S}^{S} -h\sigma c_{i \sigma}^\dagger c_{i \sigma} \]

Transvars part

\[ -\Gamma \frac{S_i^+ + S_i^-}{2} = \sum_{\sigma = -S}^{S-1} -\frac{\Gamma}{2} \sqrt{S(S+1) - \sigma(\sigma+1)} (\sigma c_{i \sigma+ 1}^\dagger c_{i \sigma} + \sigma c_{i \sigma}^\dagger c_{i \sigma+1}) \]

Parameters
[in,out]StdI
[in]S2Spin moment in \(i\) site
[in]hLongitudinal magnetic field \(h\)
[in]GammaTransvars magnetic field \(h\)
[in]isite\(i\) for \(c_{i \sigma}^\dagger\)

Definition at line 151 of file StdFace_ModelUtil.c.

References StdIntList::Gamma, StdIntList::h, StdIntList::S2, and StdFace_trans().

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), StdFace_Triangular(), and StdFace_Wannier90().

158 {
159  int ispin;
160  double S, Sz;
161 
162  S = (double)S2 * 0.5;
166  for (ispin = 0; ispin <= S2; ispin++){
173  Sz = (double)ispin - S;
174  StdFace_trans(StdI, -h * Sz, isite, ispin, isite, ispin);
185  if (ispin < S2) {
186  StdFace_trans(StdI, -0.5 * Gamma * sqrt(S*(S + 1.0) - Sz*(Sz + 1.0)),
187  isite, ispin + 1, isite, ispin);
188  StdFace_trans(StdI, -0.5 * Gamma * sqrt(S*(S + 1.0) - Sz*(Sz + 1.0)),
189  isite, ispin, isite, ispin + 1);
190  }/*if (ispin < S2)*/
191  }/*for (ispin = 0; ispin <= S2; ispin++)*/
192 }/*void StdFace_MagField*/
void StdFace_trans(struct StdIntList *StdI, double complex trans0, int isite, int ispin, int jsite, int jspin)
Add transfer to the list set StdIntList::trans and StdIntList::transindx and increment StdIntList::nt...
int S2
Total spin |S| of a local spin, input from file.
Definition: StdFace_vals.h:215
double Gamma
Transvars magnetic field, input parameter.
Definition: StdFace_vals.h:127
double h
Longitudinal magnetic field, input parameter.
Definition: StdFace_vals.h:126

◆ StdFace_MallocInteractions()

void StdFace_MallocInteractions ( struct StdIntList StdI,
int  ntransMax,
int  nintrMax 
)

Malloc Arrays for interactions.

(1) Transfer StdIntList::trans, StdIntList::transindx

(2) InterAll StdIntList::intr, StdIntList::intrindx

(3) Coulomb intra StdIntList::Cintra, StdIntList::CintraIndx

(4) Coulomb inter StdIntList::Cinter, StdIntList::CinterIndx

(5) Hund StdIntList::Hund, StdIntList::HundIndx

(6) Excahnge StdIntList::Ex, StdIntList::ExIndx

(7) PairLift StdIntList::PairLift, StdIntList::PLIndx

(7) PairHopp StdIntList::PairHopp, StdIntList::PHIndx

Parameters
[in,out]StdI
[in]ntransMaxupper limit of the number of transfer
[in]nintrMaxupper limit of the number of interaction

Definition at line 1198 of file StdFace_ModelUtil.c.

References StdIntList::Cinter, StdIntList::CinterIndx, StdIntList::Cintra, StdIntList::CintraIndx, StdIntList::Ex, StdIntList::ExIndx, StdIntList::Hund, StdIntList::HundIndx, StdIntList::intr, StdIntList::intrindx, StdIntList::Lanczos_max, StdIntList::method, StdIntList::NCinter, StdIntList::NCintra, StdIntList::NEx, StdIntList::NHund, StdIntList::nintr, StdIntList::NPairHopp, StdIntList::NPairLift, StdIntList::npump, StdIntList::ntrans, StdIntList::PairHopp, StdIntList::PairLift, StdIntList::PHIndx, StdIntList::PLIndx, StdIntList::pump, StdIntList::PumpBody, StdIntList::pumpindx, StdIntList::trans, and StdIntList::transindx.

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), StdFace_Triangular(), and StdFace_Wannier90().

1202  {
1203  int ii;
1204 #if defined(_HPhi)
1205  int it;
1206 #endif
1207 
1210  StdI->transindx = (int **)malloc(sizeof(int*) * ntransMax);
1211  StdI->trans = (double complex *)malloc(sizeof(double complex) * ntransMax);
1212  for (ii = 0; ii < ntransMax; ii++) {
1213  StdI->transindx[ii] = (int *)malloc(sizeof(int) * 4);
1214  }
1215  StdI->ntrans = 0;
1216 #if defined(_HPhi)
1217  if (strcmp(StdI->method, "timeevolution") == 0 && StdI->PumpBody == 1) {
1218  StdI->npump = (int *)malloc(sizeof(int) * StdI->Lanczos_max);
1219  StdI->pumpindx = (int ***)malloc(sizeof(int**) * StdI->Lanczos_max);
1220  StdI->pump = (double complex **)malloc(sizeof(double complex*) * StdI->Lanczos_max);
1221  for (it = 0; it < StdI->Lanczos_max; it++) {
1222  StdI->npump[it] = 0;
1223  StdI->pumpindx[it] = (int **)malloc(sizeof(int*) * ntransMax);
1224  StdI->pump[it] = (double complex *)malloc(sizeof(double complex) * ntransMax);
1225  for (ii = 0; ii < ntransMax; ii++) {
1226  StdI->pumpindx[it][ii] = (int *)malloc(sizeof(int) * 4);
1227  }
1228  }/*for (it = 0; it < StdI->Lanczos_max;)*/
1229  }/*if (strcmp(StdI->method, "timeevolution") == 0*/
1230 #endif
1231 
1234  StdI->intrindx = (int **)malloc(sizeof(int*) * nintrMax);
1235  StdI->intr = (double complex *)malloc(sizeof(double complex) * nintrMax);
1236  for (ii = 0; ii < nintrMax; ii++) {
1237  StdI->intrindx[ii] = (int *)malloc(sizeof(int) * 8);
1238  }
1239  StdI->nintr = 0;
1243  StdI->CintraIndx = (int **)malloc(sizeof(int*) * nintrMax);
1244  StdI->Cintra = (double *)malloc(sizeof(double) * nintrMax);
1245  for (ii = 0; ii < nintrMax; ii++) {
1246  StdI->CintraIndx[ii] = (int *)malloc(sizeof(int) * 1);
1247  }
1248  StdI->NCintra = 0;
1252  StdI->CinterIndx = (int **)malloc(sizeof(int*) * nintrMax);
1253  StdI->Cinter = (double *)malloc(sizeof(double) * nintrMax);
1254  for (ii = 0; ii < nintrMax; ii++) {
1255  StdI->CinterIndx[ii] = (int *)malloc(sizeof(int) * 2);
1256  }
1257  StdI->NCinter = 0;
1261  StdI->HundIndx = (int **)malloc(sizeof(int*) * nintrMax);
1262  StdI->Hund = (double *)malloc(sizeof(double) * nintrMax);
1263  for (ii = 0; ii < nintrMax; ii++) {
1264  StdI->HundIndx[ii] = (int *)malloc(sizeof(int) * 2);
1265  }
1266  StdI->NHund = 0;
1270  StdI->ExIndx = (int **)malloc(sizeof(int*) * nintrMax);
1271  StdI->Ex = (double *)malloc(sizeof(double) * nintrMax);
1272  for (ii = 0; ii < nintrMax; ii++) {
1273  StdI->ExIndx[ii] = (int *)malloc(sizeof(int) * 2);
1274  }
1275  StdI->NEx = 0;
1279  StdI->PLIndx = (int **)malloc(sizeof(int*) * nintrMax);
1280  StdI->PairLift = (double *)malloc(sizeof(double) * nintrMax);
1281  for (ii = 0; ii < nintrMax; ii++) {
1282  StdI->PLIndx[ii] = (int *)malloc(sizeof(int) * 2);
1283  }
1284  StdI->NPairLift = 0;
1288  StdI->PHIndx = (int **)malloc(sizeof(int*) * nintrMax);
1289  StdI->PairHopp = (double *)malloc(sizeof(double) * nintrMax);
1290  for (ii = 0; ii < nintrMax; ii++) {
1291  StdI->PHIndx[ii] = (int *)malloc(sizeof(int) * 2);
1292  }
1293  StdI->NPairHopp = 0;
1294 }/*void StdFace_MallocInteractions*/
int NPairHopp
Number of pair-hopping term, counted in each lattice file.
Definition: StdFace_vals.h:200
int NHund
Number of Hund term, counted in each lattice file.
Definition: StdFace_vals.h:176
double complex * intr
[StdIntList::nintr] Coefficient of general two-body term, malloc in StdFace_MallocInteractions() and ...
Definition: StdFace_vals.h:155
double * Ex
[StdIntList::NEx] Coefficient of exchange term, malloc in StdFace_MallocInteractions() and set in Std...
Definition: StdFace_vals.h:189
int ** PLIndx
[StdIntList::NPairLift][2] Site indices of pair-lift term, malloc in StdFace_MallocInteractions() and...
Definition: StdFace_vals.h:194
double complex * trans
[StdIntList::ntrans] Coefficient of one-body term, malloc in StdFace_MallocInteractions() and set in ...
Definition: StdFace_vals.h:147
int ** ExIndx
[StdIntList::NEx][2] Site indices of exchange term, malloc in StdFace_MallocInteractions() and set in...
Definition: StdFace_vals.h:186
int ntrans
Number of transfer, counted in each lattice file.
Definition: StdFace_vals.h:143
int ** CintraIndx
[StdIntList::NCintra][1] Site indices of intra-site Coulomb term, malloc in StdFace_MallocInteraction...
Definition: StdFace_vals.h:161
int NCinter
Number of inter-site Coulomb interaction, counted in each lattice file.
Definition: StdFace_vals.h:167
int NCintra
Number of intra-site Coulomb interaction, counted in each lattice file.
Definition: StdFace_vals.h:158
double * Hund
[StdIntList::NHund] Coefficient of Hund term, malloc in StdFace_MallocInteractions() and set in StdFa...
Definition: StdFace_vals.h:181
double * PairLift
[StdIntList::NPairLift] Coefficient of pair-lift term, malloc in StdFace_MallocInteractions() and set...
Definition: StdFace_vals.h:197
int *** pumpindx
[StdIntList::nt][StdIntList::npump][4] Site/spin indices of one-body term, malloc in StdFace_MallocIn...
Definition: StdFace_vals.h:279
int NPairLift
Number of pair-lift term, counted in each lattice file.
Definition: StdFace_vals.h:192
int PumpBody
one- or two-body pumping, defined from StdIntList::PumpType
Definition: StdFace_vals.h:276
int * npump
[StdIntList::nt] Number of transfer, counted in each lattice file.
Definition: StdFace_vals.h:278
int ** HundIndx
[StdIntList::NHund][2] Site indices of Hund term, malloc in StdFace_MallocInteractions() and set in S...
Definition: StdFace_vals.h:178
int Lanczos_max
The maxixmum number of iterations, input from file.
Definition: StdFace_vals.h:237
double complex ** pump
[StdIntList::nt][StdIntList::npump] Coefficient of one-body term, malloc in StdFace_MallocInteraction...
Definition: StdFace_vals.h:282
int NEx
Number of exchange term, counted in each lattice file.
Definition: StdFace_vals.h:184
char method[256]
The name of method, input from file.
Definition: StdFace_vals.h:232
int ** intrindx
[StdIntList::nintr][8] Site/spin indices of two-body term, malloc in StdFace_MallocInteractions() and...
Definition: StdFace_vals.h:152
double * PairHopp
[StdIntList::NPairLift] Coefficient of pair-hopping term, malloc in StdFace_MallocInteractions() and ...
Definition: StdFace_vals.h:205
int ** transindx
[StdIntList::ntrans][4] Site/spin indices of one-body term, malloc in StdFace_MallocInteractions() an...
Definition: StdFace_vals.h:144
double * Cintra
[StdIntList::NCintra] Coefficient of intra-site Coulomb term, malloc in StdFace_MallocInteractions() ...
Definition: StdFace_vals.h:164
int nintr
Number of InterAll, counted in each lattice file.
Definition: StdFace_vals.h:150
int ** PHIndx
[StdIntList::NPairLift][2] Site indices of pair-hopping term, malloc in StdFace_MallocInteractions() ...
Definition: StdFace_vals.h:202
double * Cinter
[StdIntList::NCinter] Coefficient of inter-site Coulomb term, malloc in StdFace_MallocInteractions() ...
Definition: StdFace_vals.h:173
int ** CinterIndx
[StdIntList::NCinter][2] Site indices of inter-site Coulomb term, malloc in StdFace_MallocInteraction...
Definition: StdFace_vals.h:170

◆ StdFace_NotUsed_c()

void StdFace_NotUsed_c ( char *  valname,
double complex  val 
)

Stop HPhi if a variable (complex) not used is specified in the input file (!=NaN).

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]valnameName of the valiable
[in]val

Definition at line 509 of file StdFace_ModelUtil.c.

References StdFace_exit().

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), and StdFace_Triangular().

513 {
514  if (isnan(creal(val)) == 0) {
515  fprintf(stdout, "\n Check ! %s is SPECIFIED but will NOT be USED. \n", valname);
516  fprintf(stdout, " Please COMMENT-OUT this line \n");
517  fprintf(stdout, " or check this input is REALLY APPROPRIATE for your purpose ! \n\n");
518  StdFace_exit(-1);
519  }
520 }/*void StdFace_NotUsed_c*/
void StdFace_exit(int errorcode)
MPI Abortation wrapper.

◆ StdFace_NotUsed_d()

void StdFace_NotUsed_d ( char *  valname,
double  val 
)

Stop HPhi if a variable (real) not used is specified in the input file (!=NaN).

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]valnameName of the valiable
[in]val

Definition at line 492 of file StdFace_ModelUtil.c.

References StdFace_exit().

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_NotUsed_J(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), StdFace_Triangular(), and StdFace_Wannier90().

496 {
497  if (isnan(val) == 0) {
498  fprintf(stdout, "\n Check ! %s is SPECIFIED but will NOT be USED. \n", valname);
499  fprintf(stdout, " Please COMMENT-OUT this line \n");
500  fprintf(stdout, " or check this input is REALLY APPROPRIATE for your purpose ! \n\n");
501  StdFace_exit(-1);
502  }
503 }/*void StdFace_NotUsed_d*/
void StdFace_exit(int errorcode)
MPI Abortation wrapper.

◆ StdFace_NotUsed_i()

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).

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]valnameName of the valiable
[in]val

Definition at line 558 of file StdFace_ModelUtil.c.

References StdIntList::NaN_i, and StdFace_exit().

Referenced by CheckModPara(), StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), and StdFace_Triangular().

562 {
563  int NaN_i = 2147483647;
564 
565  if (val != NaN_i) {
566  fprintf(stdout, "\n Check ! %s is SPECIFIED but will NOT be USED. \n", valname);
567  fprintf(stdout, " Please COMMENT-OUT this line \n");
568  fprintf(stdout, " or check this input is REALLY APPROPRIATE for your purpose ! \n\n");
569  StdFace_exit(-1);
570  }
571 }/*void StdFace_NotUsed_i*/
int NaN_i
It is used for initializing input parameter. This means that a parameter wich is not specified in inp...
Definition: StdFace_vals.h:28
void StdFace_exit(int errorcode)
MPI Abortation wrapper.

◆ StdFace_NotUsed_J()

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).

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]valnameName of the valiable*/,
[in]JAll*/,
[in]J

Definition at line 526 of file StdFace_ModelUtil.c.

References StdIntList::J, StdIntList::JAll, and StdFace_NotUsed_d().

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), and StdFace_Triangular().

531 {
532  int i1, i2;
533  char Jname[3][3][10];
534 
535  sprintf(Jname[0][0], "%sx", valname);
536  sprintf(Jname[0][1], "%sxy", valname);
537  sprintf(Jname[0][2], "%sxz", valname);
538  sprintf(Jname[1][0], "%syx", valname);
539  sprintf(Jname[1][1], "%sy", valname);
540  sprintf(Jname[1][2], "%syz", valname);
541  sprintf(Jname[2][0], "%szx", valname);
542  sprintf(Jname[2][1], "%szy", valname);
543  sprintf(Jname[2][2], "%sz", valname);
544 
545  StdFace_NotUsed_d(valname, JAll);
546 
547  for (i1 = 0; i1 < 3; i1++) {
548  for (i2 = 0; i2 < 3; i2++) {
549  StdFace_NotUsed_d(Jname[i1][i2], J[i1][i2]);
550  }/*for (j = 0; j < 3; j++)*/
551  }/*for (i = 0; i < 3; i++)*/
552 }/*void StdFace_NotUsed_J*/
double J[3][3]
Isotropic, diagonal/off-diagonal spin coupling (1st Near.), input parameter Jx, Jy, Jz, Jxy, etc.
Definition: StdFace_vals.h:100
double JAll
Isotropic, diagonal spin coupling (1st Near.), input parameter J.
Definition: StdFace_vals.h:82
void StdFace_NotUsed_d(char *valname, double val)
Stop HPhi if a variable (real) not used is specified in the input file (!=NaN).

◆ StdFace_PrintGeometry()

void StdFace_PrintGeometry ( struct StdIntList StdI)

Print geometry of sites for the pos-process of correlation function.

Definition at line 1157 of file StdFace_ModelUtil.c.

References StdIntList::box, StdIntList::Cell, StdIntList::direct, StdIntList::model, StdIntList::NCell, StdIntList::NsiteUC, and StdIntList::phase.

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), StdFace_Triangular(), and StdFace_Wannier90().

1157  {
1158  FILE *fp;
1159  int isite, iCell, ii;
1160 
1161  fp = fopen("geometry.dat", "w");
1162 
1163  for (ii = 0; ii < 3; ii++)
1164  fprintf(fp, "%25.15e %25.15e %25.15e\n",
1165  StdI->direct[ii][0], StdI->direct[ii][1], StdI->direct[ii][2]);
1166  fprintf(fp, "%25.15e %25.15e %25.15e\n",
1167  StdI->phase[0], StdI->phase[1], StdI->phase[2]);
1168  for (ii = 0; ii < 3; ii++)
1169  fprintf(fp, "%d %d %d\n",
1170  StdI->box[ii][0], StdI->box[ii][1], StdI->box[ii][2]);
1171 
1172  for (iCell = 0; iCell < StdI->NCell; iCell++) {
1173  for (isite = 0; isite < StdI->NsiteUC; isite++) {
1174  fprintf(fp, "%d %d %d %d\n",
1175  StdI->Cell[iCell][0] - StdI->Cell[0][0],
1176  StdI->Cell[iCell][1] - StdI->Cell[0][1],
1177  StdI->Cell[iCell][2] - StdI->Cell[0][2],
1178  isite);
1179  }/*for (isite = 0; isite < StdI->NsiteUC; isite++)*/
1180  }/* for (iCell = 0; iCell < StdI->NCell; iCell++)*/
1181  if (strcmp(StdI->model, "kondo") == 0) {
1182  for (iCell = 0; iCell < StdI->NCell; iCell++) {
1183  for (isite = 0; isite < StdI->NsiteUC; isite++) {
1184  fprintf(fp, "%d %d %d %d\n",
1185  StdI->Cell[iCell][0] - StdI->Cell[0][0],
1186  StdI->Cell[iCell][1] - StdI->Cell[0][1],
1187  StdI->Cell[iCell][2] - StdI->Cell[0][2],
1188  isite + StdI->NsiteUC);
1189  }/*for (isite = 0; isite < StdI->NsiteUC; isite++)*/
1190  }/* for (iCell = 0; iCell < StdI->NCell; iCell++)*/
1191  }
1192  fflush(fp);
1193  fclose(fp);
1194 }/*void StdFace_PrintGeometry()*/
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
int NsiteUC
Number of sites in the unit cell. Defined in the beginning of each lattice function.
Definition: StdFace_vals.h:53
char model[256]
Name of model, input parameter.
Definition: StdFace_vals.h:60
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:133
double direct[3][3]
The unit direct lattice vector. Set in StdFace_InitSite().
Definition: StdFace_vals.h:42
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

◆ StdFace_PrintVal_c()

void StdFace_PrintVal_c ( char *  valname,
double complex *  val,
double complex  val0 
)

Print a valiable (complex) read from the input file if it is not specified in the input file (=NaN), set the default value.

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]valnameName of the valiable
[in,out]valValiable to be set
[in]val0The default value

Definition at line 455 of file StdFace_ModelUtil.c.

Referenced by StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Orthorhombic(), StdFace_Tetragonal(), and StdFace_Triangular().

460 {
461  if (isnan(creal(*val)) == 1) {
462  *val = val0;
463  fprintf(stdout, " %15s = %-10.5f %-10.5f ###### DEFAULT VALUE IS USED ######\n", valname, creal(*val), cimag(*val));
464  }
465  else fprintf(stdout, " %15s = %-10.5f %-10.5f\n", valname, creal(*val), cimag(*val));
466 }/*void StdFace_PrintVal_c*/

◆ StdFace_PrintVal_d()

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), set the default value.

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]valnameName of the valiable
[in,out]valValiable to be set
[in]val0The default value

Definition at line 414 of file StdFace_ModelUtil.c.

Referenced by CheckModPara(), PrintExcitation(), StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_LargeValue(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), StdFace_Triangular(), StdFace_Wannier90(), and VectorPotential().

419 {
420  if (isnan(*val) == 1) {
421  *val = val0;
422  fprintf(stdout, " %15s = %-10.5f ###### DEFAULT VALUE IS USED ######\n", valname, *val);
423  }
424  else fprintf(stdout, " %15s = %-10.5f\n", valname, *val);
425 }/*void StdFace_PrintVal_d*/

◆ StdFace_PrintVal_dd()

void StdFace_PrintVal_dd ( char *  valname,
double *  val,
double  val0,
double  val1 
)

Print a valiable (real) read from the input file if it is not specified in the input file (=NaN), set the default value.

Author
Mitsuaki Kawamura (The University of Tokyo)

If the primary default value (val0) is not specified, use secondary default value

Parameters
[in]valnameName of the valiable
[in,out]valValiable to be set
[in]val0The primary default value, possible not to be specified
[in]val1The secondary default value

Definition at line 432 of file StdFace_ModelUtil.c.

438 {
439  if (isnan(*val) == 1) {
443  if (isnan(val0) == 1) *val = val1;
444  else *val = val0;
445  fprintf(stdout, " %15s = %-10.5f ###### DEFAULT VALUE IS USED ######\n", valname, *val);
446  }
447  else fprintf(stdout, " %15s = %-10.5f\n", valname, *val);
448 }/*void StdFace_PrintVal_dd*/

◆ StdFace_PrintVal_i()

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 (=2147483647, the upper limt of Int) set the default value.

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]valnameName of the valiable
[in,out]valValiable to be set
[in]val0The default value

Definition at line 473 of file StdFace_ModelUtil.c.

References StdIntList::NaN_i.

Referenced by CheckModPara(), StdFace_Chain(), StdFace_FCOrtho(), StdFace_Honeycomb(), StdFace_InitSite(), StdFace_Kagome(), StdFace_Ladder(), StdFace_main(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), StdFace_Tetragonal(), StdFace_Triangular(), StdFace_Wannier90(), and VectorPotential().

478 {
479  int NaN_i = 2147483647;/*The upper limt of Int*/
480 
481  if (*val == NaN_i) {
482  *val = val0;
483  fprintf(stdout, " %15s = %-10d ###### DEFAULT VALUE IS USED ######\n", valname, *val);
484  }
485  else fprintf(stdout, " %15s = %-10d\n", valname, *val);
486 }/*void StdFace_PrintVal_i*/
int NaN_i
It is used for initializing input parameter. This means that a parameter wich is not specified in inp...
Definition: StdFace_vals.h:28

◆ StdFace_PrintXSF()

void StdFace_PrintXSF ( struct StdIntList StdI)

Print lattice.xsf (XCrysDen format)

Definition at line 938 of file StdFace_ModelUtil.c.

References StdIntList::box, StdIntList::Cell, StdIntList::direct, StdIntList::NCell, StdIntList::NsiteUC, StdIntList::tau, and vec.

Referenced by StdFace_FCOrtho(), StdFace_Orthorhombic(), StdFace_Pyrochlore(), and StdFace_Wannier90().

938  {
939  FILE *fp;
940  int ii, jj, kk, isite, iCell;
941  double vec[3];
942 
943  fp = fopen("lattice.xsf", "w");
944  fprintf(fp, "CRYSTAL\n");
945  fprintf(fp, "PRIMVEC\n");
946  for (ii = 0; ii < 3; ii++) {
947  for (jj = 0; jj < 3; jj++) {
948  vec[jj] = 0.0;
949  for (kk = 0; kk < 3; kk++)
950  vec[jj] += (double)StdI->box[ii][kk] * StdI->direct[kk][jj];
951  }
952  fprintf(fp, "%15.5f %15.5f %15.5f\n", vec[0], vec[1], vec[2]);
953  }
954  fprintf(fp, "PRIMCOORD\n");
955  fprintf(fp, "%d 1\n", StdI->NCell * StdI->NsiteUC);
956  for (iCell = 0; iCell < StdI->NCell; iCell++) {
957  for (isite = 0; isite < StdI->NsiteUC; isite++) {
958  for (jj = 0; jj < 3; jj++) {
959  vec[jj] = 0.0;
960  for (kk = 0; kk < 3; kk++)
961  vec[jj] += ((double)StdI->Cell[iCell][kk] + StdI->tau[isite][kk])
962  * StdI->direct[kk][jj];
963  }
964  fprintf(fp, "H %15.5f %15.5f %15.5f\n", vec[0], vec[1], vec[2]);
965  }
966  }
967  fflush(fp);
968  fclose(fp);
969 }/*void StdFace_PrintXSF*/
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 ** vec
Definition: global.h:45
int NsiteUC
Number of sites in the unit cell. Defined in the beginning of each lattice function.
Definition: StdFace_vals.h:53
int ** Cell
[StdIntList][3] The cell position in the fractional coordinate. Malloc and Set in StdFace_InitSite()...
Definition: StdFace_vals.h:51
double direct[3][3]
The unit direct lattice vector. Set in StdFace_InitSite().
Definition: StdFace_vals.h:42
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 ** tau
Cell-internal site position in the fractional coordinate. Defined in the beginning of each lattice fu...
Definition: StdFace_vals.h:55

◆ StdFace_RequiredVal_i()

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, the upper limt of Int).

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]valnameName of the valiable
[in]val

Definition at line 577 of file StdFace_ModelUtil.c.

References StdIntList::NaN_i, and StdFace_exit().

Referenced by CheckModPara(), StdFace_Chain(), and StdFace_Ladder().

581 {
582  int NaN_i = 2147483647;
583 
584  if (val == NaN_i){
585  fprintf(stdout, "ERROR ! %s is NOT specified !\n", valname);
586  StdFace_exit(-1);
587  }
588  else fprintf(stdout, " %15s = %-3d\n", valname, val);
589 }/*void StdFace_RequiredVal_i*/
int NaN_i
It is used for initializing input parameter. This means that a parameter wich is not specified in inp...
Definition: StdFace_vals.h:28
void StdFace_exit(int errorcode)
MPI Abortation wrapper.

◆ StdFace_SetLabel()

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)

First print the reversed one

Then print the normal one, these are different when they cross boundary.

Parameters
[in,out]StdI
[in]fpFile pointer to lattice.gp
[in]iWposition of initial site
[in]iLposition of initial site
[in]diWTranslation from the initial site
[in]diLTranslation from the initial site
[in]isiteUCIntrinsic site index of initial site
[in]jsiteUCIntrinsic site index of final site
[out]isiteinitial site
[out]jsitefinal site
[in]connect1 for nearest neighbor, 2 for 2nd nearest
[out]CphaseBoundary phase, if it across boundary
[out]dRR_i - R_j

Definition at line 877 of file StdFace_ModelUtil.c.

References StdIntList::direct, StdFace_FindSite(), and StdIntList::tau.

Referenced by StdFace_Chain(), StdFace_Honeycomb(), StdFace_Kagome(), StdFace_Ladder(), StdFace_Tetragonal(), and StdFace_Triangular().

892 {
893  double xi, yi, xj, yj;
897  StdFace_FindSite(StdI, iW, iL, 0, -diW, -diL, 0, jsiteUC, isiteUC, isite, jsite, Cphase, dR);
898 
899  xi = StdI->direct[0][0] * ((double)iW + StdI->tau[jsiteUC][0])
900  + StdI->direct[1][0] * ((double)iL + StdI->tau[jsiteUC][1]);
901  yi = StdI->direct[0][1] * ((double)iW + StdI->tau[jsiteUC][0])
902  + StdI->direct[1][1] * ((double)iL + StdI->tau[jsiteUC][1]);
903 
904  xj = StdI->direct[0][0] * ((double)(iW - diW) + StdI->tau[isiteUC][0])
905  + StdI->direct[1][0] * ((double)(iL - diL) + StdI->tau[isiteUC][1]);
906  yj = StdI->direct[0][1] * ((double)(iW - diW) + StdI->tau[isiteUC][0])
907  + StdI->direct[1][1] * ((double)(iL - diL) + StdI->tau[isiteUC][1]);
908 
909  if (*isite < 10)fprintf(fp, "set label \"%1d\" at %f, %f center front\n", *isite, xi, yi);
910  else fprintf(fp, "set label \"%2d\" at %f, %f center front\n", *isite, xi, yi);
911  if (*jsite < 10)fprintf(fp, "set label \"%1d\" at %f, %f center front\n", *jsite, xj, yj);
912  else fprintf(fp, "set label \"%2d\" at %f, %f center front\n", *jsite, xj, yj);
913  fprintf(fp, "set arrow from %f, %f to %f, %f nohead ls %d\n", xi, yi, xj, yj, connect);
917  StdFace_FindSite(StdI, iW, iL, 0, diW, diL, 0, isiteUC, jsiteUC, isite, jsite, Cphase, dR);
918 
919  xi = StdI->direct[1][0] * ((double)iL + StdI->tau[isiteUC][1])
920  + StdI->direct[0][0] * ((double)iW + StdI->tau[isiteUC][0]);
921  yi = StdI->direct[1][1] * ((double)iL + StdI->tau[isiteUC][1])
922  + StdI->direct[0][1] * ((double)iW + StdI->tau[isiteUC][0]);
923 
924  xj = StdI->direct[0][0] * ((double)(iW + diW) + StdI->tau[jsiteUC][0])
925  + StdI->direct[1][0] * ((double)(iL + diL) + StdI->tau[jsiteUC][1]);
926  yj = StdI->direct[0][1] * ((double)(iW + diW) + StdI->tau[jsiteUC][0])
927  + StdI->direct[1][1] * ((double)(iL + diL) + StdI->tau[jsiteUC][1]);
928 
929  if (*isite < 10)fprintf(fp, "set label \"%1d\" at %f, %f center front\n", *isite, xi, yi);
930  else fprintf(fp, "set label \"%2d\" at %f, %f center front\n", *isite, xi, yi);
931  if (*jsite < 10)fprintf(fp, "set label \"%1d\" at %f, %f center front\n", *jsite, xj, yj);
932  else fprintf(fp, "set label \"%2d\" at %f, %f center front\n", *jsite, xj, yj);
933  fprintf(fp, "set arrow from %f, %f to %f, %f nohead ls %d\n", xi, yi, xj, yj, connect);
934 }/*void StdFace_SetLabel*/
double direct[3][3]
The unit direct lattice vector. Set in StdFace_InitSite().
Definition: StdFace_vals.h:42
void StdFace_FindSite(struct StdIntList *StdI, int iW, int iL, int iH, int diW, int diL, int diH, int isiteUC, int jsiteUC, int *isite, int *jsite, double complex *Cphase, double *dR)
Find the index of transfer and interaction.
double ** tau
Cell-internal site position in the fractional coordinate. Defined in the beginning of each lattice fu...
Definition: StdFace_vals.h:55

◆ StdFace_trans()

void StdFace_trans ( struct StdIntList StdI,
double complex  trans0,
int  isite,
int  ispin,
int  jsite,
int  jspin 
)

Add transfer to the list set StdIntList::trans and StdIntList::transindx and increment StdIntList::ntrans.

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in,out]StdI
[in]trans0Hopping integral \(t, mu\), etc.
[in]isite\(i\) for \(c_{i \sigma}^\dagger\)
[in]ispin\(\sigma\) for \(c_{i \sigma}^\dagger\)
[in]jsite\(j\) for \(c_{j \sigma'}\)
[in]jspin\(\sigma'\) for \(c_{j \sigma'}\)

Definition at line 55 of file StdFace_ModelUtil.c.

References StdIntList::ntrans, StdIntList::trans, and StdIntList::transindx.

Referenced by StdFace_Hopping(), StdFace_HubbardLocal(), and StdFace_MagField().

63 {
64  StdI->trans[StdI->ntrans] = trans0;
65  StdI->transindx[StdI->ntrans][0] = isite;
66  StdI->transindx[StdI->ntrans][1] = ispin;
67  StdI->transindx[StdI->ntrans][2] = jsite;
68  StdI->transindx[StdI->ntrans][3] = jspin;
69  StdI->ntrans = StdI->ntrans + 1;
70 }/*void StdFace_trans*/
double complex * trans
[StdIntList::ntrans] Coefficient of one-body term, malloc in StdFace_MallocInteractions() and set in ...
Definition: StdFace_vals.h:147
int ntrans
Number of transfer, counted in each lattice file.
Definition: StdFace_vals.h:143
int ** transindx
[StdIntList::ntrans][4] Site/spin indices of one-body term, malloc in StdFace_MallocInteractions() an...
Definition: StdFace_vals.h:144