HPhi++  3.1.0
mltply.cpp
Go to the documentation of this file.
1 /* HPhi - Quantum Lattice Model Simulator */
2 /* Copyright (C) 2015 The University of Tokyo */
3 /* This program is free software: you can redistribute it and/or modify */
4 /* it under the terms of the GNU General Public License as published by */
5 /* the Free Software Foundation, either version 3 of the License, or */
6 /* (at your option) any later version. */
7 
8 /* This program is distributed in the hope that it will be useful, */
9 /* but WITHOUT ANY WARRANTY; without even the implied warranty of */
10 /* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the */
11 /* GNU General Public License for more details. */
12 
13 /* You should have received a copy of the GNU General Public License */
14 /* along with this program. If not, see <http://www.gnu.org/licenses/>. */
15 
16 // Define Mode for mltply
17 // complex version
18 #include <bitcalc.hpp>
19 #include "mltply.hpp"
20 #include "mltplySpin.hpp"
21 #include "mltplyHubbard.hpp"
22 #include "wrapperMPI.hpp"
23 #include "CalcTime.hpp"
24 #include "mltplyCommon.hpp"
25 #include "diagonalcalc.hpp"
26 
56 int mltply(struct BindStruct *X, int nstate, std::complex<double> **tmp_v0,std::complex<double> **tmp_v1) {
57  int one = 1;
58  long int j=0;
59  long int irght=0;
60  long int ilft=0;
61  long int ihfbit=0;
62  std::complex<double> dmv;
63 
64 
65  long int i_max;
66 
67  StartTimer(1);
68  i_max = X->Check.idim_max;
69 
70  if (X->Def.iFlgGeneralSpin == FALSE) {
71  if (GetSplitBitByModel(X->Def.Nsite, X->Def.iCalcModel, &irght, &ilft, &ihfbit) != 0) {
72  return -1;
73  }
74  }
75  else {
76  if (X->Def.iCalcModel == Spin) {
77  if (GetSplitBitForGeneralSpin(X->Def.Nsite, &ihfbit, X->Def.SiteToBit) != 0) {
78  return -1;
79  }
80  }
81  }
82 
83  X->Large.i_max = i_max;
84  X->Large.irght = irght;
85  X->Large.ilft = ilft;
86  X->Large.ihfbit = ihfbit;
87  X->Large.mode = M_MLTPLY;
88 
89  StartTimer(100);
90 #pragma omp parallel for default(none) private(dmv) \
91  firstprivate(i_max) shared(tmp_v0, tmp_v1, list_Diagonal,one,nstate)
92  for (j = 1; j <= i_max; j++) {
93  dmv = list_Diagonal[j];
94  zaxpy_(&nstate, &dmv, &tmp_v1[j][0], &one, &tmp_v0[j][0], &one);
95  }
96  StopTimer(100);
97  if (X->Def.iCalcType == TimeEvolution) diagonalcalcForTE(step_i, X, &tmp_v0[0][0], &tmp_v1[0][0]);
98 
99  switch (X->Def.iCalcModel) {
100  case HubbardGC:
101  mltplyHubbardGC(X, nstate, tmp_v0, tmp_v1);
102  break;
103 
104  case KondoGC:
105  case Hubbard:
106  case Kondo:
107  mltplyHubbard(X, nstate, tmp_v0, tmp_v1);
108  break;
109 
110  case Spin:
111  mltplySpin(X, nstate, tmp_v0, tmp_v1);
112  break;
113 
114  case SpinGC:
115  mltplySpinGC(X, nstate, tmp_v0, tmp_v1);
116  break;
117 
118  default:
119  return -1;
120  }
121 
122  StopTimer(1);
123  return 0;
124 }
129  long int n,
130  std::complex<double> a,
131  std::complex<double> *x,
132  std::complex<double> *y
133 ) {
134  long int i;
135 
136 #pragma omp parallel for default(none) private(i) shared(n, a, x, y)
137  for (i = 0; i < n; i++)
138  y[i] += a * x[i];
139 }
143 void zclear(
144  long int n,
145  std::complex<double> *x
146 ) {
147  long int i;
148 
149 #pragma omp parallel for default(none) private(i) shared(n, x)
150  for (i = 0; i < n; i++)
151  x[i] = 0.0;
152 }
int mltplyHubbard(struct BindStruct *X, int nstate, std::complex< double > **tmp_v0, std::complex< double > **tmp_v1)
perform Hamiltonian vector product for (extended) Hubbard type model.
struct DefineList Def
Definision of system (Hamiltonian) etc.
Definition: struct.hpp:395
int GetSplitBitByModel(const int Nsite, const int iCalcModel, long int *irght, long int *ilft, long int *ihfbit)
function of splitting original bit into right and left spaces.
Definition: bitcalc.cpp:78
int Nsite
Number of sites in the INTRA process region.
Definition: struct.hpp:56
void zaxpy_long(long int n, std::complex< double > a, std::complex< double > *x, std::complex< double > *y)
Wrapper of zaxpy.
Definition: mltply.cpp:128
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. . First, the calculation of diago...
Definition: mltply.cpp:56
int mltplySpinGC(struct BindStruct *X, int nstate, std::complex< double > **tmp_v0, std::complex< double > **tmp_v1)
Driver function for Spin hamiltonian.
Definition: mltplySpin.cpp:381
struct LargeList Large
Variables for Matrix-Vector product.
Definition: struct.hpp:397
int mode
multiply or expectation value.
Definition: struct.hpp:331
long int irght
Used for Ogata-Lin ???
Definition: struct.hpp:344
long int ilft
Used for Ogata-Lin ???
Definition: struct.hpp:345
Bind.
Definition: struct.hpp:394
int mltplySpin(struct BindStruct *X, int nstate, std::complex< double > **tmp_v0, std::complex< double > **tmp_v1)
Driver function for Spin hamiltonian.
Definition: mltplySpin.cpp:173
int diagonalcalcForTE(const int _istep, struct BindStruct *X, std::complex< double > *tmp_v0, std::complex< double > *tmp_v1)
void zclear(long int n, std::complex< double > *x)
clear std::complex<double> array.
Definition: mltply.cpp:143
int GetSplitBitForGeneralSpin(const int Nsite, long int *ihfbit, const long int *SiteToBit)
function of getting right, left and half bits corresponding to a original space.
Definition: bitcalc.cpp:124
long int ihfbit
Used for Ogata-Lin ???
Definition: struct.hpp:346
long int i_max
Length of eigenvector.
Definition: struct.hpp:318
int step_i
Definition: global.cpp:44
int iFlgGeneralSpin
Flag for the general (Sz/=1/2) spin.
Definition: struct.hpp:86
int mltplyHubbardGC(struct BindStruct *X, int nstate, std::complex< double > **tmp_v0, std::complex< double > **tmp_v1)
perform Hamiltonian vector product for (extended) Hubbard type model (Grandcanonical).
void StopTimer(int n)
function for calculating elapse time [elapse time=StartTimer-StopTimer]
Definition: time.cpp:83
double * list_Diagonal
Definition: global.cpp:24
long int * SiteToBit
[DefineList::NsiteMPI] Similar to DefineList::Tpow. For general spin.
Definition: struct.hpp:94
int iCalcModel
Switch for model. 0:Hubbard, 1:Spin, 2:Kondo, 3:HubbardGC, 4:SpinGC, 5:KondoGC, 6:HubbardNConserved.
Definition: struct.hpp:200
struct CheckList Check
Size of the Hilbert space.
Definition: struct.hpp:396
long int idim_max
The dimension of the Hilbert space of this process.
Definition: struct.hpp:305
int iCalcType
Switch for calculation type. 0:Lanczos, 1:TPQCalc, 2:FullDiag.
Definition: struct.hpp:194
void StartTimer(int n)
function for initializing elapse time [start]
Definition: time.cpp:71