HΦ  3.2.0
mltplySpinCore.c File Reference

Functions for spin Hamiltonian (Core) More...

#include <bitcalc.h>
#include "xsetmem.h"
#include "wrapperMPI.h"
#include "mltplyCommon.h"
#include "mltplySpinCore.h"
+ Include dependency graph for mltplySpinCore.c:

Go to the source code of this file.

Functions

int child_exchange_spin_GetInfo (int iExchange, struct BindStruct *X)
 Set parameters for the bit operation of spin-exchange term. More...
 
int child_pairlift_spin_GetInfo (int iPairLift, struct BindStruct *X)
 Set parameters for the bit operation of spin-pairlift term. More...
 
int child_general_int_spin_GetInfo (struct BindStruct *X, long unsigned int isite1, long unsigned int isite2, long unsigned int sigma1, long unsigned int sigma2, long unsigned int sigma3, long unsigned int sigma4, double complex tmp_V)
 Set parameters for the bit operation of spin-general interaction term. More...
 
int X_Spin_CisAit (long unsigned int j, struct BindStruct *X, long unsigned int is1_spin, long unsigned int sigma2, long unsigned int *list_1_Org_, long unsigned int *list_2_1_, long unsigned int *list_2_2_, long unsigned int *tmp_off)
 Compute index of final wavefunction by \(c_{is}^\dagger c_{it}\) term. More...
 
int X_Spin_CisAis (long unsigned int j, struct BindStruct *X, long unsigned int is1_spin, long unsigned int sigma1)
 Compute the spin state with bit mask is1_spin. More...
 
int X_SpinGC_CisAis (long unsigned int j, struct BindStruct *X, long unsigned int is1_spin, long unsigned int sigma1)
 Compute the grandcanonical spin state with bit mask is1_spin. More...
 
int X_SpinGC_CisAit (long unsigned int j, struct BindStruct *X, long unsigned int is1_spin, long unsigned int sigma2, long unsigned int *tmp_off)
 Compute index of final wavefunction by \(c_{is}^\dagger c_{it}\) term (grandcanonical). More...
 
int X_child_exchange_spin_element (long unsigned int j, struct BindStruct *X, long unsigned int isA_up, long unsigned int isB_up, long unsigned int sigmaA, long unsigned int sigmaB, long unsigned int *tmp_off)
 Compute index of final wavefunction associated to spin-exchange term. More...
 
double complex child_exchange_spin_element (long unsigned int j, double complex *tmp_v0, double complex *tmp_v1, struct BindStruct *X, long unsigned int *tmp_off)
 Multiply Hamiltonian of exchange term of canonical spin system. More...
 
double complex GC_child_exchange_spin_element (long unsigned int j, double complex *tmp_v0, double complex *tmp_v1, struct BindStruct *X, long unsigned int *tmp_off)
 Multiply Hamiltonian of exchange term of grandcanonical spin system. More...
 
double complex GC_child_pairlift_spin_element (long unsigned int j, double complex *tmp_v0, double complex *tmp_v1, struct BindStruct *X, long unsigned int *tmp_off)
 Multiply Hamiltonian of pairlift term of grandcanonical spin system. More...
 
double complex child_CisAisCisAis_spin_element (long unsigned int j, long unsigned int isA_up, long unsigned int isB_up, long unsigned int org_sigma2, long unsigned int org_sigma4, double complex tmp_V, double complex *tmp_v0, double complex *tmp_v1, struct BindStruct *X)
 Compute \(c_{is}^\dagger c_{is} c_{is}^\dagger c_{is}\) term of canonical spsin system. More...
 
double complex GC_child_CisAisCisAis_spin_element (long unsigned int j, long unsigned int isA_up, long unsigned int isB_up, long unsigned int org_sigma2, long unsigned int org_sigma4, double complex tmp_V, double complex *tmp_v0, double complex *tmp_v1, struct BindStruct *X)
 Compute \(c_{is}^\dagger c_{is} c_{is}^\dagger c_{is}\) term of grandcanonical spsin system. More...
 
double complex GC_child_CisAisCitAiu_spin_element (long unsigned int j, long unsigned int org_sigma2, long unsigned int org_sigma4, long unsigned int isA_up, long unsigned int isB_up, double complex tmp_V, double complex *tmp_v0, double complex *tmp_v1, struct BindStruct *X, long unsigned int *tmp_off)
 Compute \(c_{is}^\dagger c_{is} c_{it}^\dagger c_{iu}\) term of grandcanonical spsin system. More...
 
double complex GC_child_CisAitCiuAiu_spin_element (long unsigned int j, long unsigned int org_sigma2, long unsigned int org_sigma4, long unsigned int isA_up, long unsigned int isB_up, double complex tmp_V, double complex *tmp_v0, double complex *tmp_v1, struct BindStruct *X, long unsigned int *tmp_off)
 Compute \(c_{is}^\dagger c_{it} c_{iu}^\dagger c_{iu}\) term of grandcanonical spsin system. More...
 
double complex GC_child_CisAitCiuAiv_spin_element (long unsigned int j, long unsigned int org_sigma2, long unsigned int org_sigma4, long unsigned int isA_up, long unsigned int isB_up, double complex tmp_V, double complex *tmp_v0, double complex *tmp_v1, struct BindStruct *X, long unsigned int *tmp_off_2)
 Compute \(c_{is}^\dagger c_{it} c_{iu}^\dagger c_{iv}\) term of grandcanonical spsin system. More...
 

Detailed Description

Functions for spin Hamiltonian (Core)

Definition in file mltplySpinCore.c.

Function Documentation

◆ child_CisAisCisAis_spin_element()

double complex child_CisAisCisAis_spin_element ( long unsigned int  j,
long unsigned int  isA_up,
long unsigned int  isB_up,
long unsigned int  org_sigma2,
long unsigned int  org_sigma4,
double complex  tmp_V,
double complex *  tmp_v0,
double complex *  tmp_v1,
struct BindStruct X 
)

Compute \(c_{is}^\dagger c_{is} c_{is}^\dagger c_{is}\) term of canonical spsin system.

Returns
Fragment of \(\langle v_1 | H_{\rm this}|v_1 \rangle\)
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]jIndex of initial wavefunction
[in]isA_upBit mask for spin 1
[in]isB_upBit mask for spin 2
[in]org_sigma2Target for spin 1
[in]org_sigma4Target for spin 2
[in]tmp_VCoupling constatnt
[in]tmp_v0Resulting wavefunction
[in]tmp_v1Wavefunction to be multiplied
[in,out]X

Definition at line 391 of file mltplySpinCore.c.

References BindStruct::Large, LargeList::mode, and X_Spin_CisAis().

Referenced by expec_cisajscktalt_SpinHalf().

401  {
402  int tmp_sgn;
403  double complex dmv;
404  double complex dam_pr = 0;
405 
406  tmp_sgn = X_Spin_CisAis(j, X, isB_up, org_sigma4);
407  tmp_sgn *= X_Spin_CisAis(j, X, isA_up, org_sigma2);
408  dmv = tmp_v1[j] * tmp_sgn * tmp_V;
409  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) { // for multply
410  tmp_v0[j] += dmv;
411  }
412  dam_pr = conj(tmp_v1[j]) * dmv;
413  return dam_pr;
414 }/*double complex child_CisAisCisAis_spin_element*/
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.h:412
int X_Spin_CisAis(long unsigned int j, struct BindStruct *X, long unsigned int is1_spin, long unsigned int sigma1)
Compute the spin state with bit mask is1_spin.
int mode
multiply or expectation value.
Definition: struct.h:330
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ child_exchange_spin_element()

double complex child_exchange_spin_element ( long unsigned int  j,
double complex *  tmp_v0,
double complex *  tmp_v1,
struct BindStruct X,
long unsigned int *  tmp_off 
)

Multiply Hamiltonian of exchange term of canonical spin system.

Returns
\(\langle v_1 | H_{\rm this}| v_1 \rangle\)
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]jIndex of initial wavefunction
[out]tmp_v0Resulting wavefunction
[in]tmp_v1Wavefunction to be multiplied
[in,out]X
[out]tmp_offIndex of final wavefunction

Definition at line 273 of file mltplySpinCore.c.

References GetOffComp(), LargeList::ihfbit, LargeList::ilft, LargeList::irght, LargeList::isA_spin, BindStruct::Large, list_1, list_2_1, list_2_2, LargeList::mode, and LargeList::tmp_J.

Referenced by child_exchange_spin(), and makeHam().

279  {
280  long unsigned int off;
281  double complex dmv;
282  long unsigned int iexchg;
283  long unsigned int is_up = X->Large.isA_spin;
284  long unsigned int irght = X->Large.irght;
285  long unsigned int ilft = X->Large.ilft;
286  long unsigned int ihfbit = X->Large.ihfbit;
287  double complex tmp_J = X->Large.tmp_J;
288  int mode = X->Large.mode;
289  double complex dam_pr = 0;
290  long unsigned int ibit_tmp;
291 
292  ibit_tmp = (list_1[j] & is_up);
293  if (ibit_tmp == 0 || ibit_tmp == is_up) {
294  return dam_pr;
295  }
296  else {
297  iexchg = list_1[j] ^ is_up;
298  GetOffComp(list_2_1, list_2_2, iexchg, irght, ilft, ihfbit, &off);
299  *tmp_off = off;
300  dmv = tmp_J * tmp_v1[j];
301  if (mode == M_MLTPLY) {
302  tmp_v0[off] += dmv;
303  }
304  dam_pr += dmv * conj(tmp_v1[off]);
305  return dam_pr;
306  }
307 }/*double complex child_exchange_spin_element*/
long unsigned int ihfbit
Used for Ogata-Lin ???
Definition: struct.h:345
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.h:412
int mode
multiply or expectation value.
Definition: struct.h:330
double complex tmp_J
Coupling constant.
Definition: struct.h:323
long unsigned int ilft
Used for Ogata-Lin ???
Definition: struct.h:344
long unsigned int * list_2_1
Definition: global.h:49
int GetOffComp(long unsigned int *_list_2_1, long unsigned int *_list_2_2, long unsigned int _ibit, const long unsigned int _irght, const long unsigned int _ilft, const long unsigned int _ihfbit, long unsigned int *_ioffComp)
function of getting off-diagonal component
Definition: bitcalc.c:195
long unsigned int * list_1
Definition: global.h:47
long unsigned int * list_2_2
Definition: global.h:50
long unsigned int isA_spin
Mask used in the bit oeration.
Definition: struct.h:346
long unsigned int irght
Used for Ogata-Lin ???
Definition: struct.h:343
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ child_exchange_spin_GetInfo()

int child_exchange_spin_GetInfo ( int  iExchange,
struct BindStruct X 
)

Set parameters for the bit operation of spin-exchange term.

Returns
Always return 0
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)

Set the exchange coupling constant (LargeList::tmp_J)

Set the bit mask for computing spin state of both site (LargeList::is1_up, LargeList::is2_up)

Set the bit mask for exchange 2 spins (LargeList::isA_spin)

Parameters
[in]iExchangeCounter of exchange interaction
[in,out]X

Definition at line 35 of file mltplySpinCore.c.

References BindStruct::Def, DefineList::ExchangeCoupling, LargeList::is1_up, LargeList::is2_up, LargeList::isA_spin, BindStruct::Large, DefineList::ParaExchangeCoupling, LargeList::tmp_J, and DefineList::Tpow.

Referenced by makeHam(), mltplyHalfSpin(), and mltplyHalfSpinGC().

38  {
39  int isite1 = X->Def.ExchangeCoupling[iExchange][0] + 1;
40  int isite2 = X->Def.ExchangeCoupling[iExchange][1] + 1;
44  X->Large.tmp_J = X->Def.ParaExchangeCoupling[iExchange];
49  X->Large.is1_up = X->Def.Tpow[isite1 - 1];
50  X->Large.is2_up = X->Def.Tpow[isite2 - 1];
54  X->Large.isA_spin = X->Large.is1_up + X->Large.is2_up;
55  return 0;
56 }/*int child_exchange_spin_GetInfo*/
struct DefineList Def
Definision of system (Hamiltonian) etc.
Definition: struct.h:410
int ** ExchangeCoupling
[DefineList::NExchangeCoupling][2] Index of exchange term. malloc in setmem_def().
Definition: struct.h:146
long unsigned int is2_up
Mask used in the bit oeration.
Definition: struct.h:327
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.h:412
double complex tmp_J
Coupling constant.
Definition: struct.h:323
double * ParaExchangeCoupling
[DefineList::NExchangeCoupling] Coupling constant of exchange term. malloc in setmem_def().
Definition: struct.h:148
long unsigned int is1_up
Mask used in the bit oeration.
Definition: struct.h:325
long unsigned int * Tpow
[2 * DefineList::NsiteMPI] malloc in setmem_def().
Definition: struct.h:90
long unsigned int isA_spin
Mask used in the bit oeration.
Definition: struct.h:346
+ Here is the caller graph for this function:

◆ child_general_int_spin_GetInfo()

int child_general_int_spin_GetInfo ( struct BindStruct X,
long unsigned int  isite1,
long unsigned int  isite2,
long unsigned int  sigma1,
long unsigned int  sigma2,
long unsigned int  sigma3,
long unsigned int  sigma4,
double complex  tmp_V 
)

Set parameters for the bit operation of spin-general interaction term.

Returns
Always return 0
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)

Set the pairlift coupling constant (LargeList::tmp_J)

Set the bit mask for computing spin state of both site (LargeList::is1_up, LargeList::is2_up)

Set the bit mask for general interaction (LargeList::is1_spin, LargeList::is2_spin, LargeList::is3_spin, LargeList::is4_spin)

Parameters
[in,out]X
[in]isite1Site 1
[in]isite2Site 2
[in]sigma1Sigma 1, final state of site 1
[in]sigma2Sigma 3, initial state of site 1
[in]sigma3Sigma 3, final state of site 2
[in]sigma4Sigma 4, initial state of site 2
[in]tmp_VGeneral interaction constatnt

Definition at line 91 of file mltplySpinCore.c.

References BindStruct::Def, LargeList::is1_spin, LargeList::is1_up, LargeList::is2_spin, LargeList::is2_up, LargeList::is3_spin, LargeList::is4_spin, LargeList::isite1, LargeList::isite2, BindStruct::Large, LargeList::tmp_V, and DefineList::Tpow.

Referenced by makeHam(), mltplyHalfSpin(), and mltplyHalfSpinGC().

100  {
104  X->Large.tmp_V = tmp_V;
105  X->Large.isite1 = isite1;
106  X->Large.isite2 = isite2;
111  X->Large.is1_up = X->Def.Tpow[isite1 - 1];
112  X->Large.is2_up = X->Def.Tpow[isite2 - 1];
117  X->Large.is1_spin = sigma1;
118  X->Large.is2_spin = sigma2;
119  X->Large.is3_spin = sigma3;
120  X->Large.is4_spin = sigma4;
121  return 0;
122 }/*int child_general_int_spin_GetInfo*/
struct DefineList Def
Definision of system (Hamiltonian) etc.
Definition: struct.h:410
long unsigned int is2_up
Mask used in the bit oeration.
Definition: struct.h:327
long unsigned int is4_spin
Mask used in the bit oeration.
Definition: struct.h:335
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.h:412
long unsigned int is1_spin
Mask used in the bit oeration.
Definition: struct.h:332
int isite2
Is it realy used ???
Definition: struct.h:337
long unsigned int is3_spin
Mask used in the bit oeration.
Definition: struct.h:334
long unsigned int is1_up
Mask used in the bit oeration.
Definition: struct.h:325
long unsigned int * Tpow
[2 * DefineList::NsiteMPI] malloc in setmem_def().
Definition: struct.h:90
long unsigned int is2_spin
Mask used in the bit oeration.
Definition: struct.h:333
double complex tmp_V
Coupling constant.
Definition: struct.h:348
int isite1
Is it realy used ???
Definition: struct.h:336
+ Here is the caller graph for this function:

◆ child_pairlift_spin_GetInfo()

int child_pairlift_spin_GetInfo ( int  iPairLift,
struct BindStruct X 
)

Set parameters for the bit operation of spin-pairlift term.

Returns
Always return 0
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)

Set the pairlift coupling constant (LargeList::tmp_J)

Set the bit mask for computing spin state of both site (LargeList::is1_up, LargeList::is2_up)

Set the bit mask for exchange 2 spins (LargeList::isA_spin)

Definition at line 63 of file mltplySpinCore.c.

References BindStruct::Def, LargeList::is1_up, LargeList::is2_up, LargeList::isA_spin, BindStruct::Large, DefineList::PairLiftCoupling, DefineList::ParaPairLiftCoupling, LargeList::tmp_J, and DefineList::Tpow.

Referenced by makeHam(), and mltplyHalfSpinGC().

66  {
67  int isite1 = X->Def.PairLiftCoupling[iPairLift][0] + 1;
68  int isite2 = X->Def.PairLiftCoupling[iPairLift][1] + 1;
72  X->Large.tmp_J = X->Def.ParaPairLiftCoupling[iPairLift];
77  X->Large.is1_up = X->Def.Tpow[isite1 - 1];
78  X->Large.is2_up = X->Def.Tpow[isite2 - 1];
82  X->Large.isA_spin = X->Large.is1_up + X->Large.is2_up;
83  return 0;
84 }/*int child_pairlift_spin_GetInfo*/
struct DefineList Def
Definision of system (Hamiltonian) etc.
Definition: struct.h:410
long unsigned int is2_up
Mask used in the bit oeration.
Definition: struct.h:327
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.h:412
int ** PairLiftCoupling
[DefineList::NPairHopping][2] Index of pair-lift term. malloc in setmem_def().
Definition: struct.h:154
double complex tmp_J
Coupling constant.
Definition: struct.h:323
long unsigned int is1_up
Mask used in the bit oeration.
Definition: struct.h:325
double * ParaPairLiftCoupling
[DefineList::NPairHopping] Coupling constant of pair-lift term. malloc in setmem_def().
Definition: struct.h:156
long unsigned int * Tpow
[2 * DefineList::NsiteMPI] malloc in setmem_def().
Definition: struct.h:90
long unsigned int isA_spin
Mask used in the bit oeration.
Definition: struct.h:346
+ Here is the caller graph for this function:

◆ GC_child_CisAisCisAis_spin_element()

double complex GC_child_CisAisCisAis_spin_element ( long unsigned int  j,
long unsigned int  isA_up,
long unsigned int  isB_up,
long unsigned int  org_sigma2,
long unsigned int  org_sigma4,
double complex  tmp_V,
double complex *  tmp_v0,
double complex *  tmp_v1,
struct BindStruct X 
)

Compute \(c_{is}^\dagger c_{is} c_{is}^\dagger c_{is}\) term of grandcanonical spsin system.

Returns
Fragment of \(\langle v_1 | H_{\rm this}|v_1 \rangle\)
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]jIndex of initial wavefunction
[in]isA_upBit mask for spin 1
[in]isB_upBit mask for spin 2
[in]org_sigma2Target for spin 1
[in]org_sigma4Target for spin 2
[in]tmp_VCoupling constatnt
[in]tmp_v0Resulting wavefunction
[in]tmp_v1Wavefunction to be multiplied
[in,out]X

Definition at line 426 of file mltplySpinCore.c.

References BindStruct::Large, LargeList::mode, and X_SpinGC_CisAis().

Referenced by expec_cisajscktalt_SpinGCHalf(), GC_child_general_int_spin(), and makeHam().

436  {
437  int tmp_sgn;
438  double complex dmv = 0;
439  double complex dam_pr = 0;
440 
441  tmp_sgn = X_SpinGC_CisAis(j, X, isB_up, org_sigma4);
442  tmp_sgn *= X_SpinGC_CisAis(j, X, isA_up, org_sigma2);
443  if (tmp_sgn != 0) {
444  dmv = tmp_v1[j] * tmp_sgn * tmp_V;
445  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) { // for multply
446  tmp_v0[j] += dmv;
447  }
448  dam_pr = conj(tmp_v1[j]) * dmv;
449  }
450  return dam_pr;
451 }/*double complex GC_child_CisAisCisAis_spin_element*/
int X_SpinGC_CisAis(long unsigned int j, struct BindStruct *X, long unsigned int is1_spin, long unsigned int sigma1)
Compute the grandcanonical spin state with bit mask is1_spin.
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.h:412
int mode
multiply or expectation value.
Definition: struct.h:330
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GC_child_CisAisCitAiu_spin_element()

double complex GC_child_CisAisCitAiu_spin_element ( long unsigned int  j,
long unsigned int  org_sigma2,
long unsigned int  org_sigma4,
long unsigned int  isA_up,
long unsigned int  isB_up,
double complex  tmp_V,
double complex *  tmp_v0,
double complex *  tmp_v1,
struct BindStruct X,
long unsigned int *  tmp_off 
)

Compute \(c_{is}^\dagger c_{is} c_{it}^\dagger c_{iu}\) term of grandcanonical spsin system.

Returns
Fragment of \(\langle v_1 | H_{\rm this}|v_1 \rangle\)
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]jIndex of initial wavefunction
[in]org_sigma2Target for spin 1
[in]org_sigma4Target for spin 2
[in]isA_upBit mask for spin 1
[in]isB_upBit mask for spin 2
[in]tmp_VCoupling constatnt
[in]tmp_v0Resulting wavefunction
[in]tmp_v1Wavefunction to be multiplied
[in,out]X
[out]tmp_offIndex of final wavefunction

Definition at line 459 of file mltplySpinCore.c.

References BindStruct::Large, LargeList::mode, X_SpinGC_CisAis(), and X_SpinGC_CisAit().

Referenced by expec_cisajscktalt_SpinGCHalf(), GC_child_general_int_spin(), and makeHam().

470  {
471  int tmp_sgn;
472  double complex dmv;
473  double complex dam_pr = 0 + 0 * I;
474  tmp_sgn = X_SpinGC_CisAit(j, X, isB_up, org_sigma4, tmp_off);
475  if (tmp_sgn != 0) {
476  tmp_sgn *= X_SpinGC_CisAis((*tmp_off + 1), X, isA_up, org_sigma2);
477  if (tmp_sgn != 0) {
478  dmv = tmp_v1[j] * tmp_sgn * tmp_V;
479  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) { // for multply
480  tmp_v0[*tmp_off + 1] += dmv;
481  }
482  dam_pr = conj(tmp_v1[*tmp_off + 1]) * dmv;
483  }/*if (tmp_sgn != 0)*/
484  }/*if (tmp_sgn != 0)*/
485  return dam_pr;
486 }/*double complex GC_child_CisAisCitAiu_spin_element*/
int X_SpinGC_CisAis(long unsigned int j, struct BindStruct *X, long unsigned int is1_spin, long unsigned int sigma1)
Compute the grandcanonical spin state with bit mask is1_spin.
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.h:412
int mode
multiply or expectation value.
Definition: struct.h:330
int X_SpinGC_CisAit(long unsigned int j, struct BindStruct *X, long unsigned int is1_spin, long unsigned int sigma2, long unsigned int *tmp_off)
Compute index of final wavefunction by term (grandcanonical).
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GC_child_CisAitCiuAiu_spin_element()

double complex GC_child_CisAitCiuAiu_spin_element ( long unsigned int  j,
long unsigned int  org_sigma2,
long unsigned int  org_sigma4,
long unsigned int  isA_up,
long unsigned int  isB_up,
double complex  tmp_V,
double complex *  tmp_v0,
double complex *  tmp_v1,
struct BindStruct X,
long unsigned int *  tmp_off 
)

Compute \(c_{is}^\dagger c_{it} c_{iu}^\dagger c_{iu}\) term of grandcanonical spsin system.

Returns
Fragment of \(\langle v_1 | H_{\rm this}|v_1 \rangle\)
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]jIndex of initial wavefunction
[in]org_sigma2Target for spin 1
[in]org_sigma4Target for spin 2
[in]isA_upBit mask for spin 1
[in]isB_upBit mask for spin 2
[in]tmp_VCoupling constatnt
[in]tmp_v0Resulting wavefunction
[in]tmp_v1Wavefunction to be multiplied
[in,out]X
[out]tmp_offIndex of final wavefunction

Definition at line 494 of file mltplySpinCore.c.

References BindStruct::Large, LargeList::mode, X_SpinGC_CisAis(), and X_SpinGC_CisAit().

Referenced by expec_cisajscktalt_SpinGCHalf(), GC_child_general_int_spin(), and makeHam().

505  {
506  int tmp_sgn;
507  double complex dmv;
508  double complex dam_pr = 0 + 0 * I;
509  tmp_sgn = X_SpinGC_CisAis(j, X, isB_up, org_sigma4);
510  if (tmp_sgn != 0) {
511  tmp_sgn *= X_SpinGC_CisAit(j, X, isA_up, org_sigma2, tmp_off);
512  if (tmp_sgn != 0) {
513  dmv = tmp_v1[j] * tmp_sgn * tmp_V;
514  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) { // for multply
515  tmp_v0[*tmp_off + 1] += dmv;
516  }
517  dam_pr = conj(tmp_v1[*tmp_off + 1]) * dmv;
518  }/*if (tmp_sgn != 0)*/
519  }/*if (tmp_sgn != 0)*/
520  return dam_pr;
521 }/*double complex GC_child_CisAitCiuAiu_spin_element*/
int X_SpinGC_CisAis(long unsigned int j, struct BindStruct *X, long unsigned int is1_spin, long unsigned int sigma1)
Compute the grandcanonical spin state with bit mask is1_spin.
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.h:412
int mode
multiply or expectation value.
Definition: struct.h:330
int X_SpinGC_CisAit(long unsigned int j, struct BindStruct *X, long unsigned int is1_spin, long unsigned int sigma2, long unsigned int *tmp_off)
Compute index of final wavefunction by term (grandcanonical).
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GC_child_CisAitCiuAiv_spin_element()

double complex GC_child_CisAitCiuAiv_spin_element ( long unsigned int  j,
long unsigned int  org_sigma2,
long unsigned int  org_sigma4,
long unsigned int  isA_up,
long unsigned int  isB_up,
double complex  tmp_V,
double complex *  tmp_v0,
double complex *  tmp_v1,
struct BindStruct X,
long unsigned int *  tmp_off_2 
)

Compute \(c_{is}^\dagger c_{it} c_{iu}^\dagger c_{iv}\) term of grandcanonical spsin system.

Returns
Fragment of \(\langle v_1 | H_{\rm this}|v_1 \rangle\)
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]jIndex of initial wavefunction
[in]org_sigma2Target for spin 1
[in]org_sigma4Target for spin 2
[in]isA_upBit mask for spin 1
[in]isB_upBit mask for spin 2
[in]tmp_VCoupling constatnt
[in]tmp_v0Resulting wavefunction
[in]tmp_v1Wavefunction to be multiplied
[in,out]X
[out]tmp_off_2Index of final wavefunction

Definition at line 529 of file mltplySpinCore.c.

References BindStruct::Large, LargeList::mode, and X_SpinGC_CisAit().

Referenced by expec_cisajscktalt_SpinGCHalf(), GC_child_general_int_spin(), and makeHam().

540  {
541  int tmp_sgn;
542  long unsigned int tmp_off_1;
543  double complex dmv;
544  double complex dam_pr = 0 + 0 * I;
545  tmp_sgn = X_SpinGC_CisAit(j, X, isB_up, org_sigma4, &tmp_off_1);
546  if (tmp_sgn != 0) {
547  tmp_sgn *= X_SpinGC_CisAit((tmp_off_1 + 1), X, isA_up, org_sigma2, tmp_off_2);
548  if (tmp_sgn != 0) {
549  dmv = tmp_v1[j] * tmp_sgn * tmp_V;
550  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) { // for multply
551  tmp_v0[*tmp_off_2 + 1] += dmv;
552  }
553  dam_pr = conj(tmp_v1[*tmp_off_2 + 1]) * dmv;
554  }/*if (tmp_sgn != 0)*/
555  }/*if (tmp_sgn != 0)*/
556  return dam_pr;
557 }/*double complex GC_child_CisAitCiuAiv_spin_element*/
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.h:412
int mode
multiply or expectation value.
Definition: struct.h:330
int X_SpinGC_CisAit(long unsigned int j, struct BindStruct *X, long unsigned int is1_spin, long unsigned int sigma2, long unsigned int *tmp_off)
Compute index of final wavefunction by term (grandcanonical).
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GC_child_exchange_spin_element()

double complex GC_child_exchange_spin_element ( long unsigned int  j,
double complex *  tmp_v0,
double complex *  tmp_v1,
struct BindStruct X,
long unsigned int *  tmp_off 
)

Multiply Hamiltonian of exchange term of grandcanonical spin system.

Returns
\(\langle v_1 | H_{\rm this}| v_1 \rangle\)
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]jIndex of initial wavefunction
[out]tmp_v0Resulting wavefunction
[in]tmp_v1Wavefunction to be multiplied
[in,out]X
[out]tmp_offIndex of final wavefunction

Definition at line 314 of file mltplySpinCore.c.

References LargeList::isA_spin, BindStruct::Large, LargeList::mode, and LargeList::tmp_J.

Referenced by GC_child_exchange_spin(), and makeHam().

320  {
321  double complex dmv;
322  long unsigned int is_up = X->Large.isA_spin;
323  double complex tmp_J = X->Large.tmp_J;
324  int mode = X->Large.mode;
325  long unsigned int list_1_j, list_1_off;
326 
327  double complex dam_pr = 0;
328  list_1_j = j - 1;
329 
330  long unsigned int ibit_tmp;
331  ibit_tmp = (list_1_j & is_up);
332  if (ibit_tmp == 0 || ibit_tmp == is_up) {
333  return dam_pr;
334  }
335  else {
336  list_1_off = list_1_j ^ is_up;
337  *tmp_off = list_1_off;
338  dmv = tmp_J * tmp_v1[j];
339  if (mode == M_MLTPLY) {
340  tmp_v0[list_1_off + 1] += dmv;
341  }
342  dam_pr += dmv * conj(tmp_v1[list_1_off + 1]);
343  return dam_pr;
344  }
345 }/*double complex GC_child_exchange_spin_element*/
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.h:412
int mode
multiply or expectation value.
Definition: struct.h:330
double complex tmp_J
Coupling constant.
Definition: struct.h:323
long unsigned int isA_spin
Mask used in the bit oeration.
Definition: struct.h:346
+ Here is the caller graph for this function:

◆ GC_child_pairlift_spin_element()

double complex GC_child_pairlift_spin_element ( long unsigned int  j,
double complex *  tmp_v0,
double complex *  tmp_v1,
struct BindStruct X,
long unsigned int *  tmp_off 
)

Multiply Hamiltonian of pairlift term of grandcanonical spin system.

Returns
\(\langle v_1 | H_{\rm this}| v_1 \rangle\)
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]jIndex of initial wavefunction
[out]tmp_v0Resulting wavefunction
[in]tmp_v1Wavefunction to be multiplied
[in,out]X
[out]tmp_offIndex of final wavefunction

Definition at line 352 of file mltplySpinCore.c.

References LargeList::isA_spin, BindStruct::Large, LargeList::mode, and LargeList::tmp_J.

Referenced by GC_child_pairlift_spin(), and makeHam().

358  {
359  double complex dmv;
360  long unsigned int is_up = X->Large.isA_spin;
361  double complex tmp_J = X->Large.tmp_J;
362  int mode = X->Large.mode;
363  double complex dam_pr = 0;
364  long unsigned int list_1_off;
365  long unsigned int list_1_j = j - 1;
366  long unsigned int ibit_tmp;
367  //ibit_tmp = ((list_1_j & is1_up) / is1_up) ^ ((list_1_j & is2_up) / is2_up);
368  ibit_tmp = (list_1_j & is_up);
369  if (ibit_tmp == 0 || ibit_tmp == is_up) {
370  list_1_off = list_1_j ^ is_up; //Change: ++ -> -- or -- -> ++
371  *tmp_off = list_1_off;
372  dmv = tmp_J * tmp_v1[j];//* ibit_tmp;
373  if (mode == M_MLTPLY) {
374  tmp_v0[list_1_off + 1] += dmv;
375  }
376  dam_pr += dmv * conj(tmp_v1[list_1_off + 1]);
377  return dam_pr;
378  }
379  else {
380  return dam_pr;
381  }
382 }/*double complex GC_child_pairlift_spin_element*/
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.h:412
int mode
multiply or expectation value.
Definition: struct.h:330
double complex tmp_J
Coupling constant.
Definition: struct.h:323
long unsigned int isA_spin
Mask used in the bit oeration.
Definition: struct.h:346
+ Here is the caller graph for this function:

◆ X_child_exchange_spin_element()

int X_child_exchange_spin_element ( long unsigned int  j,
struct BindStruct X,
long unsigned int  isA_up,
long unsigned int  isB_up,
long unsigned int  sigmaA,
long unsigned int  sigmaB,
long unsigned int *  tmp_off 
)

Compute index of final wavefunction associated to spin-exchange term.

Returns
1 if spin of site 1 is sigmaA and spin of site 2 is sigmaB. 0 if not.
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]jIndex of initial wavefunction
[in,out]X
[in]isA_upBit mask for spin 1
[in]isB_upBit mask for spin 2
[in]sigmaATarget of spin 1
[in]sigmaBTarget of spin 2
[out]tmp_offIndex of final wavefunction

Definition at line 239 of file mltplySpinCore.c.

References GetOffComp(), LargeList::ihfbit, LargeList::ilft, LargeList::irght, BindStruct::Large, list_1, list_2_1, and list_2_2.

Referenced by child_general_int_spin(), expec_cisajscktalt_SpinHalf(), and makeHam().

247  {
248  long unsigned int iexchg, off;
249  long unsigned int irght = X->Large.irght;
250  long unsigned int ilft = X->Large.ilft;
251  long unsigned int ihfbit = X->Large.ihfbit;
252  long unsigned int ibit_tmp_A, ibit_tmp_B;
253 
254  ibit_tmp_A = ((list_1[j] & isA_up) / isA_up);
255  ibit_tmp_B = ((list_1[j] & isB_up) / isB_up);
256  if (ibit_tmp_A == sigmaA && ibit_tmp_B == sigmaB) {
257  iexchg = list_1[j] ^ (isA_up + isB_up);
258  GetOffComp(list_2_1, list_2_2, iexchg, irght, ilft, ihfbit, &off);
259  *tmp_off = off;
260  return 1;
261  }
262  else {
263  *tmp_off = 0; // just tentative
264  return 0;
265  }
266 }/*int X_child_exchange_spin_element*/
long unsigned int ihfbit
Used for Ogata-Lin ???
Definition: struct.h:345
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.h:412
long unsigned int ilft
Used for Ogata-Lin ???
Definition: struct.h:344
long unsigned int * list_2_1
Definition: global.h:49
int GetOffComp(long unsigned int *_list_2_1, long unsigned int *_list_2_2, long unsigned int _ibit, const long unsigned int _irght, const long unsigned int _ilft, const long unsigned int _ihfbit, long unsigned int *_ioffComp)
function of getting off-diagonal component
Definition: bitcalc.c:195
long unsigned int * list_1
Definition: global.h:47
long unsigned int * list_2_2
Definition: global.h:50
long unsigned int irght
Used for Ogata-Lin ???
Definition: struct.h:343
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ X_Spin_CisAis()

int X_Spin_CisAis ( long unsigned int  j,
struct BindStruct X,
long unsigned int  is1_spin,
long unsigned int  sigma1 
)

Compute the spin state with bit mask is1_spin.

Returns
1 if the spin state with bit mask is1_spin is sigma1, 0 for the other.
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]jIndex of wavefunction
[in,out]X
[in]is1_spinBit mask
[in]sigma1Target spin state

Definition at line 166 of file mltplySpinCore.c.

References list_1.

Referenced by child_CisAisCisAis_spin_element(), expec_cisajs_SpinHalf(), expec_cisajscktalt_SpinHalf(), GetPairExcitedStateHalfSpin(), makeHam(), SetDiagonalInterAll(), and SetDiagonalTEInterAll().

171  {
172  int A_ibit_tmp;
173  // off = j
174  A_ibit_tmp = ((list_1[j] & is1_spin) / is1_spin) ^ (1 - sigma1);
175  return A_ibit_tmp;
176 }/*int X_Spin_CisAis*/
long unsigned int * list_1
Definition: global.h:47
+ Here is the caller graph for this function:

◆ X_Spin_CisAit()

int X_Spin_CisAit ( long unsigned int  j,
struct BindStruct X,
long unsigned int  is1_spin,
long unsigned int  sigma2,
long unsigned int *  list_1_Org_,
long unsigned int *  list_2_1_,
long unsigned int *  list_2_2_,
long unsigned int *  tmp_off 
)

Compute index of final wavefunction by \(c_{is}^\dagger c_{it}\) term.

Returns
1 if bit-mask is1_spin is sigma2, 0 for the other
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]jIndex of initial wavefunction
[in,out]X
[in]is1_spinBit mask for computing spin state
[in]sigma2Spin state at site 2
[in]list_1_Org_Similar to list_1
[in]list_2_1_Similar to list_2_1
[in]list_2_2_Similar to list_2_2
[out]tmp_offIndex of final wavefunction

Definition at line 138 of file mltplySpinCore.c.

References GetOffComp(), LargeList::ihfbit, LargeList::ilft, LargeList::irght, BindStruct::Large, and X_SpinGC_CisAit().

Referenced by GetPairExcitedStateHalfSpin().

147  {
148  long unsigned int list_1_j;
149  long unsigned int off;
150  list_1_j = list_1_Org_[j];
151  if (X_SpinGC_CisAit(list_1_j + 1, X, is1_spin, sigma2, &off) != 0) {
152  GetOffComp(list_2_1_, list_2_2_, off, X->Large.irght, X->Large.ilft, X->Large.ihfbit, tmp_off);
153  return 1;
154  }
155  else {
156  *tmp_off = 0;
157  return 0;
158  }
159 }/*int X_Spin_CisAit*/
long unsigned int ihfbit
Used for Ogata-Lin ???
Definition: struct.h:345
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.h:412
long unsigned int ilft
Used for Ogata-Lin ???
Definition: struct.h:344
int GetOffComp(long unsigned int *_list_2_1, long unsigned int *_list_2_2, long unsigned int _ibit, const long unsigned int _irght, const long unsigned int _ilft, const long unsigned int _ihfbit, long unsigned int *_ioffComp)
function of getting off-diagonal component
Definition: bitcalc.c:195
long unsigned int irght
Used for Ogata-Lin ???
Definition: struct.h:343
int X_SpinGC_CisAit(long unsigned int j, struct BindStruct *X, long unsigned int is1_spin, long unsigned int sigma2, long unsigned int *tmp_off)
Compute index of final wavefunction by term (grandcanonical).
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ X_SpinGC_CisAis()

int X_SpinGC_CisAis ( long unsigned int  j,
struct BindStruct X,
long unsigned int  is1_spin,
long unsigned int  sigma1 
)

Compute the grandcanonical spin state with bit mask is1_spin.

Returns
1 if the spin state with bit mask is1_spin is sigma1, 0 for the other.
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]jIndex of wavefunction
[in,out]X
[in]is1_spinBit mask
[in]sigma1Target spin state

Definition at line 183 of file mltplySpinCore.c.

Referenced by expec_cisajs_SpinGCHalf(), expec_cisajs_SpinHalf(), expec_cisajscktalt_SpinHalf(), GC_child_CisAisCisAis_spin_element(), GC_child_CisAisCitAiu_spin_element(), GC_child_CisAitCiuAiu_spin_element(), GetPairExcitedStateHalfSpin(), GetPairExcitedStateHalfSpinGC(), SetDiagonalInterAll(), SetDiagonalTEInterAll(), totalspin_Spin(), totalspin_SpinGC(), X_GC_child_CisAisCjuAju_spin_MPIdouble(), X_GC_child_CisAisCjuAju_spin_MPIsingle(), X_GC_child_CisAisCjuAjv_spin_MPIdouble(), and X_GC_child_CisAitCjuAju_spin_MPIdouble().

188  {
189  int A_ibit_tmp;
190  long unsigned int list_1_j;
191  // off = j
192  list_1_j = j - 1;
193  A_ibit_tmp = ((list_1_j & is1_spin) / is1_spin) ^ (1 - sigma1);
194  return A_ibit_tmp;
195 }/*int X_SpinGC_CisAis*/
+ Here is the caller graph for this function:

◆ X_SpinGC_CisAit()

int X_SpinGC_CisAit ( long unsigned int  j,
struct BindStruct X,
long unsigned int  is1_spin,
long unsigned int  sigma2,
long unsigned int *  tmp_off 
)

Compute index of final wavefunction by \(c_{is}^\dagger c_{it}\) term (grandcanonical).

Returns
1 if bit-mask is1_spin is sigma2, 0 for the other
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]jIndex of wavefunction
[in,out]X
[in]is1_spinBit mask for computing spin state
[in]sigma2Spin state at site 2
[out]tmp_offIndex of final wavefunction

Definition at line 203 of file mltplySpinCore.c.

Referenced by expec_cisajs_SpinGCHalf(), GC_child_CisAisCitAiu_spin_element(), GC_child_CisAitCiuAiu_spin_element(), GC_child_CisAitCiuAiv_spin_element(), GetPairExcitedStateHalfSpinGC(), makeHam(), mltplyHalfSpinGC(), X_GC_child_CisAitCiuAiv_spin_MPIsingle(), and X_Spin_CisAit().

209  {
210  long unsigned int list_1_j, ibit_tmp_1;
211 
212  list_1_j = j - 1;
213 
214  ibit_tmp_1 = list_1_j & is1_spin;
215  if (ibit_tmp_1 == 0 && sigma2 == 0) { // down -> up
216  *tmp_off = list_1_j + is1_spin;
217  return 1;
218  }
219  else if (ibit_tmp_1 != 0 && sigma2 == 1) { // up -> down
220  *tmp_off = list_1_j - is1_spin;
221  return 1;
222  }
223  else {
224  *tmp_off = 0;
225  return 0;
226  }
227 }/*int X_SpinGC_CisAit*/
+ Here is the caller graph for this function: