QMCPACK
DiracMatrix< T_FP > Class Template Reference

helper class to compute matrix inversion and the log value of determinant More...

+ Inheritance diagram for DiracMatrix< T_FP >:
+ Collaboration diagram for DiracMatrix< T_FP >:

Public Member Functions

 DiracMatrix ()
 
template<typename TMAT , typename ALLOC1 , typename ALLOC2 , typename TREAL , typename = std::enable_if_t<qmc_allocator_traits<ALLOC1>::is_host_accessible>, typename = std::enable_if_t<qmc_allocator_traits<ALLOC2>::is_host_accessible>>
std::enable_if_t< std::is_same< T_FP, TMAT >::value > invert_transpose (const Matrix< TMAT, ALLOC1 > &amat, Matrix< TMAT, ALLOC2 > &invMat, std::complex< TREAL > &LogDet)
 compute the inverse of the transpose of matrix A and its determinant value in log when T_FP and TMAT are the same More...
 
template<typename TMAT , typename ALLOC1 , typename ALLOC2 , typename TREAL , typename = std::enable_if_t<qmc_allocator_traits<ALLOC1>::is_host_accessible>, typename = std::enable_if_t<qmc_allocator_traits<ALLOC2>::is_host_accessible>>
std::enable_if_t<!std::is_same< T_FP, TMAT >::value > invert_transpose (const Matrix< TMAT, ALLOC1 > &amat, Matrix< TMAT, ALLOC2 > &invMat, std::complex< TREAL > &LogDet)
 compute the inverse of the transpose of matrix A and its determinant value in log when T_FP and TMAT are not the same and need scratch space psiM_fp More...
 

Private Types

using Real_FP = RealAlias< T_FP >
 

Private Member Functions

void reset (T_FP *invMat_ptr, const int lda)
 reset internal work space More...
 
template<typename TREAL >
void computeInvertAndLog (T_FP *invMat, const int n, const int lda, std::complex< TREAL > &LogDet)
 compute the inverse of invMat (in place) and the log value of determinant More...
 

Private Attributes

aligned_vector< T_FP > m_work
 
aligned_vector< int > m_pivot
 
int Lwork
 
Matrix< T_FP > psiM_fp
 scratch space used for mixed precision More...
 
aligned_vector< T_FP > LU_diag
 LU diagonal elements. More...
 

Detailed Description

template<typename T_FP>
class qmcplusplus::DiracMatrix< T_FP >

helper class to compute matrix inversion and the log value of determinant

Template Parameters
T_FPthe datatype used in the actual computation of matrix inversion

Definition at line 111 of file DiracMatrix.h.

Member Typedef Documentation

◆ Real_FP

using Real_FP = RealAlias<T_FP>
private

Definition at line 113 of file DiracMatrix.h.

Constructor & Destructor Documentation

◆ DiracMatrix()

DiracMatrix ( )
inline

Definition at line 175 of file DiracMatrix.h.

175 : Lwork(0) {}

Member Function Documentation

◆ computeInvertAndLog()

void computeInvertAndLog ( T_FP *  invMat,
const int  n,
const int  lda,
std::complex< TREAL > &  LogDet 
)
inlineprivate

compute the inverse of invMat (in place) and the log value of determinant

Template Parameters
TREALreal type
Parameters
ninvMat is n x n matrix
ldathe first dimension of invMat
LogDetlog determinant value of invMat before inversion

Definition at line 150 of file DiracMatrix.h.

Referenced by DiracMatrix< VALUE_FP >::invert_transpose().

151  {
152  BlasThreadingEnv knob(getNextLevelNumThreads());
153  if (Lwork < lda)
154  reset(invMat, lda);
155  int status = Xgetrf(n, n, invMat, lda, m_pivot.data());
156  if (status != 0)
157  {
158  std::ostringstream msg;
159  msg << "Xgetrf failed with error " << status << std::endl;
160  throw std::runtime_error(msg.str());
161  }
162  for (int i = 0; i < n; i++)
163  LU_diag[i] = invMat[i * lda + i];
164  computeLogDet(LU_diag.data(), n, m_pivot.data(), LogDet);
165  status = Xgetri(n, invMat, lda, m_pivot.data(), m_work.data(), Lwork);
166  if (status != 0)
167  {
168  std::ostringstream msg;
169  msg << "Xgetri failed with error " << status << std::endl;
170  throw std::runtime_error(msg.str());
171  }
172  }
void reset(T_FP *invMat_ptr, const int lda)
reset internal work space
Definition: DiracMatrix.h:123
aligned_vector< T_FP > LU_diag
LU diagonal elements.
Definition: DiracMatrix.h:120
int Xgetrf(int n, int m, float *restrict a, int lda, int *restrict piv)
wrappers around xgetrf lapack routines
Definition: DiracMatrix.h:31
int Xgetri(int n, float *restrict a, int lda, int *restrict piv, float *restrict work, int &lwork)
inversion of a float matrix after lu factorization
Definition: DiracMatrix.h:60
aligned_vector< int > m_pivot
Definition: DiracMatrix.h:115
void computeLogDet(const T *restrict diag, int n, const int *restrict pivot, std::complex< T_FP > &logdet)
Definition: DiracMatrix.h:100
aligned_vector< T_FP > m_work
Definition: DiracMatrix.h:114
int getNextLevelNumThreads()
get the number of threads at the next parallel level
Definition: OpenMP.h:35

◆ invert_transpose() [1/2]

std::enable_if_t<std::is_same<T_FP, TMAT>::value> invert_transpose ( const Matrix< TMAT, ALLOC1 > &  amat,
Matrix< TMAT, ALLOC2 > &  invMat,
std::complex< TREAL > &  LogDet 
)
inline

compute the inverse of the transpose of matrix A and its determinant value in log when T_FP and TMAT are the same

Template Parameters
TMATmatrix value type
TREALreal type

Definition at line 188 of file DiracMatrix.h.

Referenced by DiracMatrixComputeOMPTarget< VALUE_FP >::mw_invertTranspose(), qmcplusplus::TEST_CASE(), qmcplusplus::test_DiracDeterminant_delayed_update(), qmcplusplus::test_DiracDeterminant_second(), qmcplusplus::test_DiracDeterminantBatched_delayed_update(), qmcplusplus::test_DiracDeterminantBatched_second(), and qmcplusplus::test_inverse().

191  {
192  const int n = invMat.rows();
193  const int lda = invMat.cols();
194  simd::transpose(amat.data(), n, amat.cols(), invMat.data(), n, lda);
195  computeInvertAndLog(invMat.data(), n, lda, LogDet);
196  }
void transpose(const T *restrict A, size_t m, size_t lda, TO *restrict B, size_t n, size_t ldb)
transpose of A(m,n) to B(n,m)
Definition: algorithm.hpp:97
void computeInvertAndLog(T_FP *invMat, const int n, const int lda, std::complex< TREAL > &LogDet)
compute the inverse of invMat (in place) and the log value of determinant
Definition: DiracMatrix.h:150

◆ invert_transpose() [2/2]

std::enable_if_t<!std::is_same<T_FP, TMAT>::value> invert_transpose ( const Matrix< TMAT, ALLOC1 > &  amat,
Matrix< TMAT, ALLOC2 > &  invMat,
std::complex< TREAL > &  LogDet 
)
inline

compute the inverse of the transpose of matrix A and its determinant value in log when T_FP and TMAT are not the same and need scratch space psiM_fp

Template Parameters
TMATmatrix value type
TREALreal type

Definition at line 209 of file DiracMatrix.h.

212  {
213  const int n = invMat.rows();
214  const int lda = invMat.cols();
215  psiM_fp.resize(n, lda);
216  simd::transpose(amat.data(), n, amat.cols(), psiM_fp.data(), n, lda);
217  computeInvertAndLog(psiM_fp.data(), n, lda, LogDet);
218  invMat = psiM_fp;
219  }
void transpose(const T *restrict A, size_t m, size_t lda, TO *restrict B, size_t n, size_t ldb)
transpose of A(m,n) to B(n,m)
Definition: algorithm.hpp:97
void resize(size_type n, size_type m)
Resize the container.
Definition: OhmmsMatrix.h:99
void computeInvertAndLog(T_FP *invMat, const int n, const int lda, std::complex< TREAL > &LogDet)
compute the inverse of invMat (in place) and the log value of determinant
Definition: DiracMatrix.h:150
Matrix< T_FP > psiM_fp
scratch space used for mixed precision
Definition: DiracMatrix.h:118

◆ reset()

void reset ( T_FP *  invMat_ptr,
const int  lda 
)
inlineprivate

reset internal work space

Definition at line 123 of file DiracMatrix.h.

Referenced by DiracMatrix< VALUE_FP >::computeInvertAndLog().

124  {
125  m_pivot.resize(lda);
126  Lwork = -1;
127  T_FP tmp;
128  Real_FP lw;
129  int status = Xgetri(lda, invMat_ptr, lda, m_pivot.data(), &tmp, Lwork);
130  if (status != 0)
131  {
132  std::ostringstream msg;
133  msg << "Xgetri failed with error " << status << std::endl;
134  throw std::runtime_error(msg.str());
135  }
136 
137  lw = std::real(tmp);
138  Lwork = static_cast<int>(lw);
139  m_work.resize(Lwork);
140  LU_diag.resize(lda);
141  }
QMCTraits::RealType real
aligned_vector< T_FP > LU_diag
LU diagonal elements.
Definition: DiracMatrix.h:120
int Xgetri(int n, float *restrict a, int lda, int *restrict piv, float *restrict work, int &lwork)
inversion of a float matrix after lu factorization
Definition: DiracMatrix.h:60
aligned_vector< int > m_pivot
Definition: DiracMatrix.h:115
aligned_vector< T_FP > m_work
Definition: DiracMatrix.h:114
RealAlias< T_FP > Real_FP
Definition: DiracMatrix.h:113

Member Data Documentation

◆ LU_diag

aligned_vector<T_FP> LU_diag
private

LU diagonal elements.

Definition at line 120 of file DiracMatrix.h.

Referenced by DiracMatrix< VALUE_FP >::computeInvertAndLog(), and DiracMatrix< VALUE_FP >::reset().

◆ Lwork

int Lwork
private

◆ m_pivot

aligned_vector<int> m_pivot
private

◆ m_work

aligned_vector<T_FP> m_work
private

◆ psiM_fp

Matrix<T_FP> psiM_fp
private

scratch space used for mixed precision

Definition at line 118 of file DiracMatrix.h.

Referenced by DiracMatrix< VALUE_FP >::invert_transpose().


The documentation for this class was generated from the following file: