HPhi++  3.1.0
StdFace_ModelUtil.cpp File Reference

Various utility for constructing models. More...

#include <cstdio>
#include <cstdlib>
#include <cmath>
#include <complex>
#include "StdFace_vals.hpp"
#include <cstring>
#include <mpi.h>
#include <iostream>

Go to the source code of this file.

Functions

void StdFace_exit (int errorcode)
 MPI Abortation wrapper. More...
 
void StdFace_trans (struct StdIntList *StdI, std::complex< double > 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, std::complex< double > 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, std::complex< double > 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 (const char *valname, double *val, double val0)
 Print a valiable (real) read from the input file if it is not specified in the input file (=NaN), 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 (const char *valname, std::complex< double > *val, std::complex< double > 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 (const char *valname, int *val, int val0)
 Print a valiable (integer) read from the input file if it is not specified in the input file (=2147483647, the upper limt of Int) set the default value. More...
 
void StdFace_NotUsed_d (const char *valname, double val)
 Stop HPhi if a variable (real) not used is specified in the input file (!=NaN). More...
 
void StdFace_NotUsed_c (const char *valname, std::complex< double > val)
 Stop HPhi if a variable (complex) not used is specified in the input file (!=NaN). More...
 
void StdFace_NotUsed_J (const 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 (const char *valname, int val)
 Stop HPhi if a variable (integer) not used is specified in the input file (!=2147483647, the upper limt of Int). More...
 
void StdFace_RequiredVal_i (const char *valname, int val)
 Stop HPhi if a variable (integer) which must be specified is absent in the input file (=2147483647, 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, std::complex< double > *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, std::complex< double > *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 (double J[3][3], double JAll, double J0[3][3], double J0All, const char *J0name)
 Input nearest-neighbor spin-spin interaction. More...
 
void StdFace_InputSpin (double Jp[3][3], double JpAll, const char *Jpname)
 Input spin-spin interaction other than nearest-neighbor. More...
 
void StdFace_InputCoulombV (double V, double *V0, const 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 (std::complex< double > t, std::complex< double > *t0, const 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.cpp.

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 400 of file StdFace_ModelUtil.cpp.

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

406 {
407  StdI->Cinter[StdI->NCinter] = V;
408  StdI->CinterIndx[StdI->NCinter][0] = isite;
409  StdI->CinterIndx[StdI->NCinter][1] = jsite;
410  StdI->NCinter += 1;
411 }/*void StdFace_Coulomb*/

◆ StdFace_exit()

void StdFace_exit ( int  errorcode)

MPI Abortation wrapper.

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

Definition at line 36 of file StdFace_ModelUtil.cpp.

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_MallocInteractions(), 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().

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

◆ 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,
std::complex< double > *  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 827 of file StdFace_ModelUtil.cpp.

References StdFace_FoldSite().

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

842 {
843  int iCell, jCell, kCell, ii;
844  int nBox[3], jCellV[3];
845 
846  dR[0] = - (double)diW + StdI->tau[isiteUC][0] - StdI->tau[jsiteUC][0];
847  dR[1] = - (double)diL + StdI->tau[isiteUC][1] - StdI->tau[jsiteUC][1];
848  dR[2] = - (double)diH + StdI->tau[isiteUC][2] - StdI->tau[jsiteUC][2];
849 
850  jCellV[0] = iW + diW;
851  jCellV[1] = iL + diL;
852  jCellV[2] = iH + diH;
853  StdFace_FoldSite(StdI, jCellV, nBox, jCellV);
854  *Cphase = 1.0;
855  for (ii = 0; ii < 3; ii++) *Cphase *= pow(StdI->ExpPhase[ii], (double)nBox[ii]);
856 
857  for (kCell = 0; kCell < StdI->NCell; kCell++) {
858  if (jCellV[0] == StdI->Cell[kCell][0] &&
859  jCellV[1] == StdI->Cell[kCell][1] &&
860  jCellV[2] == StdI->Cell[kCell][2])
861  {
862  jCell = kCell;
863  }
864  if (iW == StdI->Cell[kCell][0] &&
865  iL == StdI->Cell[kCell][1] &&
866  iH == StdI->Cell[kCell][2])
867  {
868  iCell = kCell;
869  }
870  }/*for (iCell = 0; iCell < StdI->NCell; iCell++)*/
871  *isite = iCell * StdI->NsiteUC + isiteUC;
872  *jsite = jCell * StdI->NsiteUC + jsiteUC;
873  if (strcmp(StdI->model, "kondo") == 0) {
874  *isite += StdI->NCell * StdI->NsiteUC;
875  *jsite += StdI->NCell * StdI->NsiteUC;
876  }
877 }/*void StdFace_FindSite*/
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.

◆ 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 599 of file StdFace_ModelUtil.cpp.

Referenced by StdFace_FindSite(), StdFace_InitSite(), and StdFace_MallocInteractions().

606 {
607  int ii, jj, iCellV_frac[3];
611  for (ii = 0; ii < 3; ii++) {
612  iCellV_frac[ii] = 0;
613  for (jj = 0; jj < 3; jj++)iCellV_frac[ii] += StdI->rbox[ii][jj] * iCellV[jj];
614  }
618  for (ii = 0; ii < 3; ii++)
619  nBox[ii] = (iCellV_frac[ii] + StdI->NCell * 1000) / StdI->NCell - 1000;
623  for (ii = 0; ii < 3; ii++)
624  iCellV_frac[ii] -= StdI->NCell*(nBox[ii]);
625 
626  for (ii = 0; ii < 3; ii++) {
627  iCellV_fold[ii] = 0;
628  for (jj = 0; jj < 3; jj++) iCellV_fold[ii] += StdI->box[jj][ii] * iCellV_frac[jj];
629  iCellV_fold[ii] = (iCellV_fold[ii] + StdI->NCell * 1000) / StdI->NCell - 1000;
630  }/*for (ii = 0; ii < 3; ii++)*/
631 }/*static void StdFace_FoldSite*/

◆ 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 227 of file StdFace_ModelUtil.cpp.

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

235 {
236  int ispin, jspin, ZGeneral, ExGeneral;
237  double Si, Sj, Siz, Sjz;
238  std::complex<double> intr0;
242  ZGeneral = 1;
243  ExGeneral = 1;
244  if (Si2 == 1 || Sj2 == 1) {
251  ZGeneral = 0;
252 
253  StdI->Hund[StdI->NHund] = -0.5 * J[2][2];
254  StdI->HundIndx[StdI->NHund][0] = isite;
255  StdI->HundIndx[StdI->NHund][1] = jsite;
256  StdI->NHund += 1;
257 
258  StdI->Cinter[StdI->NCinter] = -0.25 * J[2][2];
259  StdI->CinterIndx[StdI->NCinter][0] = isite;
260  StdI->CinterIndx[StdI->NCinter][1] = jsite;
261  StdI->NCinter += 1;
262 
263  if (fabs(J[0][1]) < 0.000001 && fabs(J[1][0]) < 0.000001
264 #if defined(_mVMC)
265  && abs(J[0][0] - J[1][1]) < 0.000001
266 #endif
267  ) {
274  ExGeneral = 0;
275 
276 #if defined(_mVMC)
277  StdI->Ex[StdI->NEx] = - 0.25 * (J[0][0] + J[1][1]);
278 #else
279  if (strcmp(StdI->model, "kondo") == 0)
280  StdI->Ex[StdI->NEx] = -0.25 * (J[0][0] + J[1][1]);
281  else
282  StdI->Ex[StdI->NEx] = 0.25 * (J[0][0] + J[1][1]);
283 #endif
284  StdI->ExIndx[StdI->NEx][0] = isite;
285  StdI->ExIndx[StdI->NEx][1] = jsite;
286  StdI->NEx += 1;
287 
288  StdI->PairLift[StdI->NPairLift] = 0.25 * (J[0][0] - J[1][1]);
289  StdI->PLIndx[StdI->NPairLift][0] = isite;
290  StdI->PLIndx[StdI->NPairLift][1] = jsite;
291  StdI->NPairLift += 1;
292  }/*if (fabs(J[0][1]) < 0.000001 && fabs(J[1][0]) < 0.000001)*/
293  }
298  Si = 0.5 * (double)Si2;
299  Sj = 0.5 * (double)Sj2;
300 
301  for (ispin = 0; ispin <= Si2; ispin++) {
302  Siz = (double)ispin - Si;
303  for (jspin = 0; jspin <= Sj2; jspin++) {
304  Sjz = (double)jspin - Sj;
312  if (ZGeneral == 1) {
313  intr0 = J[2][2] * Siz * Sjz;
314  StdFace_intr(StdI, intr0,
315  isite, ispin, isite, ispin, jsite, jspin, jsite, jspin);
316  }
328  if ((ispin < Si2 && jspin < Sj2) && ExGeneral == 1) {
329  intr0 = 0.25 * std::complex<double>(J[0][0] + J[1][1], J[0][1] - J[1][0])
330  * sqrt(Si * (Si + 1.0) - Siz * (Siz + 1.0))
331  * sqrt(Sj * (Sj + 1.0) - Sjz * (Sjz + 1.0));
332  StdFace_intr(StdI, intr0,
333  isite, ispin + 1, isite, ispin, jsite, jspin, jsite, jspin + 1);
334  StdFace_intr(StdI, conj(intr0),
335  isite, ispin, isite, ispin + 1, jsite, jspin + 1, jsite, jspin);
336  }
348  if ((ispin < Si2 && jspin < Sj2) && ExGeneral == 1) {
349  intr0 = 0.5 * 0.5 * std::complex<double>(J[0][0] - J[1][1], - (J[0][1] + J[1][0]))
350  * sqrt(Si * (Si + 1.0) - Siz * (Siz + 1.0))
351  * sqrt(Sj * (Sj + 1.0) - Sjz * (Sjz + 1.0));
352  StdFace_intr(StdI, intr0,
353  isite, ispin + 1, isite, ispin, jsite, jspin + 1, jsite, jspin);
354  StdFace_intr(StdI, conj(intr0),
355  isite, ispin, isite, ispin + 1, jsite, jspin, jsite, jspin + 1);
356  }
367  if (ispin < Si2) {
368  intr0 = 0.5 * std::complex<double>(J[0][2], - J[1][2]) * sqrt(Si * (Si + 1.0) - Siz * (Siz + 1.0)) * Sjz;
369  StdFace_intr(StdI, intr0,
370  isite, ispin + 1, isite, ispin, jsite, jspin, jsite, jspin);
371  StdFace_intr(StdI, conj(intr0),
372  jsite, jspin, jsite, jspin, isite, ispin, isite, ispin + 1);
373  }/*if (ispin < Si2)*/
384  if (jspin < Sj2) {
385  intr0 = 0.5 * std::complex<double>(J[2][0], - J[2][1]) * Siz * sqrt(Sj * (Sj + 1.0) - Sjz * (Sjz + 1.0));
386  StdFace_intr(StdI, intr0,
387  isite, ispin, isite, ispin, jsite, jspin + 1, jsite, jspin);
388  StdFace_intr(StdI, conj(intr0),
389  jsite, jspin, jsite, jspin + 1, isite, ispin, isite, ispin);
390  }/*if (jspin < Sj2)*/
391  }/*for (jspin = 0; jspin <= Sj2; jspin++)*/
392  }/*for (ispin = 0; ispin <= Si2; ispin++)*/
393 }/*StdFace_GeneralJ*/
void StdFace_intr(struct StdIntList *StdI, std::complex< double > 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...

◆ StdFace_Hopping()

void StdFace_Hopping ( struct StdIntList *  StdI,
std::complex< double >  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 76 of file StdFace_ModelUtil.cpp.

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

83 {
84  int ispin;
90 #if defined(_HPhi)
91  int it, ii;
92  double Cphase;
93  std::complex<double> coef;
94 
95  if (strcmp(StdI->method, "timeevolution") == 0 && StdI->PumpBody == 1) {
96  for (it = 0; it < StdI->Lanczos_max; it++) {
97  Cphase = 0.0f;
98  for (ii = 0; ii < 3; ii++) Cphase += /*2.0*StdI->pi */ StdI->At[it][ii] * dR[ii];
99  coef = std::complex<double>(cos(Cphase), sin(-Cphase));
100  for (ispin = 0; ispin < 2; ispin++) {
101  StdI->pump[it][StdI->npump[it]] = coef * trans0;
102  StdI->pumpindx[it][StdI->npump[it]][0] = isite;
103  StdI->pumpindx[it][StdI->npump[it]][1] = ispin;
104  StdI->pumpindx[it][StdI->npump[it]][2] = jsite;
105  StdI->pumpindx[it][StdI->npump[it]][3] = ispin;
106  StdI->npump[it] = StdI->npump[it] + 1;
107 
108  StdI->pump[it][StdI->npump[it]] = conj(coef * trans0);
109  StdI->pumpindx[it][StdI->npump[it]][0] = jsite;
110  StdI->pumpindx[it][StdI->npump[it]][1] = ispin;
111  StdI->pumpindx[it][StdI->npump[it]][2] = isite;
112  StdI->pumpindx[it][StdI->npump[it]][3] = ispin;
113  StdI->npump[it] = StdI->npump[it] + 1;
114  }/*for (ispin = 0; ispin < 2; ispin++)*/
115  }/*for (it = 0; it < StdI->Lanczos_max; it++)*/
116  }/*if (strcmp(StdI->method, "timeevolution") == 0)*/
117  else
118 #endif
119  for (ispin = 0; ispin < 2; ispin++) {
120  StdFace_trans(StdI, trans0, jsite, ispin, isite, ispin);
121  StdFace_trans(StdI, conj(trans0), isite, ispin, jsite, ispin);
122  }/*for (ispin = 0; ispin < 2; ispin++)*/
123 }/*void StdFace_Hopping*/
void StdFace_trans(struct StdIntList *StdI, std::complex< double > trans0, int isite, int ispin, int jsite, int jspin)
Add transfer to the list set StdIntList::trans and StdIntList::transindx and increment StdIntList::nt...

◆ 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 129 of file StdFace_ModelUtil.cpp.

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

137 {
138  StdFace_trans(StdI, mu0 + 0.5 * h0, isite, 0, isite, 0);
139  StdFace_trans(StdI, mu0 - 0.5 * h0, isite, 1, isite, 1);
140  StdFace_trans(StdI, -0.5 * Gamma0, isite, 1, isite, 0);
141  StdFace_trans(StdI, -0.5 * Gamma0, isite, 0, isite, 1);
147  StdI->Cintra[StdI->NCintra] = U0;
148  StdI->CintraIndx[StdI->NCintra][0] = isite;
149  StdI->NCintra += 1;
150 }/*void StdFace_HubbardLocal*/
void StdFace_trans(struct StdIntList *StdI, std::complex< double > trans0, int isite, int ispin, int jsite, int jspin)
Add transfer to the list set StdIntList::trans and StdIntList::transindx and increment StdIntList::nt...

◆ 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 636 of file StdFace_ModelUtil.cpp.

References StdFace_exit(), StdFace_FoldSite(), and StdFace_PrintVal_i().

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

641 {
642  int bound[3][2], edge, ii, jj;
643  int ipos;
644  int nBox[3], iCellV_fold[3], iCellV[3];
645  double pos[4][2], xmin, xmax/*, offset[2], scale*/;
646 
647  fprintf(stdout, "\n @ Super-Lattice setting\n\n");
651  if (
652  (StdI->L != StdI->NaN_i || StdI->W != StdI->NaN_i || StdI->Height != StdI->NaN_i)
653  &&
654  (StdI->box[0][0] != StdI->NaN_i || StdI->box[0][1] != StdI->NaN_i || StdI->box[0][2] != StdI->NaN_i ||
655  StdI->box[1][0] != StdI->NaN_i || StdI->box[1][1] != StdI->NaN_i || StdI->box[1][2] != StdI->NaN_i ||
656  StdI->box[2][0] != StdI->NaN_i || StdI->box[2][1] != StdI->NaN_i || StdI->box[2][2] != StdI->NaN_i)
657  )
658  {
659  fprintf(stdout, "\nERROR ! (L, W, Height) and (a0W, ..., a2H) conflict !\n\n");
660  StdFace_exit(-1);
661  }
662  else if (StdI->L != StdI->NaN_i || StdI->W != StdI->NaN_i || StdI->Height != StdI->NaN_i)
663  {
664  StdFace_PrintVal_i("L", &StdI->L, 1);
665  StdFace_PrintVal_i("W", &StdI->W, 1);
666  StdFace_PrintVal_i("Height", &StdI->Height, 1);
667  for (ii = 0; ii < 3; ii++) for (jj = 0; jj < 3; jj++)
668  StdI->box[ii][jj] = 0;
669  StdI->box[0][0] = StdI->W;
670  StdI->box[1][1] = StdI->L;
671  StdI->box[2][2] = StdI->Height;
672  }
673  else
674  {
675  StdFace_PrintVal_i("a0W", &StdI->box[0][0], 1);
676  StdFace_PrintVal_i("a0L", &StdI->box[0][1], 0);
677  StdFace_PrintVal_i("a0H", &StdI->box[0][2], 0);
678  StdFace_PrintVal_i("a1W", &StdI->box[1][0], 0);
679  StdFace_PrintVal_i("a1L", &StdI->box[1][1], 1);
680  StdFace_PrintVal_i("a1H", &StdI->box[1][2], 0);
681  StdFace_PrintVal_i("a2W", &StdI->box[2][0], 0);
682  StdFace_PrintVal_i("a2L", &StdI->box[2][1], 0);
683  StdFace_PrintVal_i("a2H", &StdI->box[2][2], 1);
684  }
685  /*
686  Parameters for the 3D system will not used.
687  */
688  if (dim == 2) {
689  StdI->direct[0][2] = 0.0;
690  StdI->direct[1][2] = 0.0;
691  StdI->direct[2][0] = 0.0;
692  StdI->direct[2][1] = 0.0;
693  StdI->direct[2][2] = 1.0;
694  }
699  if (dim == 2) StdI->phase[2] = 0.0;
700  for (ii = 0; ii < 3; ii++) {
701  StdI->ExpPhase[ii] = std::complex<double>(cos(StdI->pi180 * StdI->phase[ii]), sin(StdI->pi180 * StdI->phase[ii]));
702  if (abs(StdI->ExpPhase[ii] + 1.0) < 0.000001) StdI->AntiPeriod[ii] = 1;
703  else StdI->AntiPeriod[ii] = 0;
704  }
708  StdI->tau = (double **)malloc(sizeof(double*) * StdI->NsiteUC);
709  for (ii = 0; ii < StdI->NsiteUC; ii++) {
710  StdI->tau[ii] = (double *)malloc(sizeof(double) * 3);
711  }
716  StdI->NCell = 0;
717  for (ii = 0; ii < 3; ii++) {
718  StdI->NCell += StdI->box[0][ii]
719  * StdI->box[1][(ii + 1) % 3]
720  * StdI->box[2][(ii + 2) % 3]
721  - StdI->box[0][ii]
722  * StdI->box[1][(ii + 2) % 3]
723  * StdI->box[2][(ii + 1) % 3];
724  }
725  printf(" Number of Cell = %d\n", abs(StdI->NCell));
726  if (StdI->NCell == 0) {
727  StdFace_exit(-1);
728  }
729 
730  for (ii = 0; ii < 3; ii++) {
731  for (jj = 0; jj < 3; jj++) {
732  StdI->rbox[ii][jj] = StdI->box[(ii + 1) % 3][(jj + 1) % 3] * StdI->box[(ii + 2) % 3][(jj + 2) % 3]
733  - StdI->box[(ii + 1) % 3][(jj + 2) % 3] * StdI->box[(ii + 2) % 3][(jj + 1) % 3];
734  }
735  }
736  if (StdI->NCell < 0) {
737  for (ii = 0; ii < 3; ii++)
738  for (jj = 0; jj < 3; jj++)
739  StdI->rbox[ii][jj] *= -1;
740  StdI->NCell *= -1;
741  }/*if (StdI->NCell < 0)*/
746  for (ii = 0; ii < 3; ii++) {
747  bound[ii][0] = 0;
748  bound[ii][1] = 0;
749  for (nBox[2] = 0; nBox[2] < 2; nBox[2]++) {
750  for (nBox[1] = 0; nBox[1] < 2; nBox[1]++) {
751  for (nBox[0] = 0; nBox[0] < 2; nBox[0]++) {
752  edge = 0;
753  for (jj = 0; jj < 3; jj++) edge += nBox[jj] * StdI->box[jj][ii];
754  if (edge < bound[ii][0]) bound[ii][0] = edge;
755  if (edge > bound[ii][1]) bound[ii][1] = edge;
756  }
757  }
758  }
759  }
763  StdI->Cell = (int **)malloc(sizeof(int*) * StdI->NCell);
764  for (ii = 0; ii < StdI->NCell; ii++) {
765  StdI->Cell[ii] = (int *)malloc(sizeof(int) * 3);
766  }/*for (ii = 0; ii < StdI->NCell; ii++)*/
767  jj = 0;
768  for (iCellV[2] = bound[2][0]; iCellV[2] <= bound[2][1]; iCellV[2]++) {
769  for (iCellV[1] = bound[1][0]; iCellV[1] <= bound[1][1]; iCellV[1]++) {
770  for (iCellV[0] = bound[0][0]; iCellV[0] <= bound[0][1]; iCellV[0]++) {
771  StdFace_FoldSite(StdI, iCellV, nBox, iCellV_fold);
772  if (nBox[0] == 0 && nBox[1] == 0 && nBox[2] == 0) {
773  for (ii = 0; ii < 3; ii++)
774  StdI->Cell[jj][ii] = iCellV[ii];
775  jj += 1;
776  }/*if (lUC == 1)*/
777  }/*for (iCellV[0] = bound[0][0]; iCellV[0] <= bound[0][1]; iCellV[0]++*/
778  }/*for (iCellV[1] = bound[1][0]; iCellV[1] <= bound[1][1]; iCellV[1]++)*/
779  }/*for (iCellV[2] = bound[2][0]; iCellV[2] <= bound[2][1]; iCellV[2]++)*/
783  if (dim == 2) {
784  pos[0][0] = 0.0;
785  pos[0][1] = 0.0;
786  pos[1][0] = StdI->direct[0][0] * (double)StdI->box[0][0] + StdI->direct[1][0] * (double)StdI->box[0][1];
787  pos[1][1] = StdI->direct[0][1] * (double)StdI->box[0][0] + StdI->direct[1][1] * (double)StdI->box[0][1];
788  pos[2][0] = StdI->direct[0][0] * (double)StdI->box[1][0] + StdI->direct[1][0] * (double)StdI->box[1][1];
789  pos[2][1] = StdI->direct[0][1] * (double)StdI->box[1][0] + StdI->direct[1][1] * (double)StdI->box[1][1];
790  pos[3][0] = pos[1][0] + pos[2][0];
791  pos[3][1] = pos[1][1] + pos[2][1];
792 
793  xmin = 0.0;
794  xmax = 0.0;
795  for (ipos = 0; ipos < 4; ipos++) {
796  if (pos[ipos][0] < xmin) xmin = pos[ipos][0];
797  if (pos[ipos][0] > xmax) xmax = pos[ipos][0];
798  if (pos[ipos][1] < xmin) xmin = pos[ipos][1];
799  if (pos[ipos][1] > xmax) xmax = pos[ipos][1];
800  }
801  xmin -= 2.0;
802  xmax += 2.0;
803 
804  fprintf(fp, "#set terminal pdf color enhanced \\\n");
805  fprintf(fp, "#dashed dl 1.0 size 20.0cm, 20.0cm \n");
806  fprintf(fp, "#set output \"lattice.pdf\"\n");
807  fprintf(fp, "set xrange [%f: %f]\n", xmin, xmax);
808  fprintf(fp, "set yrange [%f: %f]\n", xmin, xmax);
809  fprintf(fp, "set size square\n");
810  fprintf(fp, "unset key\n");
811  fprintf(fp, "unset tics\n");
812  fprintf(fp, "unset border\n");
813 
814  fprintf(fp, "set style line 1 lc 1 lt 1\n");
815  fprintf(fp, "set style line 2 lc 5 lt 1\n");
816  fprintf(fp, "set style line 3 lc 0 lt 1\n");
817 
818  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]);
819  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]);
820  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]);
821  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]);
822  }/*if (dim == 2)*/
823 }/*void StdFace_InitSite2D*/
void StdFace_exit(int errorcode)
MPI Abortation wrapper.
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.
void StdFace_PrintVal_i(const char *valname, int *val, int val0)
Print a valiable (integer) read from the input file if it is not specified in the input file (=214748...

◆ StdFace_InputCoulombV()

void StdFace_InputCoulombV ( double  V,
double *  V0,
const 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]V
[in]V0
[in]V0nameE.g. V1

Definition at line 1115 of file StdFace_ModelUtil.cpp.

References StdFace_exit().

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

1120 {
1121  if (std::isnan(V) == 0 && std::isnan(*V0) == 0) {
1122  fprintf(stdout, "\n ERROR! %s conflicts !\n\n", V0name);
1123  StdFace_exit(-1);
1124  }
1125  else if (std::isnan(*V0) == 0)
1126  fprintf(stdout, " %15s = %-10.5f\n", V0name, *V0);
1127  else if (std::isnan(V) == 0) {
1128  *V0 = V;
1129  fprintf(stdout, " %15s = %-10.5f\n", V0name, *V0);
1130  }
1131  else {
1132  *V0 = 0.0;
1133  }
1134 }/*void StdFace_InputCoulombV*/
void StdFace_exit(int errorcode)
MPI Abortation wrapper.

◆ StdFace_InputHopp()

void StdFace_InputHopp ( std::complex< double >  t,
std::complex< double > *  t0,
const 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]t
[in]t0
[in]t0nameE.g. t1

Definition at line 1140 of file StdFace_ModelUtil.cpp.

References StdFace_exit().

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

1145 {
1146  if (std::isnan(real(t)) == 0 && std::isnan(real(*t0)) == 0) {
1147  fprintf(stdout, "\n ERROR! %s conflicts !\n\n", t0name);
1148  StdFace_exit(-1);
1149  }
1150  else if (std::isnan(real(*t0)) == 0)
1151  fprintf(stdout, " %15s = %-10.5f %-10.5f\n", t0name, real(*t0), imag(*t0));
1152  else if (std::isnan(real(t)) == 0) {
1153  *t0 = t;
1154  fprintf(stdout, " %15s = %-10.5f %-10.5f\n", t0name, real(*t0), imag(*t0));
1155  }
1156  else {
1157  *t0 = 0.0;
1158  }
1159 }/*void StdFace_InputHopp*/
void StdFace_exit(int errorcode)
MPI Abortation wrapper.

◆ StdFace_InputSpin()

void StdFace_InputSpin ( double  Jp[3][3],
double  JpAll,
const char *  Jpname 
)

Input spin-spin interaction other than nearest-neighbor.

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

Definition at line 1067 of file StdFace_ModelUtil.cpp.

References StdFace_exit().

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

1072 {
1073  int i1, i2;
1074  char Jname[3][3][10];
1075 
1076  strcpy(Jname[0][0], "x\0");
1077  strcpy(Jname[0][1], "xy\0");
1078  strcpy(Jname[0][2], "xz\0");
1079  strcpy(Jname[1][0], "yx\0");
1080  strcpy(Jname[1][1], "y\0");
1081  strcpy(Jname[1][2], "yz\0");
1082  strcpy(Jname[2][0], "zx\0");
1083  strcpy(Jname[2][1], "zy\0");
1084  strcpy(Jname[2][2], "z\0");
1085 
1086  for (i1 = 0; i1 < 3; i1++) {
1087  for (i2 = 0; i2 < 3; i2++) {
1088  if (std::isnan(JpAll) == 0 && std::isnan(Jp[i1][i2]) == 0) {
1089  fprintf(stdout, "\n ERROR! %s and %s%s conflict !\n\n", Jpname,
1090  Jpname, Jname[i1][i2]);
1091  StdFace_exit(-1);
1092  }
1093  }/*for (j = 0; j < 3; j++)*/
1094  }/*for (i = 0; i < 3; i++)*/
1095 
1096  for (i1 = 0; i1 < 3; i1++) {
1097  for (i2 = 0; i2 < 3; i2++) {
1098  if (std::isnan(Jp[i1][i2]) == 0)
1099  fprintf(stdout, " %14s%s = %-10.5f\n", Jpname, Jname[i1][i2], Jp[i1][i2]);
1100  else if (i1 == i2 && std::isnan(JpAll) == 0) {
1101  Jp[i1][i2] = JpAll;
1102  fprintf(stdout, " %14s%s = %-10.5f\n", Jpname, Jname[i1][i2], Jp[i1][i2]);
1103  }
1104  else {
1105  Jp[i1][i2] = 0.0;
1106  }
1107  }/*for (i2 = 0; i2 < 3; i2++)*/
1108  }/*for (i = 0; i < 3; i++)*/
1109 }/*void StdFace_InputSpin*/
void StdFace_exit(int errorcode)
MPI Abortation wrapper.

◆ StdFace_InputSpinNN()

void StdFace_InputSpinNN ( double  J[3][3],
double  JAll,
double  J0[3][3],
double  J0All,
const char *  J0name 
)

Input nearest-neighbor spin-spin interaction.

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

Definition at line 979 of file StdFace_ModelUtil.cpp.

References StdFace_exit().

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

986 {
987  int i1, i2, i3, i4;
988  char Jname[3][3][10];
989 
990  strcpy(Jname[0][0], "x\0");
991  strcpy(Jname[0][1], "xy\0");
992  strcpy(Jname[0][2], "xz\0");
993  strcpy(Jname[1][0], "yx\0");
994  strcpy(Jname[1][1], "y\0");
995  strcpy(Jname[1][2], "yz\0");
996  strcpy(Jname[2][0], "zx\0");
997  strcpy(Jname[2][1], "zy\0");
998  strcpy(Jname[2][2], "z\0");
999 
1000  if (std::isnan(JAll) == 0 && std::isnan(J0All) == 0) {
1001  fprintf(stdout, "\n ERROR! %s conflict !\n\n", J0name);
1002  StdFace_exit(-1);
1003  }
1004  for (i1 = 0; i1 < 3; i1++) {
1005  for (i2 = 0; i2 < 3; i2++) {
1006  if (std::isnan(JAll) == 0 && std::isnan(J[i1][i2]) == 0) {
1007  fprintf(stdout, "\n ERROR! J%s conflict !\n\n", Jname[i1][i2]);
1008  StdFace_exit(-1);
1009  }
1010  else if (std::isnan(J0All) == 0 && std::isnan(J[i1][i2]) == 0) {
1011  fprintf(stdout, "\n ERROR! %s and J%s conflict !\n\n",
1012  J0name, Jname[i1][i2]);
1013  StdFace_exit(-1);
1014  }
1015  else if (std::isnan(J0All) == 0 && std::isnan(J0[i1][i2]) == 0) {
1016  fprintf(stdout, "\n ERROR! %s and %s%s conflict !\n\n", J0name,
1017  J0name, Jname[i1][i2]);
1018  StdFace_exit(-1);
1019  }
1020  else if (std::isnan(J0[i1][i2]) == 0 && std::isnan(JAll) == 0) {
1021  fprintf(stdout, "\n ERROR! %s%s conflict !\n\n",
1022  J0name, Jname[i1][i2]);
1023  StdFace_exit(-1);
1024  }
1025  }/*for (j = 0; j < 3; j++)*/
1026  }/*for (i = 0; i < 3; i++)*/
1027 
1028  for (i1 = 0; i1 < 3; i1++) {
1029  for (i2 = 0; i2 < 3; i2++) {
1030  for (i3 = 0; i3 < 3; i3++) {
1031  for (i4 = 0; i4 < 3; i4++) {
1032  if (std::isnan(J0[i1][i2]) == 0 && std::isnan(J[i3][i4]) == 0) {
1033  fprintf(stdout, "\n ERROR! %s%s and J%s conflict !\n\n",
1034  J0name, Jname[i1][i2], Jname[i3][i4]);
1035  StdFace_exit(-1);
1036  }
1037  }/*for (i4 = 0; i4 < 3; i4++)*/
1038  }/*for (i3 = 0; i3 < 3; i3++)*/
1039  }/*for (j = 0; j < 3; j++)*/
1040  }/*for (i = 0; i < 3; i++)*/
1041 
1042  for (i1 = 0; i1 < 3; i1++) {
1043  for (i2 = 0; i2 < 3; i2++) {
1044  if (std::isnan(J0[i1][i2]) == 0)
1045  fprintf(stdout, " %14s%s = %-10.5f\n", J0name, Jname[i1][i2], J0[i1][i2]);
1046  else if (std::isnan(J[i1][i2]) == 0) {
1047  J0[i1][i2] = J[i1][i2];
1048  fprintf(stdout, " %14s%s = %-10.5f\n", J0name, Jname[i1][i2], J0[i1][i2]);
1049  }
1050  else if (i1 == i2 && std::isnan(J0All) == 0) {
1051  J0[i1][i2] = J0All;
1052  fprintf(stdout, " %14s%s = %-10.5f\n", J0name, Jname[i1][i2], J0[i1][i2]);
1053  }
1054  else if (i1 == i2 && std::isnan(JAll) == 0) {
1055  J0[i1][i2] = JAll;
1056  fprintf(stdout, " %14s%s = %-10.5f\n", J0name, Jname[i1][i2], J0[i1][i2]);
1057  }
1058  else {
1059  J0[i1][i2] = 0.0;
1060  }
1061  }/*for (i2 = 0; i2 < 3; i2++)*/
1062  }/*for (i = 0; i < 3; i++)*/
1063 }/*void StdFace_InputSpinNN*/
void StdFace_exit(int errorcode)
MPI Abortation wrapper.

◆ StdFace_intr()

void StdFace_intr ( struct StdIntList *  StdI,
std::complex< double >  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 203 of file StdFace_ModelUtil.cpp.

Referenced by StdFace_GeneralJ().

215 {
216  StdI->intr[StdI->nintr] = intr0;
217  StdI->intrindx[StdI->nintr][0] = site1; StdI->intrindx[StdI->nintr][1] = spin1;
218  StdI->intrindx[StdI->nintr][2] = site2; StdI->intrindx[StdI->nintr][3] = spin2;
219  StdI->intrindx[StdI->nintr][4] = site3; StdI->intrindx[StdI->nintr][5] = spin3;
220  StdI->intrindx[StdI->nintr][6] = site4; StdI->intrindx[StdI->nintr][7] = spin4;
221  StdI->nintr = StdI->nintr + 1;
222 }/*void StdFace_intr*/

◆ 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 155 of file StdFace_ModelUtil.cpp.

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

162 {
163  int ispin;
164  double S, Sz;
165 
166  S = (double)S2 * 0.5;
170  for (ispin = 0; ispin <= S2; ispin++){
177  Sz = (double)ispin - S;
178  StdFace_trans(StdI, -h * Sz, isite, ispin, isite, ispin);
189  if (ispin < S2) {
190  StdFace_trans(StdI, -0.5 * Gamma * sqrt(S*(S + 1.0) - Sz*(Sz + 1.0)),
191  isite, ispin + 1, isite, ispin);
192  StdFace_trans(StdI, -0.5 * Gamma * sqrt(S*(S + 1.0) - Sz*(Sz + 1.0)),
193  isite, ispin, isite, ispin + 1);
194  }/*if (ispin < S2)*/
195  }/*for (ispin = 0; ispin <= S2; ispin++)*/
196 }/*void StdFace_MagField*/
void StdFace_trans(struct StdIntList *StdI, std::complex< double > trans0, int isite, int ispin, int jsite, int jspin)
Add transfer to the list set StdIntList::trans and StdIntList::transindx and increment StdIntList::nt...

◆ 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 1204 of file StdFace_ModelUtil.cpp.

References StdFace_exit(), StdFace_FindSite(), StdFace_FoldSite(), and StdFace_PrintVal_i().

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

1208  {
1209  int ii;
1210 #if defined(_HPhi)
1211  int it;
1212 #endif
1213 
1216  StdI->transindx = (int **)malloc(sizeof(int*) * ntransMax);
1217  StdI->trans = (std::complex<double> *)malloc(sizeof(std::complex<double>) * ntransMax);
1218  for (ii = 0; ii < ntransMax; ii++) {
1219  StdI->transindx[ii] = (int *)malloc(sizeof(int) * 4);
1220  }
1221  StdI->ntrans = 0;
1222 #if defined(_HPhi)
1223  if (strcmp(StdI->method, "timeevolution") == 0 && StdI->PumpBody == 1) {
1224  StdI->npump = (int *)malloc(sizeof(int) * StdI->Lanczos_max);
1225  StdI->pumpindx = (int ***)malloc(sizeof(int**) * StdI->Lanczos_max);
1226  StdI->pump = (std::complex<double> **)malloc(sizeof(std::complex<double>*) * StdI->Lanczos_max);
1227  for (it = 0; it < StdI->Lanczos_max; it++) {
1228  StdI->npump[it] = 0;
1229  StdI->pumpindx[it] = (int **)malloc(sizeof(int*) * ntransMax);
1230  StdI->pump[it] = (std::complex<double> *)malloc(sizeof(std::complex<double>) * ntransMax);
1231  for (ii = 0; ii < ntransMax; ii++) {
1232  StdI->pumpindx[it][ii] = (int *)malloc(sizeof(int) * 4);
1233  }
1234  }/*for (it = 0; it < StdI->Lanczos_max;)*/
1235  }/*if (strcmp(StdI->method, "timeevolution") == 0*/
1236 #endif
1237 
1240  StdI->intrindx = (int **)malloc(sizeof(int*) * nintrMax);
1241  StdI->intr = (std::complex<double> *)malloc(sizeof(std::complex<double>) * nintrMax);
1242  for (ii = 0; ii < nintrMax; ii++) {
1243  StdI->intrindx[ii] = (int *)malloc(sizeof(int) * 8);
1244  }
1245  StdI->nintr = 0;
1249  StdI->CintraIndx = (int **)malloc(sizeof(int*) * nintrMax);
1250  StdI->Cintra = (double *)malloc(sizeof(double) * nintrMax);
1251  for (ii = 0; ii < nintrMax; ii++) {
1252  StdI->CintraIndx[ii] = (int *)malloc(sizeof(int) * 1);
1253  }
1254  StdI->NCintra = 0;
1258  StdI->CinterIndx = (int **)malloc(sizeof(int*) * nintrMax);
1259  StdI->Cinter = (double *)malloc(sizeof(double) * nintrMax);
1260  for (ii = 0; ii < nintrMax; ii++) {
1261  StdI->CinterIndx[ii] = (int *)malloc(sizeof(int) * 2);
1262  }
1263  StdI->NCinter = 0;
1267  StdI->HundIndx = (int **)malloc(sizeof(int*) * nintrMax);
1268  StdI->Hund = (double *)malloc(sizeof(double) * nintrMax);
1269  for (ii = 0; ii < nintrMax; ii++) {
1270  StdI->HundIndx[ii] = (int *)malloc(sizeof(int) * 2);
1271  }
1272  StdI->NHund = 0;
1276  StdI->ExIndx = (int **)malloc(sizeof(int*) * nintrMax);
1277  StdI->Ex = (double *)malloc(sizeof(double) * nintrMax);
1278  for (ii = 0; ii < nintrMax; ii++) {
1279  StdI->ExIndx[ii] = (int *)malloc(sizeof(int) * 2);
1280  }
1281  StdI->NEx = 0;
1285  StdI->PLIndx = (int **)malloc(sizeof(int*) * nintrMax);
1286  StdI->PairLift = (double *)malloc(sizeof(double) * nintrMax);
1287  for (ii = 0; ii < nintrMax; ii++) {
1288  StdI->PLIndx[ii] = (int *)malloc(sizeof(int) * 2);
1289  }
1290  StdI->NPairLift = 0;
1294  StdI->PHIndx = (int **)malloc(sizeof(int*) * nintrMax);
1295  StdI->PairHopp = (double *)malloc(sizeof(double) * nintrMax);
1296  for (ii = 0; ii < nintrMax; ii++) {
1297  StdI->PHIndx[ii] = (int *)malloc(sizeof(int) * 2);
1298  }
1299  StdI->NPairHopp = 0;
1300 }/*void StdFace_MallocInteractions*/

◆ StdFace_NotUsed_c()

void StdFace_NotUsed_c ( const char *  valname,
std::complex< double >  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 513 of file StdFace_ModelUtil.cpp.

References StdFace_exit().

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

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

◆ StdFace_NotUsed_d()

void StdFace_NotUsed_d ( const 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 496 of file StdFace_ModelUtil.cpp.

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

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

◆ StdFace_NotUsed_i()

void StdFace_NotUsed_i ( const char *  valname,
int  val 
)

Stop HPhi if a variable (integer) not used is specified in the input file (!=2147483647, the upper limt of Int).

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

Definition at line 562 of file StdFace_ModelUtil.cpp.

References StdFace_exit().

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

566 {
567  int NaN_i = 2147483647;
568 
569  if (val != NaN_i) {
570  fprintf(stdout, "\n Check ! %s is SPECIFIED but will NOT be USED. \n", valname);
571  fprintf(stdout, " Please COMMENT-OUT this line \n");
572  fprintf(stdout, " or check this input is REALLY APPROPRIATE for your purpose ! \n\n");
573  StdFace_exit(-1);
574  }
575 }/*void StdFace_NotUsed_i*/
void StdFace_exit(int errorcode)
MPI Abortation wrapper.

◆ StdFace_NotUsed_J()

void StdFace_NotUsed_J ( const 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 530 of file StdFace_ModelUtil.cpp.

References StdFace_NotUsed_d().

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

535 {
536  int i1, i2;
537  char Jname[3][3][10];
538 
539  sprintf(Jname[0][0], "%sx", valname);
540  sprintf(Jname[0][1], "%sxy", valname);
541  sprintf(Jname[0][2], "%sxz", valname);
542  sprintf(Jname[1][0], "%syx", valname);
543  sprintf(Jname[1][1], "%sy", valname);
544  sprintf(Jname[1][2], "%syz", valname);
545  sprintf(Jname[2][0], "%szx", valname);
546  sprintf(Jname[2][1], "%szy", valname);
547  sprintf(Jname[2][2], "%sz", valname);
548 
549  StdFace_NotUsed_d(valname, JAll);
550 
551  for (i1 = 0; i1 < 3; i1++) {
552  for (i2 = 0; i2 < 3; i2++) {
553  StdFace_NotUsed_d(Jname[i1][i2], J[i1][i2]);
554  }/*for (j = 0; j < 3; j++)*/
555  }/*for (i = 0; i < 3; i++)*/
556 }/*void StdFace_NotUsed_J*/
void StdFace_NotUsed_d(const 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 1163 of file StdFace_ModelUtil.cpp.

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

1163  {
1164  FILE *fp;
1165  int isite, iCell, ii;
1166 
1167  fp = fopen("geometry.dat", "w");
1168 
1169  for (ii = 0; ii < 3; ii++)
1170  fprintf(fp, "%25.15e %25.15e %25.15e\n",
1171  StdI->direct[ii][0], StdI->direct[ii][1], StdI->direct[ii][2]);
1172  fprintf(fp, "%25.15e %25.15e %25.15e\n",
1173  StdI->phase[0], StdI->phase[1], StdI->phase[2]);
1174  for (ii = 0; ii < 3; ii++)
1175  fprintf(fp, "%d %d %d\n",
1176  StdI->box[ii][0], StdI->box[ii][1], StdI->box[ii][2]);
1177 
1178  for (iCell = 0; iCell < StdI->NCell; iCell++) {
1179  for (isite = 0; isite < StdI->NsiteUC; isite++) {
1180  fprintf(fp, "%d %d %d %d\n",
1181  StdI->Cell[iCell][0] - StdI->Cell[0][0],
1182  StdI->Cell[iCell][1] - StdI->Cell[0][1],
1183  StdI->Cell[iCell][2] - StdI->Cell[0][2],
1184  isite);
1185  }/*for (isite = 0; isite < StdI->NsiteUC; isite++)*/
1186  }/* for (iCell = 0; iCell < StdI->NCell; iCell++)*/
1187  if (strcmp(StdI->model, "kondo") == 0) {
1188  for (iCell = 0; iCell < StdI->NCell; iCell++) {
1189  for (isite = 0; isite < StdI->NsiteUC; isite++) {
1190  fprintf(fp, "%d %d %d %d\n",
1191  StdI->Cell[iCell][0] - StdI->Cell[0][0],
1192  StdI->Cell[iCell][1] - StdI->Cell[0][1],
1193  StdI->Cell[iCell][2] - StdI->Cell[0][2],
1194  isite + StdI->NsiteUC);
1195  }/*for (isite = 0; isite < StdI->NsiteUC; isite++)*/
1196  }/* for (iCell = 0; iCell < StdI->NCell; iCell++)*/
1197  }
1198  fflush(fp);
1199  fclose(fp);
1200 }/*void StdFace_PrintGeometry()*/

◆ StdFace_PrintVal_c()

void StdFace_PrintVal_c ( const char *  valname,
std::complex< double > *  val,
std::complex< double >  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 459 of file StdFace_ModelUtil.cpp.

Referenced by StdFace_Orthorhombic().

464 {
465  if (std::isnan(real(*val)) == 1) {
466  *val = val0;
467  fprintf(stdout, " %15s = %-10.5f %-10.5f ###### DEFAULT VALUE IS USED ######\n", valname, real(*val), imag(*val));
468  }
469  else fprintf(stdout, " %15s = %-10.5f %-10.5f\n", valname, real(*val), imag(*val));
470 }/*void StdFace_PrintVal_c*/

◆ StdFace_PrintVal_d()

void StdFace_PrintVal_d ( const char *  valname,
double *  val,
double  val0 
)

Print a valiable (real) read from the input file if it is not specified in the input file (=NaN), 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 418 of file StdFace_ModelUtil.cpp.

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

423 {
424  if (std::isnan(*val) == 1) {
425  *val = val0;
426  fprintf(stdout, " %15s = %-10.5f ###### DEFAULT VALUE IS USED ######\n", valname, *val);
427  }
428  else fprintf(stdout, " %15s = %-10.5f\n", valname, *val);
429 }/*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 436 of file StdFace_ModelUtil.cpp.

442 {
443  if (std::isnan(*val) == 1) {
447  if (std::isnan(val0) == 1) *val = val1;
448  else *val = val0;
449  fprintf(stdout, " %15s = %-10.5f ###### DEFAULT VALUE IS USED ######\n", valname, *val);
450  }
451  else fprintf(stdout, " %15s = %-10.5f\n", valname, *val);
452 }/*void StdFace_PrintVal_dd*/

◆ StdFace_PrintVal_i()

void StdFace_PrintVal_i ( const char *  valname,
int *  val,
int  val0 
)

Print a valiable (integer) read from the input file if it is not specified in the input file (=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 477 of file StdFace_ModelUtil.cpp.

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

482 {
483  int NaN_i = 2147483647;/*The upper limt of Int*/
484 
485  if (*val == NaN_i) {
486  *val = val0;
487  fprintf(stdout, " %15s = %-10d ###### DEFAULT VALUE IS USED ######\n", valname, *val);
488  }
489  else fprintf(stdout, " %15s = %-10d\n", valname, *val);
490 }/*void StdFace_PrintVal_i*/

◆ StdFace_PrintXSF()

void StdFace_PrintXSF ( struct StdIntList *  StdI)

Print lattice.xsf (XCrysDen format)

Definition at line 944 of file StdFace_ModelUtil.cpp.

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

944  {
945  FILE *fp;
946  int ii, jj, kk, isite, iCell;
947  double vec[3];
948 
949  fp = fopen("lattice.xsf", "w");
950  fprintf(fp, "CRYSTAL\n");
951  fprintf(fp, "PRIMVEC\n");
952  for (ii = 0; ii < 3; ii++) {
953  for (jj = 0; jj < 3; jj++) {
954  vec[jj] = 0.0;
955  for (kk = 0; kk < 3; kk++)
956  vec[jj] += (double)StdI->box[ii][kk] * StdI->direct[kk][jj];
957  }
958  fprintf(fp, "%15.5f %15.5f %15.5f\n", vec[0], vec[1], vec[2]);
959  }
960  fprintf(fp, "PRIMCOORD\n");
961  fprintf(fp, "%d 1\n", StdI->NCell * StdI->NsiteUC);
962  for (iCell = 0; iCell < StdI->NCell; iCell++) {
963  for (isite = 0; isite < StdI->NsiteUC; isite++) {
964  for (jj = 0; jj < 3; jj++) {
965  vec[jj] = 0.0;
966  for (kk = 0; kk < 3; kk++)
967  vec[jj] += ((double)StdI->Cell[iCell][kk] + StdI->tau[isite][kk])
968  * StdI->direct[kk][jj];
969  }
970  fprintf(fp, "H %15.5f %15.5f %15.5f\n", vec[0], vec[1], vec[2]);
971  }
972  }
973  fflush(fp);
974  fclose(fp);
975 }/*void StdFace_PrintXSF*/

◆ StdFace_RequiredVal_i()

void StdFace_RequiredVal_i ( const char *  valname,
int  val 
)

Stop HPhi if a variable (integer) which must be specified is absent in the input file (=2147483647, the upper limt of Int).

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

Definition at line 581 of file StdFace_ModelUtil.cpp.

References StdFace_exit().

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

585 {
586  int NaN_i = 2147483647;
587 
588  if (val == NaN_i){
589  fprintf(stdout, "ERROR ! %s is NOT specified !\n", valname);
590  StdFace_exit(-1);
591  }
592  else fprintf(stdout, " %15s = %-3d\n", valname, val);
593 }/*void StdFace_RequiredVal_i*/
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,
std::complex< double > *  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 881 of file StdFace_ModelUtil.cpp.

References StdFace_FindSite().

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

896 {
897  double xi, yi, xj, yj;
901  StdFace_FindSite(StdI, iW, iL, 0, -diW, -diL, 0, jsiteUC, isiteUC, isite, jsite, Cphase, dR);
902 
903  xi = StdI->direct[0][0] * ((double)iW + StdI->tau[jsiteUC][0])
904  + StdI->direct[1][0] * ((double)iL + StdI->tau[jsiteUC][1]);
905  yi = StdI->direct[0][1] * ((double)iW + StdI->tau[jsiteUC][0])
906  + StdI->direct[1][1] * ((double)iL + StdI->tau[jsiteUC][1]);
907 
908  xj = StdI->direct[0][0] * ((double)(iW - diW) + StdI->tau[isiteUC][0])
909  + StdI->direct[1][0] * ((double)(iL - diL) + StdI->tau[isiteUC][1]);
910  yj = StdI->direct[0][1] * ((double)(iW - diW) + StdI->tau[isiteUC][0])
911  + StdI->direct[1][1] * ((double)(iL - diL) + StdI->tau[isiteUC][1]);
912 
913  if (*isite < 10)fprintf(fp, "set label \"%1d\" at %f, %f center front\n", *isite, xi, yi);
914  else fprintf(fp, "set label \"%2d\" at %f, %f center front\n", *isite, xi, yi);
915  if (*jsite < 10)fprintf(fp, "set label \"%1d\" at %f, %f center front\n", *jsite, xj, yj);
916  else fprintf(fp, "set label \"%2d\" at %f, %f center front\n", *jsite, xj, yj);
917  if (connect < 3)
918  fprintf(fp, "set arrow from %f, %f to %f, %f nohead ls %d\n", xi, yi, xj, yj, connect);
922  StdFace_FindSite(StdI, iW, iL, 0, diW, diL, 0, isiteUC, jsiteUC, isite, jsite, Cphase, dR);
923 
924  xi = StdI->direct[1][0] * ((double)iL + StdI->tau[isiteUC][1])
925  + StdI->direct[0][0] * ((double)iW + StdI->tau[isiteUC][0]);
926  yi = StdI->direct[1][1] * ((double)iL + StdI->tau[isiteUC][1])
927  + StdI->direct[0][1] * ((double)iW + StdI->tau[isiteUC][0]);
928 
929  xj = StdI->direct[0][0] * ((double)(iW + diW) + StdI->tau[jsiteUC][0])
930  + StdI->direct[1][0] * ((double)(iL + diL) + StdI->tau[jsiteUC][1]);
931  yj = StdI->direct[0][1] * ((double)(iW + diW) + StdI->tau[jsiteUC][0])
932  + StdI->direct[1][1] * ((double)(iL + diL) + StdI->tau[jsiteUC][1]);
933 
934  if (*isite < 10)fprintf(fp, "set label \"%1d\" at %f, %f center front\n", *isite, xi, yi);
935  else fprintf(fp, "set label \"%2d\" at %f, %f center front\n", *isite, xi, yi);
936  if (*jsite < 10)fprintf(fp, "set label \"%1d\" at %f, %f center front\n", *jsite, xj, yj);
937  else fprintf(fp, "set label \"%2d\" at %f, %f center front\n", *jsite, xj, yj);
938  if (connect < 3)
939  fprintf(fp, "set arrow from %f, %f to %f, %f nohead ls %d\n", xi, yi, xj, yj, connect);
940 }/*void StdFace_SetLabel*/
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, std::complex< double > *Cphase, double *dR)
Find the index of transfer and interaction.

◆ StdFace_trans()

void StdFace_trans ( struct StdIntList *  StdI,
std::complex< double >  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 56 of file StdFace_ModelUtil.cpp.

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

64 {
65  StdI->trans[StdI->ntrans] = trans0;
66  StdI->transindx[StdI->ntrans][0] = isite;
67  StdI->transindx[StdI->ntrans][1] = ispin;
68  StdI->transindx[StdI->ntrans][2] = jsite;
69  StdI->transindx[StdI->ntrans][3] = jspin;
70  StdI->ntrans = StdI->ntrans + 1;
71 }/*void StdFace_trans*/