QMCPACK
qmcplusplus::MatrixOperators Namespace Reference

Functions

template<typename T >
void product (const Matrix< T > &A, const Matrix< T > &B, Matrix< T > &C)
 static function to perform C=AB for real matrices More...
 
template<typename T >
void product_ABt (const Matrix< T > &A, const Matrix< T > &B, Matrix< T > &C)
 
template<typename T >
void product_AtB (const Matrix< T > &A, const Matrix< T > &B, Matrix< T > &C)
 
void half_outerProduct (const Matrix< double > &M, const Vector< double > &B, int iat, Matrix< double > &C)
 
void other_half_outerProduct (const Matrix< double > &M, const Vector< double > &B, int iat, Matrix< double > &C)
 
double innerProduct (const Vector< double > &A, const Vector< double > &B)
 
template<typename T >
void transpose (Matrix< T > &A)
 
template<typename T >
void transpose (const Matrix< T > &A, Matrix< T > &B)
 
template<typename T1 , typename T2 , typename T3 >
void diag_product (const Matrix< T1 > &A, const Vector< T2 > &B, Matrix< T3 > &C)
 C = A*diag(B) More...
 
template<typename T1 , typename T2 , typename T3 >
void diag_product (const Vector< T1 > &A, const Matrix< T2 > &B, Matrix< T3 > &C)
 C = diag(A)*B. More...
 
void product (const Matrix< double > &A, const Matrix< std::complex< double >> &B, Matrix< double > &C)
 static function to perform C=AB for complex matrices More...
 
template<typename T >
void product (const Matrix< T > &A, const Vector< T > &x, Vector< T > &y)
 static function to perform y=Ax for generic matrix and vector More...
 
template<typename T >
void product_Atx (const Matrix< T > &A, const Vector< T > &x, Vector< T > &y)
 static function to perform y = A^t x for generic matrix and vector More...
 
template<typename T , unsigned D>
void product (const Matrix< T > &A, const TinyVector< T, D > *xvPtr, TinyVector< T, D > *restrict yptr)
 static function to perform y=Ax for generic matrix and vector More...
 
template<typename T , unsigned D>
void product (const Matrix< T > &A, const Tensor< T, D > *xvPtr, Tensor< T, D > *restrict yptr)
 
template<unsigned D>
void product (const Matrix< double > &A, const Vector< TinyVector< double, D >> &x, TinyVector< double, D > *restrict yptr)
 static function to perform y=Ax for generic matrix and vector More...
 
template<unsigned D>
void product (const Matrix< std::complex< double >> &A, const Vector< TinyVector< std::complex< double >, D >> &x, TinyVector< std::complex< double >, D > *restrict yptr)
 static function to perform y=Ax for generic matrix and vector More...
 
void product (const Matrix< std::complex< double >> &A, const Vector< std::complex< double >> &x, std::complex< double > *restrict yptr)
 static function to perform y=Ax for generic matrix and vector More...
 
void product (const Matrix< std::complex< double >> &A, const Vector< double > &x, std::complex< double > *restrict yptr)
 static function to perform y=Ax for generic matrix and vector More...
 
template<typename MAT1 , typename MAT2 >
void insert_columns (const MAT1 &small, MAT2 &big, int offset_c)
 copy a small matrix (N, M1) to a big matrix (N, M2), M2>M1 More...
 

Function Documentation

◆ diag_product() [1/2]

void qmcplusplus::MatrixOperators::diag_product ( const Matrix< T1 > &  A,
const Vector< T2 > &  B,
Matrix< T3 > &  C 
)
inline

C = A*diag(B)

Definition at line 114 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A, B(), and qmcplusplus::Units::charge::C.

Referenced by DensityMatrices1B::evaluate_matrix(), and OneBodyDensityMatrices::evaluateMatrix().

115 {
116  for (int i = 0; i < C.rows(); ++i)
117  for (int j = 0; j < C.cols(); ++j)
118  C(i, j) = A(i, j) * B[j];
119  //works?
120  //const int ccols = C.cols();
121  //const int ijmax = C.size();
122  //for (int ij=0; ij<ijmax; ++ij)
123  // C(ij)=A(ij)*B(ij%ccols);
124 }
double B(double x, int k, int i, const std::vector< double > &t)

◆ diag_product() [2/2]

void qmcplusplus::MatrixOperators::diag_product ( const Vector< T1 > &  A,
const Matrix< T2 > &  B,
Matrix< T3 > &  C 
)
inline

C = diag(A)*B.

Definition at line 129 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A, B(), and qmcplusplus::Units::charge::C.

130 {
131  for (int i = 0; i < C.rows(); ++i)
132  for (int j = 0; j < C.cols(); ++j)
133  C(i, j) = A[i] * B(i, j);
134 
135 
136  //const int crows = C.rows();
137  //const int ccols = C.cols();
138  //for (int i=0,ijb=0; i<crows; ++i,ijb+=ccols)
139  //{
140  // const T1 a = A(i);
141  // for (int j=0; j<ccols; ++j)
142  // {
143  // int ij = ijb+j;
144  // C(ij)=a*B(ij);
145  // }
146  //}
147  //
148  //const int crows = C.rows();
149  //const int ijmax = C.size();
150  //for (int ij=0; ij<ijmax; ++ij)
151  // C(ij)=A(ij%crows)*B(ij);
152 }
double B(double x, int k, int i, const std::vector< double > &t)

◆ half_outerProduct()

void qmcplusplus::MatrixOperators::half_outerProduct ( const Matrix< double > &  M,
const Vector< double > &  B,
int  iat,
Matrix< double > &  C 
)
inline

Definition at line 69 of file MatrixOperators.h.

References B(), and qmcplusplus::Units::charge::C.

70 {
71  for (int i = 0; i < C.rows(); i++)
72  for (int j = 0; j < C.cols(); j++)
73  C(iat, i) += M(i, j) * B[j];
74  //for (int i=0; i<C.rows(); i++)
75  // C(iat,i)+=simd::dot(M[i],B.data(),C.cols());
76 }
double B(double x, int k, int i, const std::vector< double > &t)

◆ innerProduct()

double qmcplusplus::MatrixOperators::innerProduct ( const Vector< double > &  A,
const Vector< double > &  B 
)
inline

Definition at line 87 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A, and B().

88 {
89  double tot = 0.0;
90  for (int i = 0; i < A.size(); i++)
91  tot += A[i] * B[i];
92  return tot;
93  //return simd::dot(A.data(),B.data(),A.size());
94 }
double B(double x, int k, int i, const std::vector< double > &t)

◆ insert_columns()

void qmcplusplus::MatrixOperators::insert_columns ( const MAT1 &  small,
MAT2 &  big,
int  offset_c 
)
inline

copy a small matrix (N, M1) to a big matrix (N, M2), M2>M1

Parameters
smallinput matrix
bigoutout matrix
offset_ccolumn offset
Todo:
smater and more efficient matrix, move up for others The columns [0,M1) are inserted into [offset_c,offset_c+M1).

Definition at line 34 of file CompositeSPOSet.cpp.

References copy().

Referenced by CompositeSPOSet::evaluate_notranspose().

35 {
36  const int c = small.cols();
37  for (int i = 0; i < small.rows(); ++i)
38  std::copy(small[i], small[i] + c, big[i] + offset_c);
39 }
void copy(const Array< T1, 3 > &src, Array< T2, 3 > &dest)
Definition: Blitz.h:639

◆ other_half_outerProduct()

void qmcplusplus::MatrixOperators::other_half_outerProduct ( const Matrix< double > &  M,
const Vector< double > &  B,
int  iat,
Matrix< double > &  C 
)
inline

Definition at line 78 of file MatrixOperators.h.

References B(), and qmcplusplus::Units::charge::C.

79 {
80  for (int i = 0; i < C.rows(); i++)
81  for (int j = 0; j < C.cols(); j++)
82  C(i, iat) += M(i, j) * B[j];
83  //for (int i=0; i<C.rows(); i++)
84  // C(i,iat)+=simd::dot(M[i],B.data(),C.cols());
85 }
double B(double x, int k, int i, const std::vector< double > &t)

◆ product() [1/9]

void qmcplusplus::MatrixOperators::product ( const Matrix< T > &  A,
const Matrix< T > &  B,
Matrix< T > &  C 
)
inline

static function to perform C=AB for real matrices

Call dgemm/zgemm/sgemm/cgemm via BLAS::gemm

Definition at line 36 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A, B(), qmcplusplus::Units::charge::C, and BLAS::gemm().

Referenced by DensityMatrices1B::evaluate_matrix(), PWOrbitalSet::evaluate_notranspose(), PWRealOrbitalSet::evaluate_notranspose(), AGPDeterminant::evaluateLogAndStore(), OneBodyDensityMatrices::evaluateMatrix(), DiracDeterminant< DU_TYPE >::evaluateRatiosAlltoOne(), PWOrbitalSet::evaluateValue(), PWRealOrbitalSet::evaluateValue(), LCAOrbitalSet::evaluateValue(), PWOrbitalSet::evaluateVGL(), PWRealOrbitalSet::evaluateVGL(), ForceCeperley::InitMatrix(), ForceChiesaPBCAA::InitMatrix(), MultiSlaterDetTableMethod::precomputeC_otherDs(), QMCFixedSampleLinearOptimizeBatched::previous_linear_methods_run(), QMCFixedSampleLinearOptimize::run(), Eigensolver::solveGeneralizedEigenvalues_Inv(), and QMCFixedSampleLinearOptimizeBatched::solveShiftsWithoutLMYEngine().

37 {
38  const char transa = 'N';
39  const char transb = 'N';
40  const T one(1.0);
41  const T zero(0.0);
42  BLAS::gemm(transa, transb, B.cols(), A.rows(), B.rows(), one, B.data(), B.cols(), A.data(), A.cols(), zero, C.data(),
43  C.cols());
44 }
static void gemm(char Atrans, char Btrans, int M, int N, int K, double alpha, const double *A, int lda, const double *restrict B, int ldb, double beta, double *restrict C, int ldc)
Definition: BLAS.hpp:235
double B(double x, int k, int i, const std::vector< double > &t)

◆ product() [2/9]

void qmcplusplus::MatrixOperators::product ( const Matrix< double > &  A,
const Matrix< std::complex< double >> &  B,
Matrix< double > &  C 
)
inline

static function to perform C=AB for complex matrices

Call zgemm

Definition at line 159 of file MatrixOperators.h.

160 {
161  std::cerr << " Undefined C=AB with real A and complex B " << std::endl;
162 }

◆ product() [3/9]

void qmcplusplus::MatrixOperators::product ( const Matrix< T > &  A,
const Vector< T > &  x,
Vector< T > &  y 
)
inline

static function to perform y=Ax for generic matrix and vector

Definition at line 167 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A, Vector< T, Alloc >::data(), and BLAS::gemv().

168 {
169  const char transa = 'T';
170  const T one = 1.0;
171  const T zero = 0.0;
172  BLAS::gemv(transa, A.cols(), A.rows(), one, A.data(), A.cols(), x.data(), 1, zero, y.data(), 1);
173 }
static void gemv(int n, int m, const double *restrict amat, const double *restrict x, double *restrict y)
Definition: BLAS.hpp:118

◆ product() [4/9]

void qmcplusplus::MatrixOperators::product ( const Matrix< T > &  A,
const TinyVector< T, D > *  xvPtr,
TinyVector< T, D > *restrict  yptr 
)
inline

static function to perform y=Ax for generic matrix and vector

Definition at line 189 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A, TinyVector< T, D >::begin(), and BLAS::gemm().

190 {
191  const T one = 1.0;
192  const T zero = 0.0;
193  const char transa = 'N';
194  const char transb = 'N';
195  BLAS::gemm(transa, transb, D, A.rows(), A.cols(), one, xvPtr->begin(), D, A.data(), A.cols(), zero, yptr->begin(), D);
196 }
static void gemm(char Atrans, char Btrans, int M, int N, int K, double alpha, const double *A, int lda, const double *restrict B, int ldb, double beta, double *restrict C, int ldc)
Definition: BLAS.hpp:235

◆ product() [5/9]

void qmcplusplus::MatrixOperators::product ( const Matrix< T > &  A,
const Tensor< T, D > *  xvPtr,
Tensor< T, D > *restrict  yptr 
)
inline

Definition at line 199 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A, Tensor< T, D >::begin(), and BLAS::gemm().

200 {
201  const T one = 1.0;
202  const T zero = 0.0;
203  const char transa = 'N';
204  const char transb = 'N';
205  BLAS::gemm(transa, transb, D * D, A.rows(), A.cols(), one, xvPtr->begin(), D * D, A.data(), A.cols(), zero,
206  yptr->begin(), D * D);
207 }
static void gemm(char Atrans, char Btrans, int M, int N, int K, double alpha, const double *A, int lda, const double *restrict B, int ldb, double beta, double *restrict C, int ldc)
Definition: BLAS.hpp:235

◆ product() [6/9]

void qmcplusplus::MatrixOperators::product ( const Matrix< double > &  A,
const Vector< TinyVector< double, D >> &  x,
TinyVector< double, D > *restrict  yptr 
)
inline

static function to perform y=Ax for generic matrix and vector

Definition at line 212 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A, TinyVector< T, D >::begin(), and dgemm().

215 {
216  const double one = 1.0;
217  const double zero = 0.0;
218  const char transa = 'N';
219  const char transb = 'N';
220  dgemm(transa, transb, D, A.rows(), x.size(), one, x.data()->begin(), D, A.data(), A.cols(), zero, yptr->begin(), D);
221 }
void dgemm(const char &, const char &, const int &, const int &, const int &, const double &, const double *, const int &, const double *, const int &, const double &, double *, const int &)

◆ product() [7/9]

void qmcplusplus::MatrixOperators::product ( const Matrix< std::complex< double >> &  A,
const Vector< TinyVector< std::complex< double >, D >> &  x,
TinyVector< std::complex< double >, D > *restrict  yptr 
)
inline

static function to perform y=Ax for generic matrix and vector

Definition at line 226 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A, zgemm(), and BLAS::zone.

229 {
230  const char transa = 'N';
231  const char transb = 'N';
232  const std::complex<double> zone(1.0, 0.0);
233  const std::complex<double> zero(0.0, 0.0);
234  zgemm(transa, transb, D, A.rows(), x.size(), zone, x.data()->begin(), D, A.data(), A.cols(), zero, yptr->begin(), D);
235 }
void zgemm(const char &, const char &, const int &, const int &, const int &, const std::complex< double > &, const std::complex< double > *, const int &, const std::complex< double > *, const int &, const std::complex< double > &, std::complex< double > *, const int &)
constexpr std::complex< double > zone
Definition: BLAS.hpp:52

◆ product() [8/9]

void qmcplusplus::MatrixOperators::product ( const Matrix< std::complex< double >> &  A,
const Vector< std::complex< double >> &  x,
std::complex< double > *restrict  yptr 
)
inline

static function to perform y=Ax for generic matrix and vector

Definition at line 240 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A, zgemv(), and BLAS::zone.

243 {
244  const char transa = 'T';
245  const std::complex<double> zone(1.0, 0.0);
246  const std::complex<double> zero(0.0, 0.0);
247  zgemv(transa, A.cols(), A.rows(), zone, A.data(), A.cols(), x.data(), 1, zero, yptr, 1);
248 }
constexpr std::complex< double > zone
Definition: BLAS.hpp:52
void zgemv(const char &trans, const int &nr, const int &nc, const std::complex< double > &alpha, const std::complex< double > *amat, const int &lda, const std::complex< double > *bv, const int &incx, const std::complex< double > &beta, std::complex< double > *cv, const int &incy)

◆ product() [9/9]

void qmcplusplus::MatrixOperators::product ( const Matrix< std::complex< double >> &  A,
const Vector< double > &  x,
std::complex< double > *restrict  yptr 
)
inline

static function to perform y=Ax for generic matrix and vector

Definition at line 252 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A, qmcplusplus::Units::distance::m, and qmcplusplus::n.

253 {
254  const int n = A.rows();
255  const int m = A.cols();
256  const std::complex<double>* restrict aptr = A.data();
257  for (int i = 0, ij = 0; i < n; ++i)
258  {
259  std::complex<double> t = 0.0;
260  for (int j = 0; j < m; ++j, ++ij)
261  t += aptr[ij] * x[j];
262  yptr[i] = t;
263  }
264 }

◆ product_ABt()

void qmcplusplus::MatrixOperators::product_ABt ( const Matrix< T > &  A,
const Matrix< T > &  B,
Matrix< T > &  C 
)
inline

Definition at line 47 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A, B(), qmcplusplus::Units::charge::C, and BLAS::gemm().

48 {
49  const char transa = 't';
50  const char transb = 'n';
51  const T one(1.0);
52  const T zero(0.0);
53  BLAS::gemm(transa, transb, B.rows(), A.rows(), B.cols(), one, B.data(), B.cols(), A.data(), A.cols(), zero, C.data(),
54  C.cols());
55 }
static void gemm(char Atrans, char Btrans, int M, int N, int K, double alpha, const double *A, int lda, const double *restrict B, int ldb, double beta, double *restrict C, int ldc)
Definition: BLAS.hpp:235
double B(double x, int k, int i, const std::vector< double > &t)

◆ product_AtB()

void qmcplusplus::MatrixOperators::product_AtB ( const Matrix< T > &  A,
const Matrix< T > &  B,
Matrix< T > &  C 
)
inline

Definition at line 58 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A, B(), qmcplusplus::Units::charge::C, and BLAS::gemm().

Referenced by DensityMatrices1B::evaluate_matrix(), and OneBodyDensityMatrices::evaluateMatrix().

59 {
60  const char transa = 'n';
61  const char transb = 't';
62  const T one(1.0);
63  const T zero(0.0);
64  BLAS::gemm(transa, transb, B.cols(), A.cols(), B.rows(), one, B.data(), B.cols(), A.data(), A.cols(), zero, C.data(),
65  C.cols());
66 }
static void gemm(char Atrans, char Btrans, int M, int N, int K, double alpha, const double *A, int lda, const double *restrict B, int ldb, double beta, double *restrict C, int ldc)
Definition: BLAS.hpp:235
double B(double x, int k, int i, const std::vector< double > &t)

◆ product_Atx()

void qmcplusplus::MatrixOperators::product_Atx ( const Matrix< T > &  A,
const Vector< T > &  x,
Vector< T > &  y 
)
inline

static function to perform y = A^t x for generic matrix and vector

Definition at line 178 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A, Vector< T, Alloc >::data(), and BLAS::gemv().

Referenced by LCAOrbitalSet::evaluateDetRatios().

179 {
180  const char transa = 'N';
181  const T one = 1.0;
182  const T zero = 0.0;
183  BLAS::gemv(transa, A.cols(), A.rows(), one, A.data(), A.cols(), x.data(), 1, zero, y.data(), 1);
184 }
static void gemv(int n, int m, const double *restrict amat, const double *restrict x, double *restrict y)
Definition: BLAS.hpp:118

◆ transpose() [1/2]

void qmcplusplus::MatrixOperators::transpose ( Matrix< T > &  A)
inline

Definition at line 98 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A.

99 {
100  for (int i = 0; i < A.extent(0); i++)
101  for (int j = 0; j < i; j++)
102  std::swap(A(i, j), A(j, i));
103 }

◆ transpose() [2/2]

void qmcplusplus::MatrixOperators::transpose ( const Matrix< T > &  A,
Matrix< T > &  B 
)
inline

Definition at line 106 of file MatrixOperators.h.

References qmcplusplus::Units::distance::A, B(), and qmcplusplus::simd::transpose().

107 {
108  simd::transpose(A.data(), A.rows(), A.cols(), B.data(), B.rows(), B.cols());
109 }
void transpose(const Matrix< T > &A, Matrix< T > &B)
double B(double x, int k, int i, const std::vector< double > &t)