HPhi++
3.1.0
|
Multiplying the wavefunction by the Hamiltonian. \( H v_1\). More...
#include <bitcalc.hpp>
#include "mltply.hpp"
#include "mltplySpin.hpp"
#include "mltplyHubbard.hpp"
#include "wrapperMPI.hpp"
#include "CalcTime.hpp"
#include "mltplyCommon.hpp"
#include "diagonalcalc.hpp"
Go to the source code of this file.
Functions | |
int | mltply (struct BindStruct *X, int nstate, std::complex< double > **tmp_v0, std::complex< double > **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... | |
void | zaxpy_long (long int n, std::complex< double > a, std::complex< double > *x, std::complex< double > *y) |
Wrapper of zaxpy. More... | |
void | zclear (long int n, std::complex< double > *x) |
clear std::complex<double> array. More... | |
Multiplying the wavefunction by the Hamiltonian. \( H v_1\).
add function to treat the case of generalspin
Definition in file mltply.cpp.
int mltply | ( | struct BindStruct * | X, |
int | nstate, | ||
std::complex< double > ** | tmp_v0, | ||
std::complex< double > ** | 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.
.
X | [in] Struct for getting the information of the operators. |
tmp_v0 | [in, out] |
tmp_v1 | [in] |
Definition at line 56 of file mltply.cpp.
References BindStruct::Check, BindStruct::Def, diagonalcalcForTE(), 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, DefineList::SiteToBit, StartTimer(), step_i, and StopTimer().
Referenced by CalcByFullDiag(), CalcSpectrumByBiCG(), CalcSpectrumByFullDiag(), expec_energy_flct(), LOBPCG_Main(), and MultiplyForTEM().
void zaxpy_long | ( | long int | n, |
std::complex< double > | a, | ||
std::complex< double > * | x, | ||
std::complex< double > * | y | ||
) |
Wrapper of zaxpy.
Definition at line 128 of file mltply.cpp.
Referenced by expec_cisajs_Hubbard(), expec_cisajs_HubbardGC(), expec_cisajs_SpinGeneral(), expec_cisajs_SpinHalf(), expec_cisajscktalt_SpinHalf(), GetPairExcitedStateHubbard(), GetPairExcitedStateHubbardGC(), X_child_CisAis_Hubbard_MPI(), X_child_CisAisCjtAjt_Hubbard_MPI(), X_child_CisAisCjuAju_GeneralSpin_MPIdouble(), X_GC_Ajt_MPI(), X_GC_child_AisCis_GeneralSpin_MPIdouble(), X_GC_child_AisCis_spin_MPIdouble(), X_GC_child_CisAis_GeneralSpin_MPIdouble(), X_GC_child_CisAis_Hubbard_MPI(), X_GC_child_CisAis_spin_MPIdouble(), X_GC_child_CisAisCjtAjt_Hubbard_MPI(), X_GC_child_CisAisCjuAju_GeneralSpin_MPIdouble(), X_GC_child_CisAisCjuAjv_GeneralSpin_MPIdouble(), X_GC_child_CisAisCjuAjv_spin_MPIdouble(), X_GC_child_CisAit_GeneralSpin_MPIdouble(), X_GC_child_CisAit_spin_MPIdouble(), X_GC_child_CisAitCiuAiv_spin_MPIdouble(), X_GC_child_CisAitCjuAju_GeneralSpin_MPIdouble(), X_GC_child_CisAitCjuAju_spin_MPIdouble(), X_GC_child_CisAitCjuAjv_GeneralSpin_MPIdouble(), X_GC_child_CisAjtCkuAlv_Hubbard_MPI(), X_GC_child_general_hopp_MPIdouble(), and X_GC_Cis_MPI().
void zclear | ( | long int | n, |
std::complex< double > * | x | ||
) |
clear std::complex<double> array.
Definition at line 143 of file mltply.cpp.
Referenced by CalcByFullDiag(), CalcSpectrumByBiCG(), CalcSpectrumByFullDiag(), expec_cisajs_Hubbard(), expec_cisajs_HubbardGC(), expec_cisajs_SpinGCGeneral(), expec_cisajs_SpinGCHalf(), expec_cisajs_SpinGeneral(), expec_cisajs_SpinHalf(), expec_cisajscktalt_Hubbard(), expec_cisajscktalt_HubbardGC(), expec_cisajscktalt_SpinGCGeneral(), expec_cisajscktalt_SpinGCHalf(), expec_cisajscktalt_SpinGeneral(), expec_cisajscktalt_SpinHalf(), and LOBPCG_Main().