HΦ  3.1.0
mltplyMPISpinCore.c File Reference

Functions for spin Hamiltonian + MPI (Core) More...

#include "mpi.h"
#include "Common.h"
#include "mltplyCommon.h"
#include "mltplySpinCore.h"
#include "mltplyMPISpinCore.h"
#include "bitcalc.h"
#include "wrapperMPI.h"

Go to the source code of this file.

Functions

void GC_child_CisAitCiuAiv_spin_MPIdouble (unsigned long int i_int, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Exchange and Pairlifting term in Spin model + GC When both site1 and site2 are in the inter process region. More...
 
double complex X_GC_child_CisAitCiuAiv_spin_MPIdouble (int org_isite1, int org_ispin1, int org_ispin2, int org_isite3, int org_ispin3, int org_ispin4, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 \(c_{is}^\dagger c_{it} c_{iu}^\dagger c_{iv}\) term in Spin model + GC. When both site1 and site2 are in the inter process region. More...
 
void GC_child_CisAisCjuAjv_spin_MPIdouble (unsigned long int i_int, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Wrapper for calculating CisAisCjuAjv term in Spin model + GC When both site1 and site2 are in the inter process region. More...
 
double complex X_GC_child_CisAisCjuAjv_spin_MPIdouble (int org_isite1, int org_ispin1, int org_isite3, int org_ispin3, int org_ispin4, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 CisAisCjuAjv term in Spin model + GC When both site1 and site2 are in the inter process region. More...
 
void GC_child_CisAitCjuAju_spin_MPIdouble (unsigned long int i_int, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Wrapper for calculating CisAitCjuAju term in Spin model + GC When both site1 and site2 are in the inter process region. More...
 
double complex X_GC_child_CisAitCjuAju_spin_MPIdouble (int org_isite1, int org_ispin1, int org_ispin2, int org_isite3, int org_ispin3, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 CisAisCjuAjv term in Spin model + GC When both site1 and site2 are in the inter process region. More...
 
double complex X_GC_child_CisAisCjuAju_spin_MPIdouble (int org_isite1, int org_ispin1, int org_isite3, int org_ispin3, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 CisAisCjuAjv term in Spin model + GC When both site1 and site2 are in the inter process region. More...
 
double complex X_GC_child_CisAisCjuAju_spin_MPIsingle (int org_isite1, int org_ispin1, int org_isite3, int org_ispin3, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 CisAisCjuAjv term in Spin model + GC When both site1 and site2 are in the inter process region. More...
 
void GC_child_CisAitCiuAiv_spin_MPIsingle (unsigned long int i_int, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Exchange and Pairlifting term in Spin model + GC When only site2 is in the inter process region. More...
 
double complex X_GC_child_CisAitCiuAiv_spin_MPIsingle (int org_isite1, int org_ispin1, int org_ispin2, int org_isite3, int org_ispin3, int org_ispin4, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Exchange and Pairlifting term in Spin model + GC When only site2 is in the inter process region. More...
 
void GC_child_CisAisCjuAjv_spin_MPIsingle (unsigned long int i_int, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Wrapper for CisAisCjuAjv term in Spin model + GC When only site2 is in the inter process region. More...
 
double complex X_GC_child_CisAisCjuAjv_spin_MPIsingle (int org_isite1, int org_ispin1, int org_isite3, int org_ispin3, int org_ispin4, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 CisAisCjuAjv term in Spin model + GC When only site2 is in the inter process region. More...
 
void GC_child_CisAitCjuAju_spin_MPIsingle (unsigned long int i_int, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Wrapper for CisAisCjuAjv term in Spin model + GC When only site2 is in the inter process region. More...
 
double complex X_GC_child_CisAitCjuAju_spin_MPIsingle (int org_isite1, int org_ispin1, int org_ispin2, int org_isite3, int org_ispin3, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 CisAisCjuAjv term in Spin model + GC When only site2 is in the inter process region. More...
 
double complex X_GC_child_CisAisCjuAjv_GeneralSpin_MPIdouble (int org_isite1, int org_ispin1, int org_isite3, int org_ispin3, int org_ispin4, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 \(c_{is}^\dagger c_{is} c_{ju}^\dagger c_{jv}\) term in Spin model. When both site1 and site3 are in the inter process region. More...
 
double complex X_GC_child_CisAitCjuAju_GeneralSpin_MPIdouble (int org_isite1, int org_ispin1, int org_ispin2, int org_isite3, int org_ispin3, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 \(c_{is}^\dagger c_{it} c_{ju}^\dagger c_{ju}\) term in Spin model. When both site1 and site3 are in the inter process region. More...
 
double complex X_GC_child_CisAitCjuAjv_GeneralSpin_MPIdouble (int org_isite1, int org_ispin1, int org_ispin2, int org_isite3, int org_ispin3, int org_ispin4, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Compute \(c_{is}^\dagger c_{it} c_{ju}^\dagger c_{jv}\) term in the grandcanonical general spin system when both site is in the inter process region. More...
 
double complex X_GC_child_CisAisCjuAju_GeneralSpin_MPIdouble (int org_isite1, int org_ispin1, int org_isite3, int org_ispin3, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Compute \(c_{is}^\dagger c_{is} c_{ju}^\dagger c_{ju}\) term in the grandcanonical general spin system when both site is in the inter process region. More...
 
double complex X_GC_child_CisAit_GeneralSpin_MPIdouble (int org_isite1, int org_ispin1, int org_ispin2, double complex tmp_trans, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Compute \(c_{is}^\dagger c_{it}\) term in the grandcanonical general spin system when both site is in the inter process region. More...
 
double complex X_GC_child_CisAis_GeneralSpin_MPIdouble (int org_isite1, int org_ispin1, double complex tmp_trans, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Compute \(c_{is}^\dagger c_{is}\) term in the grandcanonical general spin system when both site is in the inter process region. More...
 
double complex X_GC_child_AisCis_GeneralSpin_MPIdouble (int org_isite1, int org_ispin1, double complex tmp_trans, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Compute \(c_{is} c_{is}^\dagger\) term in the grandcanonical general spin system when both site is in the inter process region. More...
 
double complex X_child_CisAit_GeneralSpin_MPIdouble (int org_isite1, int org_ispin1, int org_ispin2, double complex tmp_trans, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1, double complex *tmp_v1buf, unsigned long int idim_max, long unsigned int *list_1_org, long unsigned int *list_1buf_org, long unsigned int _ihfbit)
 Compute \(c_{is}^\dagger c_{it}\) term in the canonical general spin system when both site is in the inter process region. More...
 
double complex X_GC_child_CisAisCjuAjv_GeneralSpin_MPIsingle (int org_isite1, int org_ispin1, int org_isite3, int org_ispin3, int org_ispin4, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Compute \(c_{is}^\dagger c_{is}c_{ju}^\dagger c_{jv}\) term in the grandcanonical general spin system when one of these site is in the inter process region. More...
 
double complex X_GC_child_CisAitCjuAju_GeneralSpin_MPIsingle (int org_isite1, int org_ispin1, int org_ispin2, int org_isite3, int org_ispin3, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Compute \(c_{is}^\dagger c_{it}c_{ju}^\dagger c_{ju}\) term in the grandcanonical general spin system when one of these site is in the inter process region. More...
 
double complex X_GC_child_CisAitCjuAjv_GeneralSpin_MPIsingle (int org_isite1, int org_ispin1, int org_ispin2, int org_isite3, int org_ispin3, int org_ispin4, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Compute \(c_{is}^\dagger c_{is}c_{ju}^\dagger c_{jv}\) term in the grandcanonical general spin system when one of these site is in the inter process region. More...
 
double complex X_GC_child_CisAisCjuAju_GeneralSpin_MPIsingle (int org_isite1, int org_ispin1, int org_isite3, int org_ispin3, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Compute \(c_{is}^\dagger c_{is}c_{ju}^\dagger c_{ju}\) term in the grandcanonical general spin system when one of these site is in the inter process region. More...
 
double complex X_child_CisAitCjuAjv_GeneralSpin_MPIdouble (int org_isite1, int org_ispin1, int org_ispin2, int org_isite3, int org_ispin3, int org_ispin4, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Compute \(c_{is}^\dagger c_{it}c_{ju}^\dagger c_{jv}\) term in the canonical general spin system when both sites are in the inter process region. More...
 
double complex X_child_CisAisCjuAju_GeneralSpin_MPIdouble (int org_isite1, int org_ispin1, int org_isite3, int org_ispin3, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Compute \(c_{is}^\dagger c_{is}c_{ju}^\dagger c_{ju}\) term in the canonical general spin system when both sites are in the inter process region. More...
 
double complex X_child_CisAisCjuAju_GeneralSpin_MPIsingle (int org_isite1, int org_ispin1, int org_isite3, int org_ispin3, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Compute \(c_{is}^\dagger c_{is}c_{ju}^\dagger c_{ju}\) term in the canonical general spin system when one of these sites is in the inter process region. More...
 
double complex X_child_CisAitCjuAjv_GeneralSpin_MPIsingle (int org_isite1, int org_ispin1, int org_ispin2, int org_isite3, int org_ispin3, int org_ispin4, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Compute \(c_{is}^\dagger c_{it}c_{ju}^\dagger c_{jv}\) term in the canonical general spin system when one of these sites is in the inter process region. More...
 
double complex X_GC_child_CisAit_spin_MPIdouble (int org_isite1, int org_ispin1, int org_ispin2, double complex tmp_trans, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Hopping term in Spin + GC When both site1 and site2 are in the inter process region. More...
 
double complex X_child_CisAit_spin_MPIdouble (int org_isite1, int org_ispin2, double complex tmp_trans, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1, double complex *tmp_v1buf, unsigned long int idim_max, long unsigned int *Tpow, long unsigned int *list_1_org, long unsigned int *list_1buf_org, long unsigned int *list_2_1_target, long unsigned int *list_2_2_target, long unsigned int _irght, long unsigned int _ilft, long unsigned int _ihfbit)
 Hopping term in Spin + Canonical for CalcSpectrum When both site1 and site2 are in the inter process region. More...
 
double complex X_GC_child_CisAis_spin_MPIdouble (int org_isite1, int org_ispin1, double complex tmp_trans, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Hopping term in Spin + GC When both site1 and site2 are in the inter process region. More...
 
double complex X_GC_child_AisCis_spin_MPIdouble (int org_isite1, int org_ispin1, double complex tmp_trans, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Hopping term in Spin + GC When both site1 and site2 are in the inter process region. More...
 

Detailed Description

Functions for spin Hamiltonian + MPI (Core)

General two body term:

1/2 spin 1/2 spin
MPI single MPI double MPI single MPI double
\(c_{is}^\dagger c_{is} c_{ju}^\dagger c_{ju}\) ::GC_child_CisAisCjuAju_spin_MPIsingle, X_GC_child_CisAisCjuAjv_spin_MPIsingle ::GC_child_CisAisCjuAju_spin_MPIdouble, X_GC_child_CisAisCjuAjv_spin_MPIdouble X_child_CisAisCjuAju_GeneralSpin_MPIsingle, X_GC_child_CisAisCjuAjv_GeneralSpin_MPIsingle X_child_CisAisCjuAju_GeneralSpin_MPIsingle, X_GC_child_CisAisCjuAjv_GeneralSpin_MPIsingle
\(c_{is}^\dagger c_{is} c_{ju}^\dagger c_{jv}\) GC_child_CisAisCjuAjv_spin_MPIsingle, X_GC_child_CisAisCjuAjv_spin_MPIsingle GC_child_CisAisCjuAjv_spin_MPIdouble, X_GC_child_CisAisCjuAjv_spin_MPIdouble ::X_child_CisAisCjuAjv_GeneralSpin_MPIsingle, X_GC_child_CisAisCjuAjv_GeneralSpin_MPIsingle ::X_child_CisAisCjuAjv_GeneralSpin_MPIsingle, X_GC_child_CisAisCjuAjv_GeneralSpin_MPIsingle
\(c_{is}^\dagger c_{it} c_{ju}^\dagger c_{ju}\) GC_child_CisAitCjuAju_spin_MPIsingle, X_GC_child_CisAisCjuAjv_spin_MPIsingle GC_child_CisAitCjuAju_spin_MPIdouble, X_GC_child_CisAisCjuAjv_spin_MPIdouble ::X_child_CisAitCjuAju_GeneralSpin_MPIsingle, X_GC_child_CisAitCjuAju_GeneralSpin_MPIsingle ::X_child_CisAitCjuAju_GeneralSpin_MPIsingle, X_GC_child_CisAitCjuAju_GeneralSpin_MPIsingle
\(c_{is}^\dagger c_{it} c_{ju}^\dagger c_{jv}\) ::GC_child_CisAitCjuAjv_spin_MPIsingle, X_GC_child_CisAisCjuAjv_spin_MPIsingle ::GC_child_CisAitCjuAjv_spin_MPIdouble, X_GC_child_CisAisCjuAjv_spin_MPIdouble X_child_CisAitCjuAjv_GeneralSpin_MPIsingle, X_GC_child_CisAitCjuAjv_GeneralSpin_MPIsingle X_child_CisAitCjuAjv_GeneralSpin_MPIsingle, X_GC_child_CisAitCjuAjv_GeneralSpin_MPIsingle

Definition in file mltplyMPISpinCore.c.

Function Documentation

◆ GC_child_CisAisCjuAjv_spin_MPIdouble()

void GC_child_CisAisCjuAjv_spin_MPIdouble ( unsigned long int  i_int,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Wrapper for calculating CisAisCjuAjv term in Spin model + GC When both site1 and site2 are in the inter process region.

Author
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]i_intInteraction ID
[in,out]X
[out]tmp_v0Result v0 = H v1
[in]tmp_v1v0 = H v1

Definition at line 191 of file mltplyMPISpinCore.c.

References X, and X_GC_child_CisAisCjuAjv_spin_MPIdouble().

Referenced by GC_child_general_int_spin_MPIdouble().

196  {
197 #ifdef MPI
198  double complex dam_pr;
200  X->Def.InterAll_OffDiagonal[i_int][0], X->Def.InterAll_OffDiagonal[i_int][1],
201  X->Def.InterAll_OffDiagonal[i_int][4], X->Def.InterAll_OffDiagonal[i_int][5],
202  X->Def.InterAll_OffDiagonal[i_int][7], X->Def.ParaInterAll_OffDiagonal[i_int], X, tmp_v0, tmp_v1);
203  X->Large.prdct += dam_pr;
204 #endif
205 }/*void GC_child_CisAitCiuAiv_spin_MPIdouble*/
double complex X_GC_child_CisAisCjuAjv_spin_MPIdouble(int org_isite1, int org_ispin1, int org_isite3, int org_ispin3, int org_ispin4, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
CisAisCjuAjv term in Spin model + GC When both site1 and site2 are in the inter process region...
struct EDMainCalStruct X
Definition: struct.h:431

◆ GC_child_CisAisCjuAjv_spin_MPIsingle()

void GC_child_CisAisCjuAjv_spin_MPIsingle ( unsigned long int  i_int,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Wrapper for CisAisCjuAjv term in Spin model + GC When only site2 is in the inter process region.

Author
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]i_intInteraction ID
[in,out]X
[out]tmp_v0Result v0 = H v1
[in]tmp_v1v0 = H v1

Definition at line 608 of file mltplyMPISpinCore.c.

References X, and X_GC_child_CisAisCjuAjv_spin_MPIsingle().

Referenced by GC_child_general_int_spin_MPIsingle().

613  {
614 #ifdef MPI
615  double complex dam_pr;
617  X->Def.InterAll_OffDiagonal[i_int][0], X->Def.InterAll_OffDiagonal[i_int][1],
618  X->Def.InterAll_OffDiagonal[i_int][4], X->Def.InterAll_OffDiagonal[i_int][5],
619  X->Def.InterAll_OffDiagonal[i_int][7], X->Def.ParaInterAll_OffDiagonal[i_int], X, tmp_v0, tmp_v1);
620  X->Large.prdct += dam_pr;
621 #endif
622 }/*void GC_child_CisAisCjuAjv_spin_MPIsingle*/
struct EDMainCalStruct X
Definition: struct.h:431
double complex X_GC_child_CisAisCjuAjv_spin_MPIsingle(int org_isite1, int org_ispin1, int org_isite3, int org_ispin3, int org_ispin4, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
CisAisCjuAjv term in Spin model + GC When only site2 is in the inter process region.

◆ GC_child_CisAitCiuAiv_spin_MPIdouble()

void GC_child_CisAitCiuAiv_spin_MPIdouble ( unsigned long int  i_int,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Exchange and Pairlifting term in Spin model + GC When both site1 and site2 are in the inter process region.

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]i_intInteraction ID
[in,out]X
[out]tmp_v0Result v0 = H v1
[in]tmp_v1v0 = H v1

Definition at line 79 of file mltplyMPISpinCore.c.

References X, and X_GC_child_CisAitCiuAiv_spin_MPIdouble().

Referenced by GC_child_general_int_spin_MPIdouble().

84 {
85 #ifdef MPI
86  double complex dam_pr;
88  X->Def.InterAll_OffDiagonal[i_int][0], X->Def.InterAll_OffDiagonal[i_int][1],
89  X->Def.InterAll_OffDiagonal[i_int][3], X->Def.InterAll_OffDiagonal[i_int][4],
90  X->Def.InterAll_OffDiagonal[i_int][5], X->Def.InterAll_OffDiagonal[i_int][7],
91  X->Def.ParaInterAll_OffDiagonal[i_int],X, tmp_v0, tmp_v1);
92  X->Large.prdct += dam_pr;
93 #endif
94 }/*void GC_child_CisAitCiuAiv_spin_MPIdouble*/
struct EDMainCalStruct X
Definition: struct.h:431
double complex X_GC_child_CisAitCiuAiv_spin_MPIdouble(int org_isite1, int org_ispin1, int org_ispin2, int org_isite3, int org_ispin3, int org_ispin4, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
term in Spin model + GC. When both site1 and site2 are in the inter process region.

◆ GC_child_CisAitCiuAiv_spin_MPIsingle()

void GC_child_CisAitCiuAiv_spin_MPIsingle ( unsigned long int  i_int,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Exchange and Pairlifting term in Spin model + GC When only site2 is in the inter process region.

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]i_intInteraction ID
[in,out]X
[out]tmp_v0Result v0 = H v1
[in]tmp_v1v0 = H v1

Definition at line 500 of file mltplyMPISpinCore.c.

References X, and X_GC_child_CisAitCiuAiv_spin_MPIsingle().

Referenced by GC_child_general_int_spin_MPIsingle().

505  {
506 #ifdef MPI
507  double complex dam_pr;
509  X->Def.InterAll_OffDiagonal[i_int][0], X->Def.InterAll_OffDiagonal[i_int][1],
510  X->Def.InterAll_OffDiagonal[i_int][3], X->Def.InterAll_OffDiagonal[i_int][4],
511  X->Def.InterAll_OffDiagonal[i_int][5], X->Def.InterAll_OffDiagonal[i_int][7],
512  X->Def.ParaInterAll_OffDiagonal[i_int], X, tmp_v0, tmp_v1);
513  X->Large.prdct += dam_pr;
514 #endif
515 }/*void GC_child_CisAitCiuAiv_spin_MPIsingle*/
double complex X_GC_child_CisAitCiuAiv_spin_MPIsingle(int org_isite1, int org_ispin1, int org_ispin2, int org_isite3, int org_ispin3, int org_ispin4, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
Exchange and Pairlifting term in Spin model + GC When only site2 is in the inter process region...
struct EDMainCalStruct X
Definition: struct.h:431

◆ GC_child_CisAitCjuAju_spin_MPIdouble()

void GC_child_CisAitCjuAju_spin_MPIdouble ( unsigned long int  i_int,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Wrapper for calculating CisAitCjuAju term in Spin model + GC When both site1 and site2 are in the inter process region.

Author
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]i_intInteraction ID
[in,out]X
[out]tmp_v0Result v0 = H v1
[in]tmp_v1v0 = H v1

Definition at line 287 of file mltplyMPISpinCore.c.

References X, and X_GC_child_CisAitCjuAju_spin_MPIdouble().

Referenced by GC_child_general_int_spin_MPIdouble().

293 {
294 #ifdef MPI
295  double complex dam_pr;
297  X->Def.InterAll_OffDiagonal[i_int][0], X->Def.InterAll_OffDiagonal[i_int][1],
298  X->Def.InterAll_OffDiagonal[i_int][3], X->Def.InterAll_OffDiagonal[i_int][4],
299  X->Def.InterAll_OffDiagonal[i_int][5], X->Def.ParaInterAll_OffDiagonal[i_int], X, tmp_v0, tmp_v1);
300  X->Large.prdct += dam_pr;
301 #endif
302 }/*void GC_child_CisAitCiuAiv_spin_MPIdouble*/
struct EDMainCalStruct X
Definition: struct.h:431
double complex X_GC_child_CisAitCjuAju_spin_MPIdouble(int org_isite1, int org_ispin1, int org_ispin2, int org_isite3, int org_ispin3, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
CisAisCjuAjv term in Spin model + GC When both site1 and site2 are in the inter process region...

◆ GC_child_CisAitCjuAju_spin_MPIsingle()

void GC_child_CisAitCjuAju_spin_MPIsingle ( unsigned long int  i_int,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Wrapper for CisAisCjuAjv term in Spin model + GC When only site2 is in the inter process region.

Author
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]i_intInteraction ID
[in,out]X
[out]tmp_v0Result v0 = H v1
[in]tmp_v1v0 = H v1

Definition at line 711 of file mltplyMPISpinCore.c.

References X, and X_GC_child_CisAitCjuAju_spin_MPIsingle().

Referenced by GC_child_general_int_spin_MPIsingle().

716  {
717 #ifdef MPI
718  double complex dam_pr;
720  X->Def.InterAll_OffDiagonal[i_int][0], X->Def.InterAll_OffDiagonal[i_int][1],
721  X->Def.InterAll_OffDiagonal[i_int][3], X->Def.InterAll_OffDiagonal[i_int][4],
722  X->Def.InterAll_OffDiagonal[i_int][5], X->Def.ParaInterAll_OffDiagonal[i_int], X, tmp_v0, tmp_v1);
723  X->Large.prdct += dam_pr;
724 #endif
725 }/*void GC_child_CisAisCjuAjv_spin_MPIsingle*/
double complex X_GC_child_CisAitCjuAju_spin_MPIsingle(int org_isite1, int org_ispin1, int org_ispin2, int org_isite3, int org_ispin3, double complex tmp_J, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
CisAisCjuAjv term in Spin model + GC When only site2 is in the inter process region.
struct EDMainCalStruct X
Definition: struct.h:431

◆ X_child_CisAisCjuAju_GeneralSpin_MPIdouble()

double complex X_child_CisAisCjuAju_GeneralSpin_MPIdouble ( int  org_isite1,
int  org_ispin1,
int  org_isite3,
int  org_ispin3,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Compute \(c_{is}^\dagger c_{is}c_{ju}^\dagger c_{ju}\) term in the canonical general spin system when both sites are in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]tmp_JCoupling constant
[in,out]X
[in,out]tmp_v0Resulting wavefunction
[in]tmp_v1Input wavefunction

Definition at line 1743 of file mltplyMPISpinCore.c.

References BitCheckGeneral(), FALSE, myrank, and X.

Referenced by expec_cisajscktalt_SpinGeneral().

1752  {
1753 #ifdef MPI
1754  unsigned long int j, num1;
1755  double complex tmp_V, dmv, dam_pr;
1756 
1757  if (org_isite1 == org_isite3 && org_ispin1 == org_ispin3) {
1758  num1 = BitCheckGeneral((unsigned long int) myrank, org_isite1 + 1, org_ispin1, X->Def.SiteToBit, X->Def.Tpow);
1759  if (num1 != FALSE) {
1760  tmp_V = tmp_J;
1761  }
1762  else {
1763  return 0.0;
1764  }
1765  }
1766  else {
1767  num1 = BitCheckGeneral((unsigned long int) myrank, org_isite1 + 1, org_ispin1, X->Def.SiteToBit, X->Def.Tpow);
1768  if (num1 != FALSE) {
1769  num1 = BitCheckGeneral((unsigned long int) myrank, org_isite3 + 1, org_ispin3, X->Def.SiteToBit,
1770  X->Def.Tpow);
1771  if (num1 != FALSE) {
1772  tmp_V = tmp_J;
1773  }
1774  else {
1775  return 0.0;
1776  }
1777  }
1778  else {
1779  return 0.0;
1780  }
1781  }
1782  dam_pr = 0.0;
1783 #pragma omp parallel default(none) reduction(+:dam_pr) firstprivate(X, tmp_V) private(j, dmv) \
1784 shared (tmp_v0, tmp_v1)
1785  {
1786  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
1787 #pragma omp for
1788  for (j = 1; j <= X->Check.idim_max; j++) {
1789  dmv = tmp_v1[j] * tmp_V;
1790  tmp_v0[j] += dmv;
1791  dam_pr += conj(tmp_v1[j]) * dmv;
1792  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1793  }
1794  else {
1795 #pragma omp for
1796  for (j = 1; j <= X->Check.idim_max; j++) {
1797  dmv = tmp_v1[j] * tmp_V;
1798  dam_pr += conj(tmp_v1[j]) * dmv;
1799  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1800  }
1801  }/*End of parallel region*/
1802  return dam_pr;
1803 #else
1804  return 0.0;
1805 #endif
1806 }/*double complex X_child_CisAisCjuAju_GeneralSpin_MPIdouble*/
#define FALSE
Definition: global.h:25
int BitCheckGeneral(const long unsigned int org_bit, const unsigned int org_isite, const unsigned int target_ispin, const long int *SiteToBit, const long unsigned int *Tpow)
bit check function for general spin
Definition: bitcalc.c:393
struct EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162

◆ X_child_CisAisCjuAju_GeneralSpin_MPIsingle()

double complex X_child_CisAisCjuAju_GeneralSpin_MPIsingle ( int  org_isite1,
int  org_ispin1,
int  org_isite3,
int  org_ispin3,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Compute \(c_{is}^\dagger c_{is}c_{ju}^\dagger c_{ju}\) term in the canonical general spin system when one of these sites is in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]tmp_JCoupling constant
[in,out]X
[in,out]tmp_v0Resulting wavefunction
[in]tmp_v1Input wavefunction

Definition at line 1812 of file mltplyMPISpinCore.c.

References BitCheckGeneral(), FALSE, list_1, myrank, and X.

Referenced by expec_cisajscktalt_SpinGeneral().

1822 {
1823 #ifdef MPI
1824  unsigned long int j, num1;
1825  double complex tmp_V, dmv, dam_pr;
1826  //MPI_Status statusMPI;
1827 
1828  num1 = BitCheckGeneral((unsigned long int) myrank, org_isite3 + 1, org_ispin3, X->Def.SiteToBit, X->Def.Tpow);
1829  if (num1 != FALSE) {
1830  tmp_V = tmp_J;
1831  }
1832  else return 0.0;
1833 
1834  dam_pr = 0.0;
1835 #pragma omp parallel default(none) reduction(+:dam_pr) \
1836 firstprivate(X, tmp_V, org_isite1, org_ispin1) private(j, dmv, num1) shared (tmp_v0, tmp_v1, list_1)
1837  {
1838  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
1839 #pragma omp for
1840  for (j = 1; j <= X->Check.idim_max; j++) {
1841  num1 = BitCheckGeneral(list_1[j], org_isite1 + 1, org_ispin1, X->Def.SiteToBit, X->Def.Tpow);
1842 
1843  dmv = tmp_v1[j] * tmp_V * num1;
1844  tmp_v0[j] += dmv;
1845  dam_pr += conj(tmp_v1[j]) * dmv;
1846  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1847  }
1848  else {
1849 #pragma omp for
1850  for (j = 1; j <= X->Check.idim_max; j++) {
1851  num1 = BitCheckGeneral(list_1[j], org_isite1 + 1, org_ispin1, X->Def.SiteToBit, X->Def.Tpow);
1852 
1853  dmv = tmp_v1[j] * tmp_V * num1;
1854  dam_pr += conj(tmp_v1[j]) * dmv;
1855  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1856  }
1857  }/*End of parallel region*/
1858  return dam_pr;
1859 #else
1860  return 0.0;
1861 #endif
1862 }/*double complex X_child_CisAisCjuAju_GeneralSpin_MPIsingle*/
#define FALSE
Definition: global.h:25
long unsigned int * list_1
Definition: global.h:47
int BitCheckGeneral(const long unsigned int org_bit, const unsigned int org_isite, const unsigned int target_ispin, const long int *SiteToBit, const long unsigned int *Tpow)
bit check function for general spin
Definition: bitcalc.c:393
struct EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162

◆ X_child_CisAit_GeneralSpin_MPIdouble()

double complex X_child_CisAit_GeneralSpin_MPIdouble ( int  org_isite1,
int  org_ispin1,
int  org_ispin2,
double complex  tmp_trans,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1,
double complex *  tmp_v1buf,
unsigned long int  idim_max,
long unsigned int *  list_1_org,
long unsigned int *  list_1buf_org,
long unsigned int  _ihfbit 
)

Compute \(c_{is}^\dagger c_{it}\) term in the canonical general spin system when both site is in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_ispin2Spin 2
[in]tmp_transCoupling constant
[in,out]X
[in,out]tmp_v0Resulting wavefunction
[in]tmp_v1Input wavefunction
[in,out]tmp_v1bufbuffer for wavefunction
[in]idim_maxSimilar to CheckList::idim_max
[in]list_1_orgSimilar to list_1
[in]list_1buf_orgSimilar to list_1buf
[in]_ihfbitSimiler to LargeList::ihfbit

Definition at line 1290 of file mltplyMPISpinCore.c.

References ConvertToList1GeneralSpin(), exitMPI(), GetOffCompGeneralSpin(), list_1_org, list_1buf_org, myrank, TRUE, v1buf, and X.

Referenced by GetPairExcitedStateGeneralSpin().

1304 {
1305 #ifdef MPI
1306  unsigned long int off, j, tmp_off,idim_max_buf;
1307  int origin, ierr;
1308  double complex tmp_V, dmv;
1309  MPI_Status statusMPI;
1310 
1311  if (GetOffCompGeneralSpin((unsigned long int) myrank, org_isite1 + 1, org_ispin1, org_ispin2,
1312  &off, X->Def.SiteToBit, X->Def.Tpow) == TRUE) {
1313  tmp_V = tmp_trans;
1314  }
1315  else if (GetOffCompGeneralSpin((unsigned long int) myrank,
1316  org_isite1 + 1, org_ispin2, org_ispin1, &off,
1317  X->Def.SiteToBit, X->Def.Tpow) == TRUE) {
1318  tmp_V = conj(tmp_trans);
1319  if (X->Large.mode == M_CORR || X->Large.mode ==M_CALCSPEC) tmp_V = 0.0;
1320  }
1321  else return 0.0;
1322 
1323  origin = (int) off;
1324 
1325  ierr = MPI_Sendrecv(&idim_max, 1, MPI_UNSIGNED_LONG, origin, 0,
1326  &idim_max_buf, 1, MPI_UNSIGNED_LONG, origin, 0,
1327  MPI_COMM_WORLD, &statusMPI);
1328  if(ierr != 0) exitMPI(-1);
1329 
1330  ierr = MPI_Sendrecv(list_1_org, idim_max + 1, MPI_UNSIGNED_LONG, origin, 0,
1331  list_1buf_org, idim_max_buf + 1, MPI_UNSIGNED_LONG, origin, 0,
1332  MPI_COMM_WORLD, &statusMPI);
1333  if (ierr != 0) exitMPI(-1);
1334 
1335  ierr = MPI_Sendrecv(tmp_v1, idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
1336  v1buf, idim_max_buf + 1, MPI_DOUBLE_COMPLEX, origin, 0,
1337  MPI_COMM_WORLD, &statusMPI);
1338  if (ierr != 0) exitMPI(-1);
1339 
1340  if (X->Large.mode == M_MLTPLY || X->Large.mode ==M_CALCSPEC) {
1341 #pragma omp parallel for default(none)\
1342 firstprivate(X, tmp_V, idim_max_buf, list_1buf_org) private(j, dmv, tmp_off) \
1343 shared (tmp_v0, tmp_v1, v1buf)
1344  for (j = 1; j <= idim_max_buf; j++) {
1345  ConvertToList1GeneralSpin(list_1buf_org[j], X->Large.ihfbit, &tmp_off);
1346  dmv = v1buf[j] * tmp_V;
1347  tmp_v0[tmp_off] += dmv;
1348  }/*for (j = 1; j <= idim_max_buf; j++)*/
1349  }
1350  else {
1351  tmp_off = 0;
1352  return 0;
1353  }
1354  return 1;
1355 #else
1356  return 0.0;
1357 #endif
1358 }/*double complex X_child_CisAit_GeneralSpin_MPIdouble*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
int ConvertToList1GeneralSpin(const long unsigned int org_ibit, const long unsigned int ihlfbit, long unsigned int *_ilist1Comp)
function of converting component to list_1
Definition: bitcalc.c:285
#define TRUE
Definition: global.h:26
long unsigned int * list_1buf_org
Definition: global.h:54
long unsigned int * list_1_org
Definition: global.h:53
struct EDMainCalStruct X
Definition: struct.h:431
int GetOffCompGeneralSpin(const long unsigned int org_ibit, const int org_isite, const int org_ispin, const int off_ispin, long unsigned int *_ioffComp, const long int *SiteToBit, const long unsigned int *Tpow)
function of getting off-diagonal component for general spin
Definition: bitcalc.c:243
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37

◆ X_child_CisAit_spin_MPIdouble()

double complex X_child_CisAit_spin_MPIdouble ( int  org_isite1,
int  org_ispin2,
double complex  tmp_trans,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1,
double complex *  tmp_v1buf,
unsigned long int  idim_max,
long unsigned int *  Tpow,
long unsigned int *  list_1_org,
long unsigned int *  list_1buf_org,
long unsigned int *  list_2_1_target,
long unsigned int *  list_2_2_target,
long unsigned int  _irght,
long unsigned int  _ilft,
long unsigned int  _ihfbit 
)

Hopping term in Spin + Canonical for CalcSpectrum When both site1 and site2 are in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Author
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]org_isite1Site 1
[in]org_ispin2Spin 2
[in]tmp_transCoupling constant
[in,out]X
[out]tmp_v0Result v0 = H v1
[in]tmp_v1v0 = H v1
[in]tmp_v1bufbuffer for wavefunction
[in]idim_maxSimilar to CheckList::idim_max
[in]TpowSimilar to DefineList::Tpow
[in]list_1_orgSimilar to list_1
[in]list_1buf_orgSimilar to list_1buf
[in]list_2_1_targetSimilar to list_2_1
[in]list_2_2_targetSimilar to list_2_2
[in]_irghtSimiler to LargeList::irght
[in]_ilftSimiler to LargeList::ilft
[in]_ihfbitSimiler to LargeList::ihfbit

Definition at line 2039 of file mltplyMPISpinCore.c.

References exitMPI(), GetOffComp(), list_1_org, list_1buf_org, myrank, v1buf, and X.

Referenced by GetPairExcitedStateHalfSpin().

2056  {
2057 #ifdef MPI
2058  int mask1, state1, ierr, origin;
2059  unsigned long int idim_max_buf, j;
2060  unsigned long int tmp_off;
2061  MPI_Status statusMPI;
2062  double complex trans, dmv;
2063 
2064  mask1 = (int)X->Def.Tpow[org_isite1];
2065  origin = myrank ^ mask1;
2066  state1 = (origin & mask1)/mask1;
2067 
2068  if(state1 == org_ispin2){
2069  trans = tmp_trans;
2070  }
2071  else{
2072  trans =0.0;
2073  }
2074 
2075  // fprintf(stdout, "Debug: myrank=%d, origin=%d, trans=%lf\n", myrank, origin, trans);
2076 
2077  ierr = MPI_Sendrecv(&idim_max, 1, MPI_UNSIGNED_LONG, origin, 0,
2078  &idim_max_buf, 1, MPI_UNSIGNED_LONG, origin, 0,
2079  MPI_COMM_WORLD, &statusMPI);
2080  if (ierr != 0) exitMPI(-1);
2081 
2082  ierr = MPI_Sendrecv(list_1_org, idim_max + 1, MPI_UNSIGNED_LONG, origin, 0,
2083  list_1buf_org, idim_max_buf + 1, MPI_UNSIGNED_LONG, origin, 0,
2084  MPI_COMM_WORLD, &statusMPI);
2085  if (ierr != 0) exitMPI(-1);
2086 
2087  ierr = MPI_Sendrecv(tmp_v1, idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
2088  v1buf, idim_max_buf + 1, MPI_DOUBLE_COMPLEX, origin, 0,
2089  MPI_COMM_WORLD, &statusMPI);
2090  if (ierr != 0) exitMPI(-1);
2091 
2092  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
2093 #pragma omp parallel for default(none) private(j, dmv, tmp_off) \
2094 firstprivate(idim_max_buf, trans, X, list_1buf_org, list_2_1_target, list_2_2_target) \
2095 shared(v1buf, tmp_v0)
2096  for (j = 1; j <= idim_max_buf; j++) {
2097  GetOffComp(list_2_1_target, list_2_2_target, list_1buf_org[j], X->Large.irght, X->Large.ilft, X->Large.ihfbit, &tmp_off);
2098  dmv = trans * v1buf[j];
2099  tmp_v0[tmp_off] += dmv;
2100  }
2101  }
2102  else {
2103  tmp_off = 0;
2104  return 0;
2105  }
2106  return 1;
2107 #else
2108  return 0.0;
2109 #endif
2110 }/*double complex X_child_CisAit_spin_MPIdouble*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
long unsigned int * list_1buf_org
Definition: global.h:54
long unsigned int * list_1_org
Definition: global.h:53
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
struct EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37

◆ X_child_CisAitCjuAjv_GeneralSpin_MPIdouble()

double complex X_child_CisAitCjuAjv_GeneralSpin_MPIdouble ( int  org_isite1,
int  org_ispin1,
int  org_ispin2,
int  org_isite3,
int  org_ispin3,
int  org_ispin4,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Compute \(c_{is}^\dagger c_{it}c_{ju}^\dagger c_{jv}\) term in the canonical general spin system when both sites are in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_ispin2Spin 2
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]org_ispin4Spin 4
[in]tmp_JCoupling constant
[in,out]X
[in,out]tmp_v0Resulting wavefunction
[in]tmp_v1Input wavefunction

Definition at line 1647 of file mltplyMPISpinCore.c.

References ConvertToList1GeneralSpin(), exitMPI(), FALSE, GetOffCompGeneralSpin(), list_1, list_1buf, myrank, TRUE, v1buf, and X.

Referenced by child_general_int_GeneralSpin_MPIdouble(), and expec_cisajscktalt_SpinGeneral().

1658  {
1659 #ifdef MPI
1660  unsigned long int tmp_off, off, j, idim_max_buf;
1661  int origin, ierr;
1662  double complex tmp_V, dmv, dam_pr;
1663  MPI_Status statusMPI;
1664  int ihermite=TRUE;
1665 
1666  if (GetOffCompGeneralSpin((unsigned long int)myrank, org_isite1 + 1, org_ispin1, org_ispin2, &tmp_off, X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1667  {
1668  if (GetOffCompGeneralSpin(tmp_off, org_isite3 + 1, org_ispin3, org_ispin4, &off, X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1669  {
1670  tmp_V = tmp_J;
1671  }
1672  else{
1673  ihermite =FALSE;
1674  }
1675  }
1676  else{
1677  ihermite=FALSE;
1678  }
1679 
1680  if(ihermite==FALSE){
1681  if(GetOffCompGeneralSpin((unsigned long int)myrank, org_isite3 + 1, org_ispin4, org_ispin3, &tmp_off, X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1682  {
1683  if (GetOffCompGeneralSpin(tmp_off, org_isite1 + 1, org_ispin2, org_ispin1, &off, X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1684  {
1685  tmp_V = conj(tmp_J);
1686  if(X->Large.mode == M_CORR|| X->Large.mode == M_CALCSPEC){
1687  tmp_V=0.0;
1688  }
1689  }
1690  else return 0.0;
1691  }
1692  else return 0.0;
1693  }
1694 
1695 
1696  origin = (int)off;
1697 
1698  ierr = MPI_Sendrecv(&X->Check.idim_max, 1, MPI_UNSIGNED_LONG, origin, 0,
1699  &idim_max_buf, 1, MPI_UNSIGNED_LONG, origin, 0,
1700  MPI_COMM_WORLD, &statusMPI);
1701  if (ierr != 0) exitMPI(-1);
1702  ierr = MPI_Sendrecv(list_1, X->Check.idim_max + 1, MPI_UNSIGNED_LONG, origin, 0,
1703  list_1buf, idim_max_buf + 1, MPI_UNSIGNED_LONG, origin, 0,
1704  MPI_COMM_WORLD, &statusMPI);
1705  if (ierr != 0) exitMPI(-1);
1706  ierr = MPI_Sendrecv(tmp_v1, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
1707  v1buf, idim_max_buf + 1, MPI_DOUBLE_COMPLEX, origin, 0,
1708  MPI_COMM_WORLD, &statusMPI);
1709  if (ierr != 0) exitMPI(-1);
1710 
1711  dam_pr = 0.0;
1712 #pragma omp parallel default(none) reduction(+:dam_pr) firstprivate(X, tmp_V, idim_max_buf) \
1713 private(j, dmv, off) shared (tmp_v0, tmp_v1, list_1buf, v1buf)
1714  {
1715  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
1716 #pragma omp for
1717  for (j = 1; j <= idim_max_buf; j++) {
1718  ConvertToList1GeneralSpin(list_1buf[j], X->Check.sdim, &off);
1719  dmv = v1buf[j] * tmp_V;
1720  tmp_v0[off] += dmv;
1721  dam_pr += conj(tmp_v1[off]) * dmv;
1722  }/*for (j = 1; j <= idim_max_buf; j++)*/
1723  }
1724  else {
1725 #pragma omp for
1726  for (j = 1; j <= idim_max_buf; j++) {
1727  ConvertToList1GeneralSpin(list_1buf[j], X->Check.sdim, &off);
1728  dmv = v1buf[j] * tmp_V;
1729  dam_pr += conj(tmp_v1[off]) * dmv;
1730  }/*for (j = 1; j <= idim_max_buf; j++)*/
1731  }
1732  }/*End of parallel region*/
1733  return dam_pr;
1734 #else
1735  return 0.0;
1736 #endif
1737 }/*double complex X_child_CisAitCjuAjv_GeneralSpin_MPIdouble*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
long unsigned int * list_1buf
Definition: global.h:48
int ConvertToList1GeneralSpin(const long unsigned int org_ibit, const long unsigned int ihlfbit, long unsigned int *_ilist1Comp)
function of converting component to list_1
Definition: bitcalc.c:285
#define TRUE
Definition: global.h:26
#define FALSE
Definition: global.h:25
long unsigned int * list_1
Definition: global.h:47
struct EDMainCalStruct X
Definition: struct.h:431
int GetOffCompGeneralSpin(const long unsigned int org_ibit, const int org_isite, const int org_ispin, const int off_ispin, long unsigned int *_ioffComp, const long int *SiteToBit, const long unsigned int *Tpow)
function of getting off-diagonal component for general spin
Definition: bitcalc.c:243
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37

◆ X_child_CisAitCjuAjv_GeneralSpin_MPIsingle()

double complex X_child_CisAitCjuAjv_GeneralSpin_MPIsingle ( int  org_isite1,
int  org_ispin1,
int  org_ispin2,
int  org_isite3,
int  org_ispin3,
int  org_ispin4,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Compute \(c_{is}^\dagger c_{it}c_{ju}^\dagger c_{jv}\) term in the canonical general spin system when one of these sites is in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_ispin2Spin 2
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]org_ispin4Spin 4
[in]tmp_JCoupling constant
[in,out]X
[in,out]tmp_v0Resulting wavefunction
[in]tmp_v1Input wavefunction

Definition at line 1868 of file mltplyMPISpinCore.c.

References ConvertToList1GeneralSpin(), exitMPI(), GetOffCompGeneralSpin(), list_1, list_1buf, myrank, TRUE, v1buf, and X.

Referenced by child_general_int_GeneralSpin_MPIsingle(), and expec_cisajscktalt_SpinGeneral().

1879  {
1880 #ifdef MPI
1881  unsigned long int tmp_off, off, j, idim_max_buf;
1882  int origin, ierr, isite, IniSpin, FinSpin;
1883  double complex tmp_V, dmv, dam_pr;
1884  MPI_Status statusMPI;
1885 
1886  if (GetOffCompGeneralSpin((unsigned long int)myrank,
1887  org_isite3 + 1, org_ispin3, org_ispin4, &off,
1888  X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1889  {
1890  tmp_V = tmp_J;
1891  isite = org_isite1 + 1;
1892  IniSpin = org_ispin2;
1893  FinSpin = org_ispin1;
1894  }
1895  else if (GetOffCompGeneralSpin((unsigned long int)myrank,
1896  org_isite3 + 1, org_ispin4, org_ispin3, &off, X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1897  {
1898  tmp_V = conj(tmp_J);
1899  if (X->Large.mode == M_CORR || X->Large.mode == M_CALCSPEC) tmp_V = 0.0;
1900  isite = org_isite1 + 1;
1901  IniSpin = org_ispin1;
1902  FinSpin = org_ispin2;
1903  }
1904  else return 0.0;
1905 
1906  origin = (int)off;
1907 
1908  ierr = MPI_Sendrecv(&X->Check.idim_max, 1, MPI_UNSIGNED_LONG, origin, 0,
1909  &idim_max_buf, 1, MPI_UNSIGNED_LONG, origin, 0,
1910  MPI_COMM_WORLD, &statusMPI);
1911  if (ierr != 0) exitMPI(-1);
1912  ierr = MPI_Sendrecv(list_1, X->Check.idim_max + 1, MPI_UNSIGNED_LONG, origin, 0,
1913  list_1buf, idim_max_buf + 1, MPI_UNSIGNED_LONG, origin, 0,
1914  MPI_COMM_WORLD, &statusMPI);
1915  if (ierr != 0) exitMPI(-1);
1916  ierr = MPI_Sendrecv(tmp_v1, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
1917  v1buf, idim_max_buf + 1, MPI_DOUBLE_COMPLEX, origin, 0,
1918  MPI_COMM_WORLD, &statusMPI);
1919  if (ierr != 0) exitMPI(-1);
1920 
1921  dam_pr = 0.0;
1922 #pragma omp parallel default(none) reduction(+:dam_pr) \
1923 firstprivate(X, tmp_V, idim_max_buf, IniSpin, FinSpin, isite) \
1924 private(j, dmv, off, tmp_off) shared (tmp_v0, tmp_v1, list_1buf, v1buf)
1925  {
1926  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
1927 #pragma omp for
1928  for (j = 1; j <= idim_max_buf; j++) {
1929 
1930  if (GetOffCompGeneralSpin(list_1buf[j], isite, IniSpin, FinSpin, &tmp_off,
1931  X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1932  {
1933  ConvertToList1GeneralSpin(tmp_off, X->Check.sdim, &off);
1934  dmv = v1buf[j] * tmp_V;
1935  tmp_v0[off] += dmv;
1936  dam_pr += conj(tmp_v1[off]) * dmv;
1937  }
1938  }/*for (j = 1; j <= idim_max_buf; j++)*/
1939  }
1940  else {
1941 #pragma omp for
1942  for (j = 1; j <= idim_max_buf; j++) {
1943 
1944  if (GetOffCompGeneralSpin(list_1buf[j], isite, IniSpin, FinSpin, &tmp_off,
1945  X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1946  {
1947  ConvertToList1GeneralSpin(tmp_off, X->Check.sdim, &off);
1948  dmv = v1buf[j] * tmp_V;
1949  dam_pr += conj(tmp_v1[off]) * dmv;
1950  }
1951  }/*for (j = 1; j <= idim_max_buf; j++)*/
1952  }
1953  }/*End of parallel region*/
1954  return dam_pr;
1955 #else
1956  return 0.0;
1957 #endif
1958 }/*double complex X_child_CisAitCjuAjv_GeneralSpin_MPIsingle*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
long unsigned int * list_1buf
Definition: global.h:48
int ConvertToList1GeneralSpin(const long unsigned int org_ibit, const long unsigned int ihlfbit, long unsigned int *_ilist1Comp)
function of converting component to list_1
Definition: bitcalc.c:285
#define TRUE
Definition: global.h:26
long unsigned int * list_1
Definition: global.h:47
struct EDMainCalStruct X
Definition: struct.h:431
int GetOffCompGeneralSpin(const long unsigned int org_ibit, const int org_isite, const int org_ispin, const int off_ispin, long unsigned int *_ioffComp, const long int *SiteToBit, const long unsigned int *Tpow)
function of getting off-diagonal component for general spin
Definition: bitcalc.c:243
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37

◆ X_GC_child_AisCis_GeneralSpin_MPIdouble()

double complex X_GC_child_AisCis_GeneralSpin_MPIdouble ( int  org_isite1,
int  org_ispin1,
double complex  tmp_trans,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Compute \(c_{is} c_{is}^\dagger\) term in the grandcanonical general spin system when both site is in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]tmp_transCoupling constant
[in,out]X
[in,out]tmp_v0Resulting wavefunction
[in]tmp_v1Input wavefunction

Definition at line 1240 of file mltplyMPISpinCore.c.

References BitCheckGeneral(), myrank, and X.

Referenced by GetPairExcitedStateGeneralSpinGC().

1247  {
1248 #ifdef MPI
1249  unsigned long int j, num1;
1250  double complex tmp_V, dmv, dam_pr;
1251  //MPI_Status statusMPI;
1252 
1253  num1 = BitCheckGeneral((unsigned long int) myrank,
1254  org_isite1 + 1, org_ispin1, X->Def.SiteToBit, X->Def.Tpow);
1255  if (num1 == 0) {
1256  tmp_V = tmp_trans;
1257  }
1258  else return 0.0;
1259 
1260  dam_pr = 0.0;
1261 #pragma omp parallel default(none) reduction(+:dam_pr) firstprivate(X, tmp_V) private(j, dmv) \
1262 shared (tmp_v0, tmp_v1)
1263  {
1264  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
1265 #pragma omp for
1266  for (j = 1; j <= X->Check.idim_max; j++) {
1267  dmv = tmp_v1[j] * tmp_V;
1268  tmp_v0[j] += dmv;
1269  dam_pr += conj(tmp_v1[j]) * dmv;
1270  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1271  }
1272  else {
1273 #pragma omp for
1274  for (j = 1; j <= X->Check.idim_max; j++) {
1275  dmv = tmp_v1[j] * tmp_V;
1276  dam_pr += conj(tmp_v1[j]) * dmv;
1277  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1278  }
1279  }/*End of Parallel region*/
1280  return dam_pr;
1281 #else
1282  return 0.0;
1283 #endif
1284 }/*double complex X_GC_child_AisCis_GeneralSpin_MPIdouble*/
int BitCheckGeneral(const long unsigned int org_bit, const unsigned int org_isite, const unsigned int target_ispin, const long int *SiteToBit, const long unsigned int *Tpow)
bit check function for general spin
Definition: bitcalc.c:393
struct EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162

◆ X_GC_child_AisCis_spin_MPIdouble()

double complex X_GC_child_AisCis_spin_MPIdouble ( int  org_isite1,
int  org_ispin1,
double complex  tmp_trans,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Hopping term in Spin + GC When both site1 and site2 are in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Author
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]tmp_transCoupling constant
[in,out]X
[out]tmp_v0Result v0 = H v1
[in]tmp_v1v0 = H v1

Definition at line 2164 of file mltplyMPISpinCore.c.

References myrank, and X.

Referenced by GetPairExcitedStateHalfSpinGC().

2171  {
2172 #ifdef MPI
2173  long unsigned int j;
2174  int mask1;
2175  int ibit1;
2176  double complex dam_pr;
2177  mask1 = (int)X->Def.Tpow[org_isite1];
2178  ibit1 = (((unsigned long int)myrank& mask1) / mask1) ^ (1 - org_ispin1);
2179 
2180  dam_pr = 0.0;
2181 #pragma omp parallel reduction(+:dam_pr)default(none) shared(tmp_v1, tmp_v0, ibit1) \
2182  firstprivate(X, tmp_trans) private(j)
2183  {
2184  if (ibit1 == 0) {
2185  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) { // for multply
2186 #pragma omp for
2187  for (j = 1; j <= X->Check.idim_max; j++) {
2188  tmp_v0[j] += tmp_v1[j] * tmp_trans;
2189  dam_pr += tmp_trans * conj(tmp_v1[j]) * tmp_v1[j];
2190  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
2191  }
2192  else {
2193 #pragma omp for
2194  for (j = 1; j <= X->Check.idim_max; j++) {
2195  dam_pr += tmp_trans * conj(tmp_v1[j]) * tmp_v1[j];
2196  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
2197  }
2198  }/*if (ibit1 == 0)*/
2199  }/*End of parallel region*/
2200  return dam_pr;
2201 #else
2202  return 0.0;
2203 #endif
2204 }/*double complex X_GC_child_AisCis_spin_MPIdouble*/
struct EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162

◆ X_GC_child_CisAis_GeneralSpin_MPIdouble()

double complex X_GC_child_CisAis_GeneralSpin_MPIdouble ( int  org_isite1,
int  org_ispin1,
double complex  tmp_trans,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Compute \(c_{is}^\dagger c_{is}\) term in the grandcanonical general spin system when both site is in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]tmp_transCoupling constant
[in,out]X
[in,out]tmp_v0Resulting wavefunction
[in]tmp_v1Input wavefunction

Definition at line 1190 of file mltplyMPISpinCore.c.

References BitCheckGeneral(), myrank, and X.

Referenced by expec_cisajs_SpinGCGeneral(), GetPairExcitedStateGeneralSpinGC(), and X_GC_child_CisAitCjuAjv_GeneralSpin_MPIdouble().

1197  {
1198 #ifdef MPI
1199  unsigned long int j, num1;
1200  double complex tmp_V, dmv, dam_pr;
1201  //MPI_Status statusMPI;
1202 
1203  num1 = BitCheckGeneral((unsigned long int) myrank,
1204  org_isite1 + 1, org_ispin1, X->Def.SiteToBit, X->Def.Tpow);
1205  if (num1 != 0) {
1206  tmp_V = tmp_trans;
1207  }
1208  else return 0.0;
1209 
1210  dam_pr = 0.0;
1211 #pragma omp parallel default(none) reduction(+:dam_pr) firstprivate(X, tmp_V) private(j, dmv) \
1212 shared (tmp_v0, tmp_v1)
1213  {
1214  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
1215 #pragma omp for
1216  for (j = 1; j <= X->Check.idim_max; j++) {
1217  dmv = tmp_v1[j] * tmp_V;
1218  tmp_v0[j] += dmv;
1219  dam_pr += conj(tmp_v1[j]) * dmv;
1220  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1221  }
1222  else {
1223 #pragma omp for
1224  for (j = 1; j <= X->Check.idim_max; j++) {
1225  dmv = tmp_v1[j] * tmp_V;
1226  dam_pr += conj(tmp_v1[j]) * dmv;
1227  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1228  }
1229  }/*End of parallel region*/
1230  return dam_pr;
1231 #else
1232  return 0.0;
1233 #endif
1234 }/*double complex X_GC_child_CisAis_GeneralSpin_MPIdouble*/
int BitCheckGeneral(const long unsigned int org_bit, const unsigned int org_isite, const unsigned int target_ispin, const long int *SiteToBit, const long unsigned int *Tpow)
bit check function for general spin
Definition: bitcalc.c:393
struct EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162

◆ X_GC_child_CisAis_spin_MPIdouble()

double complex X_GC_child_CisAis_spin_MPIdouble ( int  org_isite1,
int  org_ispin1,
double complex  tmp_trans,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Hopping term in Spin + GC When both site1 and site2 are in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Author
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]tmp_transCoupling constant
[in,out]X
[out]tmp_v0Result v0 = H v1
[in]tmp_v1v0 = H v1

Definition at line 2117 of file mltplyMPISpinCore.c.

References myrank, and X.

Referenced by expec_cisajs_SpinGCHalf(), expec_cisajscktalt_SpinGCHalf(), GetPairExcitedStateHalfSpinGC(), and X_GC_child_CisAitCiuAiv_spin_MPIdouble().

2124  {
2125 #ifdef MPI
2126  long unsigned int j;
2127  int mask1;
2128  int ibit1;
2129  double complex dam_pr;
2130  mask1 = (int)X->Def.Tpow[org_isite1];
2131  ibit1 = (((unsigned long int)myrank& mask1)/mask1)^(1-org_ispin1);
2132 
2133  dam_pr = 0.0;
2134 #pragma omp parallel reduction(+:dam_pr)default(none) shared(tmp_v1, tmp_v0, ibit1) \
2135  firstprivate(X, tmp_trans) private(j)
2136  {
2137  if (ibit1 != 0) {
2138  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) { // for multply
2139 #pragma omp for
2140  for (j = 1; j <= X->Check.idim_max; j++) {
2141  tmp_v0[j] += tmp_v1[j] * tmp_trans;
2142  dam_pr += tmp_trans * conj(tmp_v1[j]) * tmp_v1[j];
2143  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
2144  }
2145  else {
2146 #pragma omp for
2147  for (j = 1; j <= X->Check.idim_max; j++) {
2148  dam_pr += tmp_trans * conj(tmp_v1[j]) * tmp_v1[j];
2149  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
2150  }
2151  }/*if (ibit1 != 0)*/
2152  }/*End of parallel region*/
2153  return dam_pr;
2154 #else
2155  return 0.0;
2156 #endif
2157 }/*double complex X_GC_child_CisAis_spin_MPIdouble*/
struct EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162

◆ X_GC_child_CisAisCjuAju_GeneralSpin_MPIdouble()

double complex X_GC_child_CisAisCjuAju_GeneralSpin_MPIdouble ( int  org_isite1,
int  org_ispin1,
int  org_isite3,
int  org_ispin3,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Compute \(c_{is}^\dagger c_{is} c_{ju}^\dagger c_{ju}\) term in the grandcanonical general spin system when both site is in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]tmp_JCoupling constant
[in,out]X
[in,out]tmp_v0Resulting wavefunction
[in]tmp_v1Input wavefunction

Definition at line 1070 of file mltplyMPISpinCore.c.

References BitCheckGeneral(), myrank, TRUE, and X.

Referenced by expec_cisajscktalt_SpinGCGeneral().

1079  {
1080 #ifdef MPI
1081  unsigned long int j, num1;
1082  double complex tmp_V, dmv, dam_pr;
1083  //MPI_Status statusMPI;
1084 
1085  num1 = BitCheckGeneral((unsigned long int) myrank, org_isite1 + 1, org_ispin1, X->Def.SiteToBit, X->Def.Tpow);
1086 
1087  if (num1 == TRUE) {
1088  num1 = BitCheckGeneral((unsigned long int) myrank, org_isite3 + 1, org_ispin3, X->Def.SiteToBit, X->Def.Tpow);
1089  if (num1 == TRUE) {
1090  tmp_V = tmp_J;
1091  }
1092  else return 0.0;
1093  }
1094  else return 0.0;
1095 
1096  dam_pr = 0.0;
1097 #pragma omp parallel default(none) reduction(+:dam_pr) firstprivate(X, tmp_V) private(j, dmv) \
1098 shared (tmp_v0, tmp_v1)
1099  {
1100  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
1101 #pragma omp for
1102  for (j = 1; j <= X->Check.idim_max; j++) {
1103  dmv = tmp_v1[j] * tmp_V;
1104  tmp_v0[j] += dmv;
1105  dam_pr += conj(tmp_v1[j]) * dmv;
1106  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1107  }
1108  else {
1109 #pragma omp for
1110  for (j = 1; j <= X->Check.idim_max; j++) {
1111  dmv = tmp_v1[j] * tmp_V;
1112  dam_pr += conj(tmp_v1[j]) * dmv;
1113  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1114  }
1115  }/*End of parallel region*/
1116  return dam_pr;
1117 #else
1118  return 0.0;
1119 #endif
1120 }/*double complex X_GC_child_CisAisCjuAju_GeneralSpin_MPIdouble*/
#define TRUE
Definition: global.h:26
int BitCheckGeneral(const long unsigned int org_bit, const unsigned int org_isite, const unsigned int target_ispin, const long int *SiteToBit, const long unsigned int *Tpow)
bit check function for general spin
Definition: bitcalc.c:393
struct EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162

◆ X_GC_child_CisAisCjuAju_GeneralSpin_MPIsingle()

double complex X_GC_child_CisAisCjuAju_GeneralSpin_MPIsingle ( int  org_isite1,
int  org_ispin1,
int  org_isite3,
int  org_ispin3,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Compute \(c_{is}^\dagger c_{is}c_{ju}^\dagger c_{ju}\) term in the grandcanonical general spin system when one of these site is in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]tmp_JCoupling constant
[in,out]X
[in,out]tmp_v0Resulting wavefunction
[in]tmp_v1Input wavefunction

Definition at line 1593 of file mltplyMPISpinCore.c.

References BitCheckGeneral(), FALSE, myrank, and X.

Referenced by expec_cisajscktalt_SpinGCGeneral().

1602  {
1603 #ifdef MPI
1604  unsigned long int j, num1;
1605  double complex tmp_V, dmv, dam_pr;
1606  //MPI_Status statusMPI;
1607 
1608  num1 = BitCheckGeneral((unsigned long int)myrank, org_isite3+1, org_ispin3, X->Def.SiteToBit, X->Def.Tpow);
1609  if (num1 != FALSE) {
1610  tmp_V = tmp_J;
1611  }
1612  else return 0.0;
1613 
1614  dam_pr = 0.0;
1615 #pragma omp parallel default(none) reduction(+:dam_pr) \
1616 firstprivate(X, tmp_V, org_isite1, org_ispin1) private(j, dmv, num1) shared (tmp_v0, tmp_v1)
1617  {
1618  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
1619 #pragma omp for
1620  for (j = 1; j <= X->Check.idim_max; j++) {
1621  num1 = BitCheckGeneral(j - 1, org_isite1 + 1, org_ispin1, X->Def.SiteToBit, X->Def.Tpow);
1622 
1623  dmv = tmp_v1[j] * tmp_V * num1;
1624  tmp_v0[j] += dmv;
1625  dam_pr += conj(tmp_v1[j]) * dmv;
1626  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1627  }
1628  else {
1629 #pragma omp for
1630  for (j = 1; j <= X->Check.idim_max; j++) {
1631  num1 = BitCheckGeneral(j - 1, org_isite1 + 1, org_ispin1, X->Def.SiteToBit, X->Def.Tpow);
1632  dmv = tmp_v1[j] * tmp_V * num1;
1633  dam_pr += conj(tmp_v1[j]) * dmv;
1634  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1635  }
1636  }/*End of parallel region*/
1637  return dam_pr;
1638 #else
1639  return 0.0;
1640 #endif
1641 }/*double complex X_GC_child_CisAisCjuAju_GeneralSpin_MPIsingle*/
#define FALSE
Definition: global.h:25
int BitCheckGeneral(const long unsigned int org_bit, const unsigned int org_isite, const unsigned int target_ispin, const long int *SiteToBit, const long unsigned int *Tpow)
bit check function for general spin
Definition: bitcalc.c:393
struct EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162

◆ X_GC_child_CisAisCjuAju_spin_MPIdouble()

double complex X_GC_child_CisAisCjuAju_spin_MPIdouble ( int  org_isite1,
int  org_ispin1,
int  org_isite3,
int  org_ispin3,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

CisAisCjuAjv term in Spin model + GC When both site1 and site2 are in the inter process region.

Returns
\(\langle v_1 | H_{\rm this} | v_1 \rangle\)
Author
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]tmp_JCopupling constatnt
[in,out]X
[in,out]tmp_v0\({\bf v}_0=H {\bf v}_1\)
[in]tmp_v1Vector to be producted

Definition at line 397 of file mltplyMPISpinCore.c.

References myrank, X, and X_SpinGC_CisAis().

Referenced by expec_cisajscktalt_SpinGCHalf().

406  {
407 #ifdef MPI
408  long unsigned int mask1, mask2, num1,num2;
409  unsigned long int j;
410 // MPI_Status statusMPI;
411  double complex dmv, dam_pr;
412  mask1 = (int)X->Def.Tpow[org_isite1];
413  mask2 = (int)X->Def.Tpow[org_isite3];
414  num1 = X_SpinGC_CisAis((unsigned long int)myrank + 1, X, mask1, org_ispin1);
415  num2 = X_SpinGC_CisAis((unsigned long int)myrank + 1, X, mask2, org_ispin3);
416 
417  dam_pr = 0.0;
418 #pragma omp parallel default(none) reduction(+:dam_pr) private(j, dmv) \
419  firstprivate(tmp_J, X, num1, num2) shared(tmp_v1, tmp_v0)
420  {
421  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
422 #pragma omp for
423  for (j = 1; j <= X->Check.idim_max; j++) {
424  dmv = num1*num2*tmp_v1[j] * tmp_J;
425  tmp_v0[j] += dmv;
426  dam_pr += conj(tmp_v1[j]) * dmv;
427  }/*for (j = 1; j <= X->Check.idim_max; j++) */
428  }
429  else {
430 #pragma omp for
431  for (j = 1; j <= X->Check.idim_max; j++) {
432  dmv = num1 * num2 * tmp_v1[j] * tmp_J;
433  dam_pr += conj(tmp_v1[j]) * dmv;
434  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
435  }
436  }/*End of parallel region*/
437  return(dam_pr);
438 #else
439  return 0.0;
440 #endif
441 }/*double complex X_GC_child_CisAisCjuAju_spin_MPIdouble*/
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 EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162

◆ X_GC_child_CisAisCjuAju_spin_MPIsingle()

double complex X_GC_child_CisAisCjuAju_spin_MPIsingle ( int  org_isite1,
int  org_ispin1,
int  org_isite3,
int  org_ispin3,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

CisAisCjuAjv term in Spin model + GC When both site1 and site2 are in the inter process region.

Returns
\(\langle v_1 | H_{\rm this} | v_1 \rangle\)
Author
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]tmp_JCopupling constatnt
[in,out]X
[in,out]tmp_v0\({\bf v}_0=H {\bf v}_1\)
[in]tmp_v1Vector to be producted

Definition at line 448 of file mltplyMPISpinCore.c.

References myrank, X, and X_SpinGC_CisAis().

Referenced by expec_cisajscktalt_SpinGCHalf().

457  {
458 #ifdef MPI
459  long unsigned int mask1, mask2, num1, num2;
460  unsigned long int j;
461 // MPI_Status statusMPI;
462  double complex Jint, dmv, dam_pr;
463  Jint = tmp_J;
464  mask1 = (int)X->Def.Tpow[org_isite1];
465  mask2 = (int)X->Def.Tpow[org_isite3];
466  num2 = X_SpinGC_CisAis((unsigned long int) myrank + 1, X, mask2, org_ispin3);
467 
468  dam_pr = 0.0;
469 #pragma omp parallel default(none) reduction(+:dam_pr) private(j, dmv, num1) \
470  firstprivate(Jint, X, num2, mask1, org_ispin1) shared(tmp_v1, tmp_v0)
471  {
472  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
473 #pragma omp for
474  for (j = 1; j <= X->Check.idim_max; j++) {
475  num1 = X_SpinGC_CisAis(j, X, mask1, org_ispin1);
476  dmv = Jint * num1 * num2 * tmp_v1[j];
477  tmp_v0[j] += dmv;
478  dam_pr += conj(tmp_v1[j]) * dmv;
479  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
480  }
481  else {
482 #pragma omp for
483  for (j = 1; j <= X->Check.idim_max; j++) {
484  num1 = X_SpinGC_CisAis(j, X, mask1, org_ispin1);
485  dmv = Jint * num1 * num2 * tmp_v1[j];
486  dam_pr += conj(tmp_v1[j]) * dmv;
487  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
488  }
489  }/*End of parallel region*/
490  return (dam_pr);
491 #else
492  return 0.0;
493 #endif
494 }/*double complex X_GC_child_CisAisCjuAju_spin_MPIdouble*/
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 EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162

◆ X_GC_child_CisAisCjuAjv_GeneralSpin_MPIdouble()

double complex X_GC_child_CisAisCjuAjv_GeneralSpin_MPIdouble ( int  org_isite1,
int  org_ispin1,
int  org_isite3,
int  org_ispin3,
int  org_ispin4,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

\(c_{is}^\dagger c_{is} c_{ju}^\dagger c_{jv}\) term in Spin model. When both site1 and site3 are in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]org_ispin4Spin 4
[in]tmp_JCopupling constatnt
[in,out]X
[in,out]tmp_v0\({\bf v}_0=H {\bf v}_1\)
[in]tmp_v1Vector to be producted

Definition at line 824 of file mltplyMPISpinCore.c.

References BitCheckGeneral(), exitMPI(), FALSE, GetOffCompGeneralSpin(), myrank, TRUE, v1buf, and X.

Referenced by expec_cisajscktalt_SpinGCGeneral(), and GC_child_general_int_GeneralSpin_MPIdouble().

834  {
835 #ifdef MPI
836  unsigned long int off, j;
837  int origin, ierr;
838  double complex tmp_V, dmv, dam_pr;
839  MPI_Status statusMPI;
840  int ihermite = TRUE;
841  if (org_isite1 == org_isite3 && org_ispin1 == org_ispin4) {//cisaisciuais=0 && cisaiucisais=0
842  return 0.0;
843  }
844 
845  if (GetOffCompGeneralSpin((unsigned long int) myrank, org_isite3 + 1, org_ispin4, org_ispin3,
846  &off, X->Def.SiteToBit, X->Def.Tpow) == TRUE) {
847  if (BitCheckGeneral(off, org_isite1 + 1, org_ispin1, X->Def.SiteToBit, X->Def.Tpow) == TRUE) {
848  tmp_V = tmp_J;
849  }
850  else {
851  ihermite = FALSE;
852  }
853  }
854  else {
855  ihermite = FALSE;
856  }
857 
858  if (ihermite == FALSE) {
859  if (BitCheckGeneral((unsigned long int) myrank, org_isite1 + 1, org_ispin1, X->Def.SiteToBit, X->Def.Tpow) ==
860  TRUE &&
861  GetOffCompGeneralSpin((unsigned long int) myrank, org_isite3 + 1, org_ispin3, org_ispin4, &off,
862  X->Def.SiteToBit, X->Def.Tpow) == TRUE) {
863  tmp_V = conj(tmp_J);
864  if (X->Large.mode == M_CORR || X->Large.mode == M_CALCSPEC) tmp_V = 0.0;
865  }
866  else {
867  return 0.0;
868  }
869  }
870  origin = (int)off;
871  ierr = MPI_Sendrecv(tmp_v1, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
872  v1buf, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
873  MPI_COMM_WORLD, &statusMPI);
874  if (ierr != 0) exitMPI(-1);
875 
876  dam_pr = 0.0;
877 #pragma omp parallel default(none) reduction(+:dam_pr) firstprivate(X, tmp_V) \
878 private(j, dmv) shared (tmp_v0, tmp_v1, v1buf)
879  {
880  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
881 #pragma omp for
882  for (j = 1; j <= X->Check.idim_max; j++) {
883  dmv = v1buf[j] * tmp_V;
884  tmp_v0[j] += dmv;
885  dam_pr += conj(tmp_v1[j]) * dmv;
886  }
887  }
888  else {
889 #pragma omp for
890  for (j = 1; j <= X->Check.idim_max; j++) {
891  dmv = v1buf[j] * tmp_V;
892  dam_pr += conj(tmp_v1[j]) * dmv;
893  }
894  }
895  }/*End of parallel region*/
896  return dam_pr;
897 #else
898  return 0.0;
899 #endif
900 }/*double complex X_GC_child_CisAisCjuAjv_GeneralSpin_MPIdouble*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
#define TRUE
Definition: global.h:26
#define FALSE
Definition: global.h:25
int BitCheckGeneral(const long unsigned int org_bit, const unsigned int org_isite, const unsigned int target_ispin, const long int *SiteToBit, const long unsigned int *Tpow)
bit check function for general spin
Definition: bitcalc.c:393
struct EDMainCalStruct X
Definition: struct.h:431
int GetOffCompGeneralSpin(const long unsigned int org_ibit, const int org_isite, const int org_ispin, const int off_ispin, long unsigned int *_ioffComp, const long int *SiteToBit, const long unsigned int *Tpow)
function of getting off-diagonal component for general spin
Definition: bitcalc.c:243
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37

◆ X_GC_child_CisAisCjuAjv_GeneralSpin_MPIsingle()

double complex X_GC_child_CisAisCjuAjv_GeneralSpin_MPIsingle ( int  org_isite1,
int  org_ispin1,
int  org_isite3,
int  org_ispin3,
int  org_ispin4,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Compute \(c_{is}^\dagger c_{is}c_{ju}^\dagger c_{jv}\) term in the grandcanonical general spin system when one of these site is in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]org_ispin4Spin 4
[in]tmp_JCoupling constant
[in,out]X
[in,out]tmp_v0Resulting wavefunction
[in]tmp_v1Input wavefunction

Definition at line 1364 of file mltplyMPISpinCore.c.

References BitCheckGeneral(), exitMPI(), GetOffCompGeneralSpin(), myrank, TRUE, v1buf, and X.

Referenced by expec_cisajscktalt_SpinGCGeneral(), and GC_child_general_int_GeneralSpin_MPIsingle().

1374  {
1375 #ifdef MPI
1376  unsigned long int off, j, num1;
1377  int origin, ierr, isite, IniSpin;
1378  double complex tmp_V, dmv, dam_pr;
1379  MPI_Status statusMPI;
1380 
1381  if (GetOffCompGeneralSpin((unsigned long int)myrank,
1382  org_isite3 + 1, org_ispin3, org_ispin4, &off,
1383  X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1384  {
1385  tmp_V = tmp_J;
1386  isite = org_isite1 + 1;
1387  IniSpin = org_ispin1;
1388  }
1389  else if (GetOffCompGeneralSpin((unsigned long int)myrank,
1390  org_isite3 + 1, org_ispin4, org_ispin3, &off,
1391  X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1392  {
1393  tmp_V = conj(tmp_J);
1394  if (X->Large.mode == M_CORR || X->Large.mode == M_CALCSPEC) tmp_V = 0.0;
1395  isite = org_isite1 + 1;
1396  IniSpin = org_ispin1;
1397  }
1398  else return 0.0;
1399 
1400  origin = (int)off;
1401 
1402  ierr = MPI_Sendrecv(tmp_v1, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
1403  v1buf, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
1404  MPI_COMM_WORLD, &statusMPI);
1405  if (ierr != 0) exitMPI(-1);
1406 
1407  dam_pr = 0.0;
1408 #pragma omp parallel default(none) reduction(+:dam_pr) firstprivate(X, tmp_V, isite, IniSpin) \
1409 private(j, dmv, num1) shared (tmp_v0, tmp_v1, v1buf)
1410  {
1411  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
1412 #pragma omp for
1413  for (j = 1; j <= X->Check.idim_max; j++) {
1414  num1 = BitCheckGeneral(j - 1, isite, IniSpin, X->Def.SiteToBit, X->Def.Tpow);
1415  if (num1 != 0) {
1416  dmv = v1buf[j] * tmp_V;
1417  tmp_v0[j] += dmv;
1418  dam_pr += conj(tmp_v1[j]) * dmv;
1419  }/*if (num1 != 0)*/
1420  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1421  }
1422  else {
1423 #pragma omp for
1424  for (j = 1; j <= X->Check.idim_max; j++) {
1425  num1 = BitCheckGeneral(j - 1, isite, IniSpin, X->Def.SiteToBit, X->Def.Tpow);
1426  if (num1 != 0) {
1427  dmv = v1buf[j] * tmp_V;
1428  dam_pr += conj(tmp_v1[j]) * dmv;
1429  }/*if (num1 != 0)*/
1430  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1431  }
1432  }/*End of parallel region*/
1433  return dam_pr;
1434 #else
1435  return 0.0;
1436 #endif
1437 }/*double complex X_GC_child_CisAisCjuAjv_GeneralSpin_MPIsingle*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
#define TRUE
Definition: global.h:26
int BitCheckGeneral(const long unsigned int org_bit, const unsigned int org_isite, const unsigned int target_ispin, const long int *SiteToBit, const long unsigned int *Tpow)
bit check function for general spin
Definition: bitcalc.c:393
struct EDMainCalStruct X
Definition: struct.h:431
int GetOffCompGeneralSpin(const long unsigned int org_ibit, const int org_isite, const int org_ispin, const int off_ispin, long unsigned int *_ioffComp, const long int *SiteToBit, const long unsigned int *Tpow)
function of getting off-diagonal component for general spin
Definition: bitcalc.c:243
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37

◆ X_GC_child_CisAisCjuAjv_spin_MPIdouble()

double complex X_GC_child_CisAisCjuAjv_spin_MPIdouble ( int  org_isite1,
int  org_ispin1,
int  org_isite3,
int  org_ispin3,
int  org_ispin4,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

CisAisCjuAjv term in Spin model + GC When both site1 and site2 are in the inter process region.

Returns
\(\langle v_1 | H_{\rm this} | v_1 \rangle\)
Author
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]org_ispin4Spin 4
[in]tmp_JCopupling constatnt
[in,out]X
[in,out]tmp_v0\({\bf v}_0=H {\bf v}_1\)
[in]tmp_v1Vector to be producted

Definition at line 212 of file mltplyMPISpinCore.c.

References exitMPI(), myrank, TRUE, v1buf, X, and X_SpinGC_CisAis().

Referenced by expec_cisajscktalt_SpinGCHalf(), and GC_child_CisAisCjuAjv_spin_MPIdouble().

222  {
223 #ifdef MPI
224  int mask1, mask2, state2, ierr;
225  long int origin, num1;
226  unsigned long int idim_max_buf, j;
227  MPI_Status statusMPI;
228  double complex Jint, dmv, dam_pr;
229 
230  if (org_isite1 == org_isite3 && org_ispin1 == org_ispin4) {//CisAisCitAis
231  return 0.0;
232  }
233 
234  mask1 = (int)X->Def.Tpow[org_isite1];
235  mask2 = (int)X->Def.Tpow[org_isite3];
236  origin = myrank ^ mask2;
237  state2 = (origin & mask2) / mask2;
238  num1 = X_SpinGC_CisAis((unsigned long int) myrank + 1, X, mask1, org_ispin1);
239  if (num1 != 0 && state2 == org_ispin4) {
240  Jint = tmp_J;
241  }
242  else if (X_SpinGC_CisAis(origin + 1, X, mask1, org_ispin1) == TRUE && state2 == org_ispin3) {
243  Jint = conj(tmp_J);
244  if (X->Large.mode == M_CORR || X->Large.mode == M_CALCSPEC) Jint = 0;
245  }
246  else {
247  return 0.0;
248  }
249 
250  ierr = MPI_Sendrecv(&X->Check.idim_max, 1, MPI_UNSIGNED_LONG, origin, 0,
251  &idim_max_buf, 1, MPI_UNSIGNED_LONG, origin, 0,
252  MPI_COMM_WORLD, &statusMPI);
253  if (ierr != 0) exitMPI(-1);
254  ierr = MPI_Sendrecv(tmp_v1, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
255  v1buf, idim_max_buf + 1, MPI_DOUBLE_COMPLEX, origin, 0,
256  MPI_COMM_WORLD, &statusMPI);
257  if (ierr != 0) exitMPI(-1);
258 
259  dam_pr = 0.0;
260  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
261 #pragma omp parallel for default(none) reduction(+:dam_pr) private(j, dmv) \
262  firstprivate(idim_max_buf, Jint, X) shared(v1buf, tmp_v1, tmp_v0)
263  for (j = 1; j <= idim_max_buf; j++) {
264  dmv = Jint * v1buf[j];
265  tmp_v0[j] += dmv;
266  dam_pr += conj(tmp_v1[j]) * dmv;
267  }
268  }
269  else {
270 #pragma omp parallel for default(none) reduction(+:dam_pr) private(j, dmv) \
271  firstprivate(idim_max_buf, Jint, X) shared(v1buf, tmp_v1, tmp_v0)
272  for (j = 1; j <= idim_max_buf; j++) {
273  dmv = Jint * v1buf[j];
274  dam_pr += conj(tmp_v1[j]) * dmv;
275  }
276  }
277  return (dam_pr);
278 #else
279  return 0.0;
280 #endif
281 }/*double complex X_GC_child_CisAisCjuAjv_spin_MPIdouble*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
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.
#define TRUE
Definition: global.h:26
struct EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37

◆ X_GC_child_CisAisCjuAjv_spin_MPIsingle()

double complex X_GC_child_CisAisCjuAjv_spin_MPIsingle ( int  org_isite1,
int  org_ispin1,
int  org_isite3,
int  org_ispin3,
int  org_ispin4,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

CisAisCjuAjv term in Spin model + GC When only site2 is in the inter process region.

Returns
\(\langle v_1 | H_{\rm this} | v_1 \rangle\)
Author
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 2
[in]org_isite3Site 1
[in]org_ispin3Spin 2
[in]org_ispin4Spin 2
[in]tmp_JCopupling constatnt
[in,out]X
[in,out]tmp_v0\({\bf v}_0=H {\bf v}_1\)
[in]tmp_v1Vector to be producted

Definition at line 629 of file mltplyMPISpinCore.c.

References exitMPI(), myrank, v1buf, and X.

Referenced by expec_cisajscktalt_SpinGCHalf(), and GC_child_CisAisCjuAjv_spin_MPIsingle().

639  {
640 #ifdef MPI
641  int mask2, state2, ierr, origin;
642  unsigned long int mask1, idim_max_buf, j, state1, state1check;
643  MPI_Status statusMPI;
644  double complex Jint, dmv, dam_pr;
645  /*
646  Prepare index in the inter PE
647  */
648  mask2 = (int)X->Def.Tpow[org_isite3];
649  origin = myrank ^ mask2;
650  state2 = (origin & mask2) / mask2;
651  if (state2 == org_ispin4) {
652  state1check = (unsigned long int) org_ispin1;
653  Jint = tmp_J;
654  }
655  else if (state2 == org_ispin3) {
656  state1check = (unsigned long int) org_ispin1;
657  Jint = conj(tmp_J);
658  if (X->Large.mode == M_CORR || X->Large.mode == M_CALCSPEC) {
659  Jint = 0;
660  }
661  }
662  else return 0.0;
663 
664  ierr = MPI_Sendrecv(&X->Check.idim_max, 1, MPI_UNSIGNED_LONG, origin, 0,
665  &idim_max_buf, 1, MPI_UNSIGNED_LONG, origin, 0,
666  MPI_COMM_WORLD, &statusMPI);
667  if (ierr != 0) exitMPI(-1);
668  ierr = MPI_Sendrecv(tmp_v1, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
669  v1buf, idim_max_buf + 1, MPI_DOUBLE_COMPLEX, origin, 0,
670  MPI_COMM_WORLD, &statusMPI);
671  if (ierr != 0) exitMPI(-1);
672  /*
673  Index in the intra PE
674  */
675  mask1 = X->Def.Tpow[org_isite1];
676 
677  dam_pr = 0.0;
678 #pragma omp parallel default(none) reduction(+:dam_pr) private(j, dmv, state1) \
679  firstprivate(idim_max_buf, Jint, X, state1check, mask1) shared(v1buf, tmp_v1, tmp_v0)
680  {
681  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
682  for (j = 0; j < idim_max_buf; j++) {
683  state1 = (j & mask1) / mask1;
684  if (state1 == state1check) {
685  dmv = Jint * v1buf[j + 1];
686  tmp_v0[j + 1] += dmv;
687  dam_pr += conj(tmp_v1[j + 1]) * dmv;
688  }/*if (state1 == state1check)*/
689  }/*for (j = 0; j < idim_max_buf; j++)*/
690  }
691  else {
692  for (j = 0; j < idim_max_buf; j++) {
693  state1 = (j & mask1) / mask1;
694  if (state1 == state1check) {
695  dmv = Jint * v1buf[j + 1];
696  dam_pr += conj(tmp_v1[j + 1]) * dmv;
697  }/*if (state1 == state1check)*/
698  }/*for (j = 0; j < idim_max_buf; j++)*/
699  }
700  }/*End of parallel region*/
701  return (dam_pr);
702 #else
703  return 0.0;
704 #endif
705 }/*void GC_child_CisAitCiuAiv_spin_MPIsingle*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
struct EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37

◆ X_GC_child_CisAit_GeneralSpin_MPIdouble()

double complex X_GC_child_CisAit_GeneralSpin_MPIdouble ( int  org_isite1,
int  org_ispin1,
int  org_ispin2,
double complex  tmp_trans,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Compute \(c_{is}^\dagger c_{it}\) term in the grandcanonical general spin system when both site is in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_ispin2Spin 2
[in]tmp_transCoupling constant
[in,out]X
[in,out]tmp_v0Resulting wavefunction
[in]tmp_v1Input wavefunction

Definition at line 1126 of file mltplyMPISpinCore.c.

References exitMPI(), GetOffCompGeneralSpin(), myrank, TRUE, v1buf, and X.

Referenced by expec_cisajs_SpinGCGeneral(), GetPairExcitedStateGeneralSpinGC(), and mltplyGeneralSpinGC().

1134  {
1135 #ifdef MPI
1136  unsigned long int off, j;
1137  int origin, ierr;
1138  double complex tmp_V, dmv, dam_pr;
1139  MPI_Status statusMPI;
1140 
1141  if (GetOffCompGeneralSpin((unsigned long int) myrank, org_isite1 + 1, org_ispin1, org_ispin2,
1142  &off, X->Def.SiteToBit, X->Def.Tpow) == TRUE) {
1143  tmp_V = tmp_trans;
1144  }
1145  else if (GetOffCompGeneralSpin((unsigned long int) myrank,
1146  org_isite1 + 1, org_ispin2, org_ispin1, &off,
1147  X->Def.SiteToBit, X->Def.Tpow) == TRUE) {
1148  tmp_V = conj(tmp_trans);
1149  if (X->Large.mode == M_CORR || X->Large.mode == M_CALCSPEC) tmp_V = 0.0;
1150  }
1151  else return 0.0;
1152 
1153  origin = (int)off;
1154 
1155  ierr = MPI_Sendrecv(tmp_v1, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
1156  v1buf, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
1157  MPI_COMM_WORLD, &statusMPI);
1158  if (ierr != 0) exitMPI(-1);
1159 
1160  dam_pr = 0.0;
1161 #pragma omp parallel default(none) reduction(+:dam_pr) firstprivate(X, tmp_V) private(j, dmv) \
1162 shared (tmp_v0, tmp_v1, v1buf)
1163  {
1164  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
1165 #pragma omp for
1166  for (j = 1; j <= X->Check.idim_max; j++) {
1167  dmv = v1buf[j] * tmp_V;
1168  tmp_v0[j] += dmv;
1169  dam_pr += conj(tmp_v1[j]) * dmv;
1170  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1171  }
1172  else {
1173 #pragma omp for
1174  for (j = 1; j <= X->Check.idim_max; j++) {
1175  dmv = v1buf[j] * tmp_V;
1176  dam_pr += conj(tmp_v1[j]) * dmv;
1177  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1178  }
1179  }/*End of parallel region*/
1180  return dam_pr;
1181 #else
1182  return 0.0;
1183 #endif
1184 }/*double complex X_GC_child_CisAit_GeneralSpin_MPIdouble*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
#define TRUE
Definition: global.h:26
struct EDMainCalStruct X
Definition: struct.h:431
int GetOffCompGeneralSpin(const long unsigned int org_ibit, const int org_isite, const int org_ispin, const int off_ispin, long unsigned int *_ioffComp, const long int *SiteToBit, const long unsigned int *Tpow)
function of getting off-diagonal component for general spin
Definition: bitcalc.c:243
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37

◆ X_GC_child_CisAit_spin_MPIdouble()

double complex X_GC_child_CisAit_spin_MPIdouble ( int  org_isite1,
int  org_ispin1,
int  org_ispin2,
double complex  tmp_trans,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Hopping term in Spin + GC When both site1 and site2 are in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Author
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_ispin2Spin 2
[in]tmp_transCoupling constant
[in,out]X
[out]tmp_v0Result v0 = H v1
[in]tmp_v1v0 = H v1

Definition at line 1965 of file mltplyMPISpinCore.c.

References exitMPI(), myrank, v1buf, and X.

Referenced by expec_cisajs_SpinGCHalf(), GetPairExcitedStateHalfSpinGC(), and mltplyHalfSpinGC().

1973 {
1974 #ifdef MPI
1975  int mask1, state1, ierr, origin;
1976  unsigned long int idim_max_buf, j;
1977  MPI_Status statusMPI;
1978  double complex trans, dmv, dam_pr;
1979 
1980  mask1 = (int)X->Def.Tpow[org_isite1];
1981  origin = myrank ^ mask1;
1982  state1 = (origin & mask1)/mask1;
1983 
1984  //fprintf(stdout, "Debug: myrank=%d, origin=%d, state1=%d\n", myrank, origin, state1);
1985 
1986  if(state1 == org_ispin2){
1987  trans = tmp_trans;
1988  }
1989  else if(state1 == org_ispin1) {
1990  trans = conj(tmp_trans);
1991  if(X->Large.mode == M_CORR|| X->Large.mode ==M_CALCSPEC){
1992  trans = 0.0;
1993  }
1994  }
1995  else{
1996  return 0.0;
1997  }
1998 
1999  ierr = MPI_Sendrecv(&X->Check.idim_max, 1, MPI_UNSIGNED_LONG, origin, 0,
2000  &idim_max_buf, 1, MPI_UNSIGNED_LONG, origin, 0,
2001  MPI_COMM_WORLD, &statusMPI);
2002  if (ierr != 0) exitMPI(-1);
2003  ierr = MPI_Sendrecv(tmp_v1, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
2004  v1buf, idim_max_buf + 1, MPI_DOUBLE_COMPLEX, origin, 0,
2005  MPI_COMM_WORLD, &statusMPI);
2006  if (ierr != 0) exitMPI(-1);
2007 
2008  dam_pr = 0.0;
2009 #pragma omp parallel default(none) reduction(+:dam_pr) private(j, dmv) \
2010 firstprivate(idim_max_buf, trans, X) shared(v1buf, tmp_v1, tmp_v0)
2011  {
2012  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
2013 #pragma omp for
2014  for (j = 1; j <= X->Check.idim_max; j++) {
2015  dmv = trans * v1buf[j];
2016  tmp_v0[j] += dmv;
2017  dam_pr += conj(tmp_v1[j]) * dmv;
2018  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
2019  }
2020  else {
2021 #pragma omp for
2022  for (j = 1; j <= X->Check.idim_max; j++) {
2023  dmv = trans * v1buf[j];
2024  dam_pr += conj(tmp_v1[j]) * dmv;
2025  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
2026  }
2027  }/*End of parallel region*/
2028  return (dam_pr);
2029 #else
2030  return 0.0;
2031 #endif
2032 }/*double complex X_GC_child_CisAit_spin_MPIdouble*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
struct EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37

◆ X_GC_child_CisAitCiuAiv_spin_MPIdouble()

double complex X_GC_child_CisAitCiuAiv_spin_MPIdouble ( int  org_isite1,
int  org_ispin1,
int  org_ispin2,
int  org_isite3,
int  org_ispin3,
int  org_ispin4,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

\(c_{is}^\dagger c_{it} c_{iu}^\dagger c_{iv}\) term in Spin model + GC. When both site1 and site2 are in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Author
Kazuyoshi Yoshimi (The University of Tokyo)
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]org_isite1site i
[in]org_ispin1spin s
[in]org_ispin2spin t
[in]org_isite3site i?
[in]org_ispin3spin u
[in]org_ispin4spin v
[in]tmp_JCopupling constatnt
[in,out]X
[in,out]tmp_v0\({\bf v}_0=H {\bf v}_1\)
[in]tmp_v1Vector to be producted

Definition at line 103 of file mltplyMPISpinCore.c.

References exitMPI(), myrank, v1buf, X, and X_GC_child_CisAis_spin_MPIdouble().

Referenced by expec_cisajscktalt_SpinGCHalf(), GC_child_CisAitCiuAiv_spin_MPIdouble(), mltplyHalfSpinGC(), and totalspin_SpinGC().

114  {
115 #ifdef MPI
116  int mask1, mask2, state1, state2, ierr, origin;
117  unsigned long int idim_max_buf, j;
118  MPI_Status statusMPI;
119  double complex Jint, dmv, dam_pr;
120 
121  mask1 = (int)X->Def.Tpow[org_isite1];
122  mask2 = (int)X->Def.Tpow[org_isite3];
123  if (org_isite1 != org_isite3) {
124  origin = myrank ^ (mask1 + mask2);
125  }
126  else {
127  if (org_ispin1 == org_ispin4 && org_ispin2 == org_ispin3) { //CisAitCitAis=CisAis
128  dam_pr = X_GC_child_CisAis_spin_MPIdouble(org_isite1, org_ispin1, tmp_J, X, tmp_v0, tmp_v1);
129  return (dam_pr);
130  }
131  else { //CisAitCisAit=0
132  return 0.0;
133  }
134  }
135 
136  state1 = (origin & mask1) / mask1;
137  state2 = (origin & mask2) / mask2;
138 
139  if (state1 == org_ispin2 && state2 == org_ispin4) {
140  Jint = tmp_J;
141  }
142  else if (state1 == org_ispin1 && state2 == org_ispin3) {
143  Jint = conj(tmp_J);
144  if (X->Large.mode == M_CORR || X->Large.mode == M_CALCSPEC) {
145  Jint = 0;
146  }
147  }
148  else {
149  return 0;
150  }
151 
152  ierr = MPI_Sendrecv(&X->Check.idim_max, 1, MPI_UNSIGNED_LONG, origin, 0,
153  &idim_max_buf, 1, MPI_UNSIGNED_LONG, origin, 0,
154  MPI_COMM_WORLD, &statusMPI);
155  if (ierr != 0) exitMPI(-1);
156  ierr = MPI_Sendrecv(tmp_v1, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
157  v1buf, idim_max_buf + 1, MPI_DOUBLE_COMPLEX, origin, 0,
158  MPI_COMM_WORLD, &statusMPI);
159  if (ierr != 0) exitMPI(-1);
160 
161  dam_pr = 0.0;
162 #pragma omp parallel default(none) reduction(+:dam_pr) private(j, dmv) \
163  firstprivate(idim_max_buf, Jint, X) shared(v1buf, tmp_v1, tmp_v0)
164  {
165  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
166 #pragma omp for
167  for (j = 1; j <= idim_max_buf; j++) {
168  dmv = Jint * v1buf[j];
169  tmp_v0[j] += dmv;
170  dam_pr += conj(tmp_v1[j]) * dmv;
171  }/*for (j = 1; j <= idim_max_buf; j++)*/
172  }
173  else {
174 #pragma omp for
175  for (j = 1; j <= idim_max_buf; j++) {
176  dmv = Jint * v1buf[j];
177  dam_pr += conj(tmp_v1[j]) * dmv;
178  }/*for (j = 1; j <= idim_max_buf; j++)*/
179  }
180  }/*End of parallel region*/
181  return dam_pr;
182 #else
183  return 0.0;
184 #endif
185 }/*void GC_child_CisAitCiuAiv_spin_MPIdouble*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
double complex X_GC_child_CisAis_spin_MPIdouble(int org_isite1, int org_ispin1, double complex tmp_trans, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
Hopping term in Spin + GC When both site1 and site2 are in the inter process region.
struct EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37

◆ X_GC_child_CisAitCiuAiv_spin_MPIsingle()

double complex X_GC_child_CisAitCiuAiv_spin_MPIsingle ( int  org_isite1,
int  org_ispin1,
int  org_ispin2,
int  org_isite3,
int  org_ispin3,
int  org_ispin4,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Exchange and Pairlifting term in Spin model + GC When only site2 is in the inter process region.

Returns
\(\langle v_1 | H_{\rm this} | v_1 \rangle\)
Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_ispin2Spin 2
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]org_ispin4Spin 4
[in]tmp_JCopupling constatnt
[in,out]X
[in,out]tmp_v0\({\bf v}_0=H {\bf v}_1\)
[in]tmp_v1Vector to be producted

Definition at line 522 of file mltplyMPISpinCore.c.

References exitMPI(), myrank, v1buf, X, and X_SpinGC_CisAit().

Referenced by expec_cisajscktalt_SpinGCHalf(), GC_child_CisAitCiuAiv_spin_MPIsingle(), mltplyHalfSpinGC(), and totalspin_SpinGC().

533  {
534 #ifdef MPI
535  int mask2, state2, ierr, origin;
536  unsigned long int mask1, idim_max_buf, j, ioff, state1, state1check;
537  MPI_Status statusMPI;
538  double complex Jint, dmv, dam_pr;
539  /*
540  Prepare index in the inter PE
541  */
542  mask2 = (int)X->Def.Tpow[org_isite3];
543  origin = myrank ^ mask2;
544  state2 = (origin & mask2) / mask2;
545 
546  if (state2 == org_ispin4) {
547  state1check = (unsigned long int) org_ispin2;
548  Jint = tmp_J;
549  }
550  else if (state2 == org_ispin3) {
551  state1check = (unsigned long int) org_ispin1;
552  Jint = conj(tmp_J);
553  if (X->Large.mode == M_CORR || X->Large.mode == M_CALCSPEC) {
554  Jint = 0;
555  }
556  }
557  else return 0.0;
558 
559  ierr = MPI_Sendrecv(&X->Check.idim_max, 1, MPI_UNSIGNED_LONG, origin, 0,
560  &idim_max_buf, 1, MPI_UNSIGNED_LONG, origin, 0,
561  MPI_COMM_WORLD, &statusMPI);
562  if (ierr != 0) exitMPI(-1);
563  ierr = MPI_Sendrecv(tmp_v1, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
564  v1buf, idim_max_buf + 1, MPI_DOUBLE_COMPLEX, origin, 0,
565  MPI_COMM_WORLD, &statusMPI);
566  if (ierr != 0) exitMPI(-1);
567  /*
568  Index in the intra PE
569  */
570  mask1 = X->Def.Tpow[org_isite1];
571 
572  dam_pr = 0.0;
573 #pragma omp parallel default(none) reduction(+:dam_pr) private(j, dmv, state1, ioff) \
574  firstprivate(idim_max_buf, Jint, X, state1check, mask1) shared(v1buf, tmp_v1, tmp_v0)
575  {
576  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
577 #pragma omp for
578  for (j = 0; j < idim_max_buf; j++) {
579  state1 = X_SpinGC_CisAit(j + 1, X, mask1, state1check, &ioff);
580  if (state1 != 0) {
581  dmv = Jint * v1buf[j + 1];
582  tmp_v0[ioff + 1] += dmv;
583  dam_pr += conj(tmp_v1[ioff + 1]) * dmv;
584  }/*if (state1 != 0)*/
585  }/*for (j = 0; j < idim_max_buf; j++)*/
586  }
587  else {
588 #pragma omp for
589  for (j = 0; j < idim_max_buf; j++) {
590  state1 = X_SpinGC_CisAit(j + 1, X, mask1, state1check, &ioff);
591  if (state1 != 0) {
592  dmv = Jint * v1buf[j + 1];
593  dam_pr += conj(tmp_v1[ioff + 1]) * dmv;
594  }/*if (state1 != 0)*/
595  }/*for (j = 0; j < idim_max_buf; j++)*/
596  }
597  }/*End of parallel region*/
598  return (dam_pr);
599 #else
600  return 0.0;
601 #endif
602 }/*void GC_child_CisAitCiuAiv_spin_MPIsingle*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
struct EDMainCalStruct X
Definition: struct.h:431
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).
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37

◆ X_GC_child_CisAitCjuAju_GeneralSpin_MPIdouble()

double complex X_GC_child_CisAitCjuAju_GeneralSpin_MPIdouble ( int  org_isite1,
int  org_ispin1,
int  org_ispin2,
int  org_isite3,
int  org_ispin3,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

\(c_{is}^\dagger c_{it} c_{ju}^\dagger c_{ju}\) term in Spin model. When both site1 and site3 are in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_ispin2Spin 2
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]tmp_JCopupling constatnt
[in,out]X
[in,out]tmp_v0\({\bf v}_0=H {\bf v}_1\)
[in]tmp_v1Vector to be producted

Definition at line 906 of file mltplyMPISpinCore.c.

References BitCheckGeneral(), exitMPI(), GetOffCompGeneralSpin(), myrank, TRUE, v1buf, and X.

Referenced by expec_cisajscktalt_SpinGCGeneral().

916  {
917 #ifdef MPI
918  unsigned long int j, off;
919  int origin, ierr;
920  double complex tmp_V, dmv, dam_pr;
921  MPI_Status statusMPI;
922 
923  if (org_isite1 == org_isite3 && org_ispin1 == org_ispin3) {//cisaitcisais=0 && cisaiscitais=0
924  return 0.0;
925  }
926 
927  if (BitCheckGeneral((unsigned long int) myrank, org_isite3 + 1, org_ispin3, X->Def.SiteToBit, X->Def.Tpow) == TRUE
928  && GetOffCompGeneralSpin((unsigned long int) myrank, org_isite1 + 1, org_ispin2, org_ispin1,
929  &off, X->Def.SiteToBit, X->Def.Tpow) == TRUE) {
930  tmp_V = tmp_J;
931  }
932  else if (GetOffCompGeneralSpin((unsigned long int) myrank, org_isite1 + 1, org_ispin1, org_ispin2, &off,
933  X->Def.SiteToBit, X->Def.Tpow) == TRUE) {
934  if (BitCheckGeneral(off, org_isite3 + 1, org_ispin3, X->Def.SiteToBit, X->Def.Tpow) == TRUE) {
935  tmp_V = conj(tmp_J);
936  if (X->Large.mode == M_CORR || X->Large.mode == M_CALCSPEC) tmp_V = 0.0;
937  }
938  else return 0.0;
939  }
940  else return 0.0;
941 
942  origin = (int)off;
943 
944  ierr = MPI_Sendrecv(tmp_v1, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
945  v1buf, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
946  MPI_COMM_WORLD, &statusMPI);
947  if (ierr != 0) exitMPI(-1);
948 
949  dam_pr = 0.0;
950 #pragma omp parallel default(none) reduction(+:dam_pr) firstprivate(X, tmp_V) private(j, dmv) \
951 shared (tmp_v0, tmp_v1, v1buf)
952  {
953  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
954 #pragma omp for
955  for (j = 1; j <= X->Check.idim_max; j++) {
956  dmv = v1buf[j] * tmp_V;
957  tmp_v0[j] += dmv;
958  dam_pr += conj(tmp_v1[j]) * dmv;
959  }
960  }
961  else {
962 #pragma omp for
963  for (j = 1; j <= X->Check.idim_max; j++) {
964  dmv = v1buf[j] * tmp_V;
965  dam_pr += conj(tmp_v1[j]) * dmv;
966  }
967  }
968  }/*End of parallel region*/
969  return dam_pr;
970 #else
971  return 0.0;
972 #endif
973 }/*double complex X_GC_child_CisAitCjuAju_GeneralSpin_MPIdouble*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
#define TRUE
Definition: global.h:26
int BitCheckGeneral(const long unsigned int org_bit, const unsigned int org_isite, const unsigned int target_ispin, const long int *SiteToBit, const long unsigned int *Tpow)
bit check function for general spin
Definition: bitcalc.c:393
struct EDMainCalStruct X
Definition: struct.h:431
int GetOffCompGeneralSpin(const long unsigned int org_ibit, const int org_isite, const int org_ispin, const int off_ispin, long unsigned int *_ioffComp, const long int *SiteToBit, const long unsigned int *Tpow)
function of getting off-diagonal component for general spin
Definition: bitcalc.c:243
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37

◆ X_GC_child_CisAitCjuAju_GeneralSpin_MPIsingle()

double complex X_GC_child_CisAitCjuAju_GeneralSpin_MPIsingle ( int  org_isite1,
int  org_ispin1,
int  org_ispin2,
int  org_isite3,
int  org_ispin3,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Compute \(c_{is}^\dagger c_{it}c_{ju}^\dagger c_{ju}\) term in the grandcanonical general spin system when one of these site is in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_ispin2Spin 2
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]tmp_JCoupling constant
[in,out]X
[in,out]tmp_v0Resulting wavefunction
[in]tmp_v1Input wavefunction

Definition at line 1443 of file mltplyMPISpinCore.c.

References BitCheckGeneral(), GetOffCompGeneralSpin(), myrank, TRUE, and X.

Referenced by expec_cisajscktalt_SpinGCGeneral(), GC_child_general_int_GeneralSpin_MPIdouble(), and GC_child_general_int_GeneralSpin_MPIsingle().

1453  {
1454 #ifdef MPI
1455  unsigned long int num1, j, off;
1456  int isite, IniSpin, FinSpin;
1457  double complex tmp_V, dmv, dam_pr;
1458  //MPI_Status statusMPI;
1459 
1460  num1 = BitCheckGeneral((unsigned long int)myrank,
1461  org_isite3+1, org_ispin3, X->Def.SiteToBit, X->Def.Tpow);
1462  if(num1 != 0){
1463  tmp_V = tmp_J;
1464  isite = org_isite1 + 1;
1465  IniSpin = org_ispin2;
1466  FinSpin = org_ispin1;
1467  }
1468  else return 0.0;
1469 
1470  dam_pr = 0.0;
1471 #pragma omp parallel default(none) reduction(+:dam_pr) \
1472 firstprivate(X, tmp_V, isite, IniSpin, FinSpin) private(j, dmv, num1, off) \
1473 shared (tmp_v0, tmp_v1, v1buf)
1474  {
1475  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
1476 #pragma omp for
1477  for (j = 1; j <= X->Check.idim_max; j++) {
1478  if (GetOffCompGeneralSpin(j - 1, isite, IniSpin, FinSpin, &off,
1479  X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1480  {
1481  dmv = tmp_v1[j] * tmp_V;
1482  tmp_v0[off + 1] += dmv;
1483  dam_pr += conj(tmp_v1[off + 1]) * dmv;
1484  }
1485  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1486  }
1487  else {
1488 #pragma omp for
1489  for (j = 1; j <= X->Check.idim_max; j++) {
1490  if (GetOffCompGeneralSpin(j - 1, isite, IniSpin, FinSpin, &off,
1491  X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1492  {
1493  dmv = tmp_v1[j] * tmp_V;
1494  dam_pr += conj(tmp_v1[off + 1]) * dmv;
1495  }
1496  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1497  }
1498  }/*End of parallel region*/
1499  return dam_pr;
1500 #else
1501  return 0.0;
1502 #endif
1503 }/*double complex X_GC_child_CisAitCjuAju_GeneralSpin_MPIsingle*/
#define TRUE
Definition: global.h:26
int BitCheckGeneral(const long unsigned int org_bit, const unsigned int org_isite, const unsigned int target_ispin, const long int *SiteToBit, const long unsigned int *Tpow)
bit check function for general spin
Definition: bitcalc.c:393
struct EDMainCalStruct X
Definition: struct.h:431
int GetOffCompGeneralSpin(const long unsigned int org_ibit, const int org_isite, const int org_ispin, const int off_ispin, long unsigned int *_ioffComp, const long int *SiteToBit, const long unsigned int *Tpow)
function of getting off-diagonal component for general spin
Definition: bitcalc.c:243
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162

◆ X_GC_child_CisAitCjuAju_spin_MPIdouble()

double complex X_GC_child_CisAitCjuAju_spin_MPIdouble ( int  org_isite1,
int  org_ispin1,
int  org_ispin2,
int  org_isite3,
int  org_ispin3,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

CisAisCjuAjv term in Spin model + GC When both site1 and site2 are in the inter process region.

Returns
\(\langle v_1 | H_{\rm this} | v_1 \rangle\)
Author
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_ispin2Spin 2
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]tmp_JCopupling constatnt
[in,out]X
[in,out]tmp_v0\({\bf v}_0=H {\bf v}_1\)
[in]tmp_v1Vector to be producted

Definition at line 309 of file mltplyMPISpinCore.c.

References exitMPI(), myrank, v1buf, X, and X_SpinGC_CisAis().

Referenced by expec_cisajscktalt_SpinGCHalf(), and GC_child_CisAitCjuAju_spin_MPIdouble().

319  {
320 #ifdef MPI
321  int mask1, mask2, state1, ierr, num1;
322  long int origin;
323  unsigned long int idim_max_buf, j;
324  MPI_Status statusMPI;
325  double complex Jint, dmv, dam_pr;
326 
327  if (org_isite1 == org_isite3 && org_ispin1 == org_ispin3) {//cisaitcisais
328  return 0.0;
329  }
330 
331  mask1 = (int)X->Def.Tpow[org_isite1];
332  origin = myrank ^ mask1;
333  state1 = (origin & mask1) / mask1;
334  mask2 = (int)X->Def.Tpow[org_isite3];
335  num1 = X_SpinGC_CisAis(origin + 1, X, mask2, org_ispin3);
336  if (state1 == org_ispin2) {
337  if (num1 != 0) {
338  Jint = tmp_J;
339  }
340  else {
341  return 0.0;
342  }
343  }/*if (state1 == org_ispin2)*/
344  else {//state1 = org_ispin1
345  num1 = X_SpinGC_CisAis((unsigned long int) myrank + 1, X, mask2, org_ispin3);
346  if (num1 != 0) {
347  Jint = conj(tmp_J);
348  if (X->Large.mode == M_CORR || X->Large.mode == M_CALCSPEC) {
349  Jint = 0;
350  }
351  }
352  else {
353  return 0.0;
354  }
355  }
356 
357  ierr = MPI_Sendrecv(&X->Check.idim_max, 1, MPI_UNSIGNED_LONG, origin, 0,
358  &idim_max_buf, 1, MPI_UNSIGNED_LONG, origin, 0,
359  MPI_COMM_WORLD, &statusMPI);
360  if (ierr != 0) exitMPI(-1);
361  ierr = MPI_Sendrecv(tmp_v1, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
362  v1buf, idim_max_buf + 1, MPI_DOUBLE_COMPLEX, origin, 0,
363  MPI_COMM_WORLD, &statusMPI);
364  if (ierr != 0) exitMPI(-1);
365 
366  dam_pr = 0.0;
367 #pragma omp parallel default(none) reduction(+:dam_pr) private(j, dmv) \
368  firstprivate(idim_max_buf, Jint, X) shared(v1buf, tmp_v1, tmp_v0)
369  {
370  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
371 #pragma omp for
372  for (j = 1; j <= idim_max_buf; j++) {
373  dmv = Jint * v1buf[j];
374  tmp_v0[j] += dmv;
375  dam_pr += conj(tmp_v1[j]) * dmv;
376  }/*for (j = 1; j <= idim_max_buf; j++)*/
377  }
378  else {
379 #pragma omp for
380  for (j = 1; j <= idim_max_buf; j++) {
381  dmv = Jint * v1buf[j];
382  dam_pr += conj(tmp_v1[j]) * dmv;
383  }/*for (j = 1; j <= idim_max_buf; j++)*/
384  }
385  }/*End of parallel region*/
386  return (dam_pr);
387 #else
388  return 0.0;
389 #endif
390 }/*double complex X_GC_child_CisAisCjuAjv_spin_MPIdouble*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
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 EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37

◆ X_GC_child_CisAitCjuAju_spin_MPIsingle()

double complex X_GC_child_CisAitCjuAju_spin_MPIsingle ( int  org_isite1,
int  org_ispin1,
int  org_ispin2,
int  org_isite3,
int  org_ispin3,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

CisAisCjuAjv term in Spin model + GC When only site2 is in the inter process region.

Returns
\(\langle v_1 | H_{\rm this} | v_1 \rangle\)
Author
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_ispin2Spin 2
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]tmp_JCopupling constatnt
[in,out]X
[in,out]tmp_v0\({\bf v}_0=H {\bf v}_1\)
[in]tmp_v1Vector to be producted

Definition at line 732 of file mltplyMPISpinCore.c.

References myrank, and X.

Referenced by expec_cisajscktalt_SpinGCHalf(), and GC_child_CisAitCjuAju_spin_MPIsingle().

742  {
743 #ifdef MPI
744  int mask2, state2;
745  unsigned long int mask1, j, ioff, state1, state1check;
746  //MPI_Status statusMPI;
747  double complex Jint, dmv, dam_pr;
748  /*
749  Prepare index in the inter PE
750  */
751  mask2 = (int)X->Def.Tpow[org_isite3];
752  state2 = (myrank & mask2) / mask2;
753 
754  if (state2 == org_ispin3) {
755  state1check = org_ispin2;
756  Jint = tmp_J;
757  }
758  else {
759  return 0.0;
760  }
761 
762  mask1 = (int)X->Def.Tpow[org_isite1];
763 
764  dam_pr = 0.0;
765 #pragma omp parallel default(none) reduction(+:dam_pr) private(j, dmv, state1, ioff) \
766  firstprivate(Jint, X, state1check, mask1) shared(tmp_v1, tmp_v0)
767  {
768  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
769 #pragma omp for
770  for (j = 0; j < X->Check.idim_max; j++) {
771 
772  state1 = (j & mask1) / mask1;
773  ioff = j ^ mask1;
774  if (state1 == state1check) {
775  dmv = Jint * tmp_v1[j + 1];
776  }
777  else {
778  dmv = conj(Jint) * tmp_v1[j + 1];
779  }
780  tmp_v0[ioff + 1] += dmv;
781  dam_pr += conj(tmp_v1[ioff + 1]) * dmv;
782  }/*for (j = 0; j < X->Check.idim_max; j++)*/
783  }
784  else if (X->Large.mode == M_CORR) {
785 #pragma omp for
786  for (j = 0; j < X->Check.idim_max; j++) {
787 
788  state1 = (j & mask1) / mask1;
789  ioff = j ^ mask1;
790  if (state1 == state1check) {
791  dmv = Jint * tmp_v1[j + 1];
792  }
793  else {
794  dmv = 0.0;
795  }
796  dam_pr += conj(tmp_v1[ioff + 1]) * dmv;
797  }/*for (j = 0; j < X->Check.idim_max; j++)*/
798  }
799  else {
800 #pragma omp for
801  for (j = 0; j < X->Check.idim_max; j++) {
802  state1 = (j & mask1) / mask1;
803  ioff = j ^ mask1;
804  if (state1 == state1check) {
805  dmv = Jint * tmp_v1[j + 1];
806  }
807  else {
808  dmv = conj(Jint) * tmp_v1[j + 1];
809  }
810  dam_pr += conj(tmp_v1[ioff + 1]) * dmv;
811  }/*for (j = 0; j < X->Check.idim_max; j++)*/
812  }
813  }/*End of parallel region*/
814  return (dam_pr);
815 #else
816  return 0.0;
817 #endif
818 }/*void GC_child_CisAitCiuAiv_spin_MPIsingle*/
struct EDMainCalStruct X
Definition: struct.h:431
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162

◆ X_GC_child_CisAitCjuAjv_GeneralSpin_MPIdouble()

double complex X_GC_child_CisAitCjuAjv_GeneralSpin_MPIdouble ( int  org_isite1,
int  org_ispin1,
int  org_ispin2,
int  org_isite3,
int  org_ispin3,
int  org_ispin4,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Compute \(c_{is}^\dagger c_{it} c_{ju}^\dagger c_{jv}\) term in the grandcanonical general spin system when both site is in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_ispin2Spin 2
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]org_ispin4Spin 4
[in]tmp_JCoupling constant
[in,out]X
[in,out]tmp_v0Resulting wavefunction
[in]tmp_v1Input wavefunction

Definition at line 979 of file mltplyMPISpinCore.c.

References exitMPI(), FALSE, GetOffCompGeneralSpin(), myrank, TRUE, v1buf, X, and X_GC_child_CisAis_GeneralSpin_MPIdouble().

Referenced by expec_cisajscktalt_SpinGCGeneral(), and GC_child_general_int_GeneralSpin_MPIdouble().

990  {
991 #ifdef MPI
992  unsigned long int tmp_off, off, j;
993  int origin, ierr, ihermite;
994  double complex tmp_V, dmv, dam_pr;
995  MPI_Status statusMPI;
996 
997  ihermite = TRUE;
998 
999  if (org_isite1 == org_isite3 && org_ispin1 == org_ispin4 &&
1000  org_ispin2 == org_ispin3) { //cisaitcitais=cisais && cisaitcitais =cisais
1001  dam_pr = X_GC_child_CisAis_GeneralSpin_MPIdouble(org_isite1, org_ispin1, tmp_J, X, tmp_v0, tmp_v1);
1002  return (dam_pr);
1003  }
1004  //cisaitcisait
1005  if (GetOffCompGeneralSpin((unsigned long int) myrank, org_isite1 + 1, org_ispin1, org_ispin2,
1006  &tmp_off, X->Def.SiteToBit, X->Def.Tpow) == TRUE) {
1007 
1008  if (GetOffCompGeneralSpin(tmp_off, org_isite3 + 1, org_ispin3, org_ispin4,
1009  &off, X->Def.SiteToBit, X->Def.Tpow) == TRUE) {
1010 
1011  tmp_V = tmp_J;
1012  }
1013  else ihermite = FALSE;
1014  }
1015  else {
1016  ihermite = FALSE;
1017  }
1018 
1019  if (ihermite == FALSE) {
1020  if (GetOffCompGeneralSpin((unsigned long int) myrank, org_isite3 + 1, org_ispin4, org_ispin3, &tmp_off,
1021  X->Def.SiteToBit, X->Def.Tpow) == TRUE) {
1022 
1023  if (GetOffCompGeneralSpin(tmp_off, org_isite1 + 1, org_ispin2, org_ispin1, &off, X->Def.SiteToBit,
1024  X->Def.Tpow) == TRUE) {
1025  tmp_V = conj(tmp_J);
1026  if (X->Large.mode == M_CORR || X->Large.mode == M_CALCSPEC) tmp_V = 0.0;
1027  }
1028  else return 0.0;
1029  }
1030  else return 0.0;
1031  }
1032 
1033  origin = (int)off;
1034 
1035  ierr = MPI_Sendrecv(tmp_v1, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
1036  v1buf, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
1037  MPI_COMM_WORLD, &statusMPI);
1038  if (ierr != 0) exitMPI(-1);
1039 
1040  dam_pr = 0.0;
1041 #pragma omp parallel default(none) reduction(+:dam_pr) firstprivate(X, tmp_V) private(j, dmv) \
1042  shared (tmp_v0, tmp_v1, v1buf)
1043  {
1044  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
1045 #pragma omp for
1046  for (j = 1; j <= X->Check.idim_max; j++) {
1047  dmv = v1buf[j] * tmp_V;
1048  tmp_v0[j] += dmv;
1049  dam_pr += conj(tmp_v1[j]) * dmv;
1050  }
1051  }
1052  else {
1053 #pragma omp for
1054  for (j = 1; j <= X->Check.idim_max; j++) {
1055  dmv = v1buf[j] * tmp_V;
1056  dam_pr += conj(tmp_v1[j]) * dmv;
1057  }
1058  }
1059  }/*End of parallel region*/
1060  return dam_pr;
1061 #else
1062  return 0.0;
1063 #endif
1064 }/*double complex X_GC_child_CisAitCjuAjv_GeneralSpin_MPIdouble*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
double complex X_GC_child_CisAis_GeneralSpin_MPIdouble(int org_isite1, int org_ispin1, double complex tmp_trans, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
Compute term in the grandcanonical general spin system when both site is in the inter process region...
#define TRUE
Definition: global.h:26
#define FALSE
Definition: global.h:25
struct EDMainCalStruct X
Definition: struct.h:431
int GetOffCompGeneralSpin(const long unsigned int org_ibit, const int org_isite, const int org_ispin, const int off_ispin, long unsigned int *_ioffComp, const long int *SiteToBit, const long unsigned int *Tpow)
function of getting off-diagonal component for general spin
Definition: bitcalc.c:243
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37

◆ X_GC_child_CisAitCjuAjv_GeneralSpin_MPIsingle()

double complex X_GC_child_CisAitCjuAjv_GeneralSpin_MPIsingle ( int  org_isite1,
int  org_ispin1,
int  org_ispin2,
int  org_isite3,
int  org_ispin3,
int  org_ispin4,
double complex  tmp_J,
struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Compute \(c_{is}^\dagger c_{is}c_{ju}^\dagger c_{jv}\) term in the grandcanonical general spin system when one of these site is in the inter process region.

Returns
\(\langle v_1| H_{\rm this} | v_1 \rangle\)
Parameters
[in]org_isite1Site 1
[in]org_ispin1Spin 1
[in]org_ispin2Spin 2
[in]org_isite3Site 3
[in]org_ispin3Spin 3
[in]org_ispin4Spin 4
[in]tmp_JCoupling constant
[in,out]X
[in,out]tmp_v0Resulting wavefunction
[in]tmp_v1Input wavefunction

Definition at line 1509 of file mltplyMPISpinCore.c.

References exitMPI(), GetOffCompGeneralSpin(), myrank, TRUE, v1buf, and X.

Referenced by expec_cisajscktalt_SpinGCGeneral(), and GC_child_general_int_GeneralSpin_MPIsingle().

1520  {
1521 #ifdef MPI
1522  unsigned long int off, j;
1523  int origin, ierr, isite, IniSpin, FinSpin;
1524  double complex tmp_V, dmv, dam_pr;
1525  MPI_Status statusMPI;
1526 
1527  if (GetOffCompGeneralSpin((unsigned long int)myrank,
1528  org_isite3 + 1, org_ispin3, org_ispin4, &off,
1529  X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1530  {
1531  tmp_V = tmp_J;
1532  isite = org_isite1 + 1;
1533  IniSpin = org_ispin2;
1534  FinSpin = org_ispin1;
1535  }
1536  else if (GetOffCompGeneralSpin((unsigned long int)myrank,
1537  org_isite3 + 1, org_ispin4, org_ispin3, &off,
1538  X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1539  {
1540  tmp_V = conj(tmp_J);
1541  if (X->Large.mode == M_CORR || X->Large.mode == M_CALCSPEC) tmp_V = 0.0;
1542  isite = org_isite1 + 1;
1543  IniSpin = org_ispin1;
1544  FinSpin = org_ispin2;
1545  }
1546  else return 0.0;
1547 
1548  origin = (int)off;
1549 
1550  ierr = MPI_Sendrecv(tmp_v1, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
1551  v1buf, X->Check.idim_max + 1, MPI_DOUBLE_COMPLEX, origin, 0,
1552  MPI_COMM_WORLD, &statusMPI);
1553  if (ierr != 0) exitMPI(-1);
1554 
1555  dam_pr = 0.0;
1556 #pragma omp parallel default(none) reduction(+:dam_pr) \
1557 firstprivate(X, tmp_V, isite, IniSpin, FinSpin) private(j, dmv, off) shared (tmp_v0, tmp_v1, v1buf)
1558  {
1559  if (X->Large.mode == M_MLTPLY || X->Large.mode == M_CALCSPEC) {
1560 #pragma omp for
1561  for (j = 1; j <= X->Check.idim_max; j++) {
1562  if (GetOffCompGeneralSpin(j - 1, isite, IniSpin, FinSpin, &off,
1563  X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1564  {
1565  dmv = v1buf[j] * tmp_V;
1566  tmp_v0[off + 1] += dmv;
1567  dam_pr += conj(tmp_v1[off + 1]) * dmv;
1568  }
1569  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1570  }
1571  else {
1572 #pragma omp for
1573  for (j = 1; j <= X->Check.idim_max; j++) {
1574  if (GetOffCompGeneralSpin(j - 1, isite, IniSpin, FinSpin, &off,
1575  X->Def.SiteToBit, X->Def.Tpow) == TRUE)
1576  {
1577  dmv = v1buf[j] * tmp_V;
1578  dam_pr += conj(tmp_v1[off + 1]) * dmv;
1579  }
1580  }/*for (j = 1; j <= X->Check.idim_max; j++)*/
1581  }
1582  }/*End of parallel region*/
1583  return dam_pr;
1584 #else
1585  return 0.0;
1586 #endif
1587 }/*double complex X_GC_child_CisAitCjuAjv_GeneralSpin_MPIsingle*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
#define TRUE
Definition: global.h:26
struct EDMainCalStruct X
Definition: struct.h:431
int GetOffCompGeneralSpin(const long unsigned int org_ibit, const int org_isite, const int org_ispin, const int off_ispin, long unsigned int *_ioffComp, const long int *SiteToBit, const long unsigned int *Tpow)
function of getting off-diagonal component for general spin
Definition: bitcalc.c:243
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
double complex * v1buf
Definition: global.h:37