HΦ  3.2.0
mltply.c File Reference

Multiplying the wavefunction by the Hamiltonian. \( H v_1\). More...

#include <bitcalc.h>
#include "mltply.h"
#include "mltplySpin.h"
#include "mltplyHubbard.h"
#include "wrapperMPI.h"
#include "CalcTime.h"
#include "mltplyCommon.h"
#include "diagonalcalc.h"
+ Include dependency graph for mltply.c:

Go to the source code of this file.

Functions

int mltply (struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
 Parent function of multiplying the wavefunction by the Hamiltonian. \( H v_1\).
First, the calculation of diagonal term is done by using the list \( \verb|list_diaognal| \).
Next, the calculation of off-diagonal term is done.
. More...
 

Detailed Description

Multiplying the wavefunction by the Hamiltonian. \( H v_1\).

Version
0.2

add function to treat the case of generalspin

Version
0.1
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)

Definition in file mltply.c.

Function Documentation

◆ mltply()

int mltply ( struct BindStruct X,
double complex *  tmp_v0,
double complex *  tmp_v1 
)

Parent function of multiplying the wavefunction by the Hamiltonian. \( H v_1\).
First, the calculation of diagonal term is done by using the list \( \verb|list_diaognal| \).
Next, the calculation of off-diagonal term is done.
.

Note
If \( \verb|mode| \) in BindStruct X is \( \verb|M_CORR| \), the wave function is not updated. The expected values are only calculated.
Otherwise, the wavefunction \( v_0 \) is updated as \( v_0 += H v_1\).
Parameters
X[in] Struct for getting the information of the operators.
tmp_v0[in, out]
tmp_v1[in]
Returns
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)

Definition at line 56 of file mltply.c.

References BindStruct::Check, BindStruct::Def, diagonalcalcForTE(), FALSE, GetSplitBitByModel(), GetSplitBitForGeneralSpin(), LargeList::i_max, DefineList::iCalcModel, DefineList::iCalcType, CheckList::idim_max, DefineList::iFlgGeneralSpin, LargeList::ihfbit, LargeList::ilft, LargeList::irght, BindStruct::Large, list_Diagonal, mltplyHubbard(), mltplyHubbardGC(), mltplySpin(), mltplySpinGC(), LargeList::mode, DefineList::Nsite, LargeList::prdct, DefineList::SiteToBit, StartTimer(), step_i, StopTimer(), and SumMPI_dc().

Referenced by CalcSpectrumByBiCG(), CG_EigenVector(), expec_energy_flct(), Lanczos_EigenValue(), Lanczos_EigenVector(), Lanczos_GetTridiagonalMatrixComponents(), LOBPCG_Main(), MultiplyForTEM(), and PowerLanczos().

56  {
57  long unsigned int j=0;
58  long unsigned int irght=0;
59  long unsigned int ilft=0;
60  long unsigned int ihfbit=0;
61 
62  double complex dam_pr;
63 
64  long unsigned int i_max;
65 
66  StartTimer(1);
67  i_max = X->Check.idim_max;
68  X->Large.prdct = 0.0;
69  dam_pr = 0.0;
70 
71  if(i_max!=0){
72  if (X->Def.iFlgGeneralSpin == FALSE) {
73  if (GetSplitBitByModel(X->Def.Nsite, X->Def.iCalcModel, &irght, &ilft, &ihfbit) != 0) {
74  return -1;
75  }
76  }
77  else{
78  if(X->Def.iCalcModel==Spin){
79  if (GetSplitBitForGeneralSpin(X->Def.Nsite, &ihfbit, X->Def.SiteToBit) != 0) {
80  return -1;
81  }
82  }
83  }
84  }
85  else{
86  irght=0;
87  ilft=0;
88  ihfbit=0;
89  }
90  X->Large.i_max = i_max;
91  X->Large.irght = irght;
92  X->Large.ilft = ilft;
93  X->Large.ihfbit = ihfbit;
94  X->Large.mode = M_MLTPLY;
95 
96  StartTimer(100);
97 #pragma omp parallel for default(none) reduction(+:dam_pr) firstprivate(i_max) shared(tmp_v0, tmp_v1, list_Diagonal)
98  for (j = 1; j <= i_max; j++) {
99  tmp_v0[j] += (list_Diagonal[j]) * tmp_v1[j];
100  dam_pr += (list_Diagonal[j]) * conj(tmp_v1[j]) * tmp_v1[j];
101  }
102  X->Large.prdct += dam_pr;
103  StopTimer(100);
104  if (X->Def.iCalcType == TimeEvolution) diagonalcalcForTE(step_i, X, tmp_v0, tmp_v1);
105 
106  switch (X->Def.iCalcModel) {
107  case HubbardGC:
108  mltplyHubbardGC(X, tmp_v0, tmp_v1);
109  break;
110 
111  case KondoGC:
112  case Hubbard:
113  case Kondo:
114  mltplyHubbard(X, tmp_v0, tmp_v1);
115  break;
116 
117  case Spin:
118  mltplySpin(X, tmp_v0, tmp_v1);
119  break;
120 
121  case SpinGC:
122  mltplySpinGC(X, tmp_v0, tmp_v1);
123  break;
124 
125  default:
126  return -1;
127  }
128 
129  X->Large.prdct = SumMPI_dc(X->Large.prdct);
130  StopTimer(1);
131  return 0;
132 }
struct DefineList Def
Definision of system (Hamiltonian) etc.
Definition: struct.h:410
void StartTimer(int n)
function for initializing elapse time [start]
Definition: time.c:71
int GetSplitBitForGeneralSpin(const int Nsite, long unsigned int *ihfbit, const long int *SiteToBit)
function of getting right, left and half bits corresponding to a original space.
Definition: bitcalc.c:124
double complex SumMPI_dc(double complex norm)
MPI wrapper function to obtain sum of Double complex across processes.
Definition: wrapperMPI.c:205
long unsigned int ihfbit
Used for Ogata-Lin ???
Definition: struct.h:345
unsigned long int idim_max
The dimension of the Hilbert space of this process.
Definition: struct.h:303
void StopTimer(int n)
function for calculating elapse time [elapse time=StartTimer-StopTimer]
Definition: time.c:83
int mltplySpinGC(struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
Driver function for Spin hamiltonian.
Definition: mltplySpin.c:393
double complex prdct
The expectation value of the energy.
Definition: struct.h:314
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.h:412
int mode
multiply or expectation value.
Definition: struct.h:330
unsigned int Nsite
Number of sites in the INTRA process region.
Definition: struct.h:56
long unsigned int ilft
Used for Ogata-Lin ???
Definition: struct.h:344
int GetSplitBitByModel(const int Nsite, const int iCalcModel, long unsigned int *irght, long unsigned int *ilft, long unsigned int *ihfbit)
function of splitting original bit into right and left spaces.
Definition: bitcalc.c:78
long int i_max
Length of eigenvector.
Definition: struct.h:317
int iFlgGeneralSpin
Flag for the general (Sz/=1/2) spin.
Definition: struct.h:86
#define FALSE
Definition: global.h:25
long int * SiteToBit
[DefineList::NsiteMPI] Similar to DefineList::Tpow. For general spin.
Definition: struct.h:94
int mltplyHubbardGC(struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
perform Hamiltonian vector product for (extended) Hubbard type model (Grandcanonical).
long unsigned int irght
Used for Ogata-Lin ???
Definition: struct.h:343
double * list_Diagonal
Definition: global.h:46
int iCalcModel
Switch for model. 0:Hubbard, 1:Spin, 2:Kondo, 3:HubbardGC, 4:SpinGC, 5:KondoGC, 6:HubbardNConserved.
Definition: struct.h:198
int diagonalcalcForTE(const int _istep, struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
Definition: diagonalcalc.c:197
int mltplySpin(struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
Driver function for Spin hamiltonian.
Definition: mltplySpin.c:175
struct CheckList Check
Size of the Hilbert space.
Definition: struct.h:411
int step_i
Definition: global.h:66
int mltplyHubbard(struct BindStruct *X, double complex *tmp_v0, double complex *tmp_v1)
perform Hamiltonian vector product for (extended) Hubbard type model.
int iCalcType
Switch for calculation type. 0:Lanczos, 1:TPQCalc, 2:FullDiag.
Definition: struct.h:192
+ Here is the call graph for this function:
+ Here is the caller graph for this function: