![]() |
QMCPACK
|
Interfaces to blas library. More...
Functions | |
static void | axpy (int n, double x, const double *a, double *b) |
static void | axpy (int n, double x, const double *a, int incx, double *b, int incy) |
static void | axpy (int n, const double *a, double *b) |
static void | axpy (int n, float x, const float *a, int incx, float *b, int incy) |
static void | axpy (int n, float x, const float *a, float *b) |
static void | axpy (int n, const std::complex< float > x, const std::complex< float > *a, int incx, std::complex< float > *b, int incy) |
static void | axpy (int n, const std::complex< double > x, const std::complex< double > *a, int incx, std::complex< double > *b, int incy) |
static float | norm2 (int n, const float *a, int incx=1) |
static float | norm2 (int n, const std::complex< float > *a, int incx=1) |
static double | norm2 (int n, const double *a, int incx=1) |
static double | norm2 (int n, const std::complex< double > *a, int incx=1) |
static void | scal (int n, float alpha, float *x, int incx=1) |
static void | scal (int n, std::complex< float > alpha, std::complex< float > *x, int incx=1) |
static void | scal (int n, double alpha, double *x, int incx=1) |
static void | scal (int n, std::complex< double > alpha, std::complex< double > *x, int incx=1) |
static void | scal (int n, double alpha, std::complex< double > *x, int incx=1) |
static void | scal (int n, float alpha, std::complex< float > *x, int incx=1) |
static void | gemv (int n, int m, const double *restrict amat, const double *restrict x, double *restrict y) |
static void | gemv (int n, int m, const float *restrict amat, const float *restrict x, float *restrict y) |
static void | gemv (int n, int m, const std::complex< double > *restrict amat, const std::complex< double > *restrict x, std::complex< double > *restrict y) |
static void | gemv (int n, int m, const std::complex< float > *restrict amat, const std::complex< float > *restrict x, std::complex< float > *restrict y) |
static void | gemv_trans (int n, int m, const double *restrict amat, const double *restrict x, double *restrict y) |
static void | gemv_trans (int n, int m, const float *restrict amat, const float *restrict x, float *restrict y) |
static void | gemv_trans (int n, int m, const std::complex< double > *restrict amat, const std::complex< double > *restrict x, std::complex< double > *restrict y) |
static void | gemv_trans (int n, int m, const std::complex< float > *restrict amat, const std::complex< float > *restrict x, std::complex< float > *restrict y) |
static void | gemv (char trans_in, int n, int m, double alpha, const double *restrict amat, int lda, const double *x, int incx, double beta, double *y, int incy) |
static void | gemv (char trans_in, int n, int m, float alpha, const float *restrict amat, int lda, const float *x, int incx, float beta, float *y, int incy) |
static void | gemv (char trans_in, int n, int m, const std::complex< double > &alpha, const std::complex< double > *restrict amat, int lda, const std::complex< double > *restrict x, int incx, const std::complex< double > &beta, std::complex< double > *y, int incy) |
static void | gemv (char trans_in, int n, int m, const std::complex< float > &alpha, const std::complex< float > *restrict amat, int lda, const std::complex< float > *restrict x, int incx, const std::complex< float > &beta, std::complex< float > *y, int incy) |
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) |
static void | gemm (char Atrans, char Btrans, int M, int N, int K, float alpha, const float *A, int lda, const float *restrict B, int ldb, float beta, float *restrict C, int ldc) |
static void | gemm (char Atrans, char Btrans, int M, int N, int K, std::complex< double > alpha, const std::complex< double > *A, int lda, const std::complex< double > *restrict B, int ldb, std::complex< double > beta, std::complex< double > *restrict C, int ldc) |
static void | gemm (char Atrans, char Btrans, int M, int N, int K, std::complex< float > alpha, const std::complex< float > *A, int lda, const std::complex< float > *restrict B, int ldb, std::complex< float > beta, std::complex< float > *restrict C, int ldc) |
template<typename T > | |
static T | dot (int n, const T *restrict a, const T *restrict b) |
template<typename T > | |
static std::complex< T > | dot (int n, const std::complex< T > *restrict a, const T *restrict b) |
template<typename T > | |
static std::complex< T > | dot (int n, const std::complex< T > *restrict a, const std::complex< T > *restrict b) |
template<typename T > | |
static std::complex< T > | dot (int n, const T *restrict a, const std::complex< T > *restrict b) |
template<typename T > | |
static T | dot (int n, const T *restrict a, int incx, const T *restrict b, int incy) |
template<typename T > | |
static std::complex< T > | dot (int n, const std::complex< T > *restrict a, int incx, const T *restrict b, int incy) |
template<typename T > | |
static std::complex< T > | dot (int n, const T *restrict a, int incx, const std::complex< T > *restrict b, int incy) |
template<typename T > | |
static std::complex< T > | dot (int n, const std::complex< T > *restrict a, int incx, const std::complex< T > *restrict b, int incy) |
template<typename T > | |
static void | copy (int n, const T *restrict a, T *restrict b) |
template<typename T > | |
static void | copy (T *restrict target, const T *restrict source, int n) |
copy using memcpy(target,source,size) More... | |
template<typename T > | |
static void | copy (int n, const std::complex< T > *restrict a, T *restrict b) |
template<typename T > | |
static void | copy (int n, const T *restrict a, std::complex< T > *restrict b) |
template<typename T > | |
static void | copy (int n, const T *restrict x, int incx, T *restrict y, int incy) |
static void | ger (int m, int n, double alpha, const double *x, int incx, const double *y, int incy, double *a, int lda) |
static void | ger (int m, int n, float alpha, const float *x, int incx, const float *y, int incy, float *a, int lda) |
static void | ger (int m, int n, const std::complex< double > &alpha, const std::complex< double > *x, int incx, const std::complex< double > *y, int incy, std::complex< double > *a, int lda) |
static void | ger (int m, int n, const std::complex< float > &alpha, const std::complex< float > *x, int incx, const std::complex< float > *y, int incy, std::complex< float > *a, int lda) |
Variables | |
constexpr int | INCX = 1 |
constexpr int | INCY = 1 |
constexpr char | UPLO = 'L' |
constexpr char | TRANS = 'T' |
constexpr char | NOTRANS = 'N' |
constexpr float | sone = 1.0e0 |
constexpr float | szero = 0.0e0 |
constexpr double | done = 1.0e0 |
constexpr double | dzero = 0.0e0 |
constexpr std::complex< float > | cone = 1.0e0 |
constexpr std::complex< float > | czero = 0.0e0 |
constexpr std::complex< double > | zone = 1.0e0 |
constexpr std::complex< double > | zzero = 0.0e0 |
Interfaces to blas library.
static data members to facilitate /Fortran blas interface static member functions to use blas functions
Arguments (float/double/complex<float>/complex<double>) determine which BLAS routines are actually used. Note that symv can be call in many ways.
|
inlinestatic |
Definition at line 55 of file BLAS.hpp.
References daxpy(), INCX, INCY, and qmcplusplus::n.
Referenced by RotatedSPOs::table_method_eval(), and RotatedSPOs::table_method_evalWF().
|
inlinestatic |
Definition at line 57 of file BLAS.hpp.
References daxpy(), and qmcplusplus::n.
|
inlinestatic |
|
inlinestatic |
Definition at line 64 of file BLAS.hpp.
References qmcplusplus::n, and saxpy().
|
inlinestatic |
|
inlinestatic |
Definition at line 71 of file BLAS.hpp.
References caxpy(), and qmcplusplus::n.
|
inlinestatic |
Definition at line 81 of file BLAS.hpp.
References qmcplusplus::n, and zaxpy().
|
inlinestatic |
Definition at line 381 of file BLAS.hpp.
References qmcplusplus::n.
Referenced by qmcplusplus::det_col_update(), LCAOrbitalSet::mw_evaluateVGL(), RotatedSPOs::table_method_eval(), and RotatedSPOs::table_method_evalWF().
|
inlinestatic |
copy using memcpy(target,source,size)
target | starting address of the targe |
source | starting address of the source |
number | of elements to copy |
Definition at line 392 of file BLAS.hpp.
References qmcplusplus::n.
|
inlinestatic |
|
inlinestatic |
|
inlinestatic |
|
inlinestatic |
Definition at line 304 of file BLAS.hpp.
References qmcplusplus::n.
Referenced by PWOrbitalSet::evaluate(), PWRealOrbitalSet::evaluate(), qmcplusplus::getRatioByColSubstitution(), qmcplusplus::getRatiosByRowSubstitution(), and qmcplusplus::getRatiosByRowSubstitution_dummy().
|
inlinestatic |
|
inlinestatic |
|
inlinestatic |
|
inlinestatic |
|
inlinestatic |
|
inlinestatic |
|
inlinestatic |
|
inlinestatic |
Definition at line 235 of file BLAS.hpp.
References qmcplusplus::Units::distance::A, B(), qmcplusplus::Units::charge::C, dgemm(), qmcplusplus::Units::energy::K, qmcplusplus::lda, and qmcplusplus::Units::force::N.
Referenced by LCAOrbitalSet::applyRotation(), SplineR2R< ST >::applyRotation(), SplineC2C< ST >::applyRotation(), RotatedSPOs::constructDeltaRotation(), RotatedSPOs::evaluateDerivatives(), RotatedSPOs::evaluateDerivativesWF(), RotatedSPOs::evaluateDerivRatios(), RotatedSPOs::exponentiate_antisym_matrix(), RotatedSPOs::log_antisym_matrix(), qmcplusplus::testing::makeRngSpdMatrix(), LCAOrbitalSet::mw_evaluateValueImplGEMM(), LCAOrbitalSet::mw_evaluateValueVPsImplGEMM(), LCAOrbitalSet::mw_evaluateVGLImplGEMM(), qmcplusplus::MatrixOperators::product(), qmcplusplus::MatrixOperators::product_ABt(), qmcplusplus::Product_ABt(), qmcplusplus::MatrixOperators::product_AtB(), RotatedSPOs::table_method_eval(), RotatedSPOs::table_method_evalWF(), qmcplusplus::test_gemm(), qmcplusplus::test_one_gemm(), and DelayedUpdate< T, T_FP >::updateInvMat().
|
inlinestatic |
Definition at line 252 of file BLAS.hpp.
References qmcplusplus::Units::distance::A, B(), qmcplusplus::Units::charge::C, qmcplusplus::Units::energy::K, qmcplusplus::lda, qmcplusplus::Units::force::N, and sgemm().
|
inlinestatic |
Definition at line 269 of file BLAS.hpp.
References qmcplusplus::Units::distance::A, B(), qmcplusplus::Units::charge::C, qmcplusplus::Units::energy::K, qmcplusplus::lda, qmcplusplus::Units::force::N, and zgemm().
|
inlinestatic |
Definition at line 286 of file BLAS.hpp.
References qmcplusplus::Units::distance::A, B(), qmcplusplus::Units::charge::C, cgemm(), qmcplusplus::Units::energy::K, qmcplusplus::lda, and qmcplusplus::Units::force::N.
|
inlinestatic |
Definition at line 118 of file BLAS.hpp.
References dgemv(), done, dzero, INCX, INCY, qmcplusplus::Units::distance::m, qmcplusplus::n, and NOTRANS.
Referenced by DelayedUpdate< T, T_FP >::acceptRow(), DelayedUpdateSYCL< T, T_FP >::acceptRow(), DelayedUpdateCUDA< T, T_FP >::acceptRow(), qmcplusplus::det_col_update(), qmcplusplus::det_row_update(), SOECPComponent::evaluateValueAndDerivatives(), NonLocalECPComponent::evaluateValueAndDerivatives(), DelayedUpdate< T, T_FP >::getInvRow(), DelayedUpdateSYCL< T, T_FP >::getInvRow(), DelayedUpdateCUDA< T, T_FP >::getInvRow(), qmcplusplus::getRatiosByRowSubstitution(), qmcplusplus::MatrixOperators::product(), qmcplusplus::MatrixOperators::product_Atx(), qmcplusplus::test_gemv(), qmcplusplus::test_gemv_batched(), qmcplusplus::test_one_gemv(), and DelayedUpdate< T, T_FP >::updateInvMat().
|
inlinestatic |
Definition at line 123 of file BLAS.hpp.
References done, dzero, INCX, INCY, qmcplusplus::Units::distance::m, qmcplusplus::n, NOTRANS, and sgemv().
|
inlinestatic |
Definition at line 128 of file BLAS.hpp.
References INCX, INCY, qmcplusplus::Units::distance::m, qmcplusplus::n, NOTRANS, zgemv(), zone, and zzero.
|
inlinestatic |
Definition at line 137 of file BLAS.hpp.
References cgemv(), cone, czero, INCX, INCY, qmcplusplus::Units::distance::m, qmcplusplus::n, and NOTRANS.
|
inlinestatic |
Definition at line 175 of file BLAS.hpp.
References dgemv(), qmcplusplus::lda, qmcplusplus::Units::distance::m, and qmcplusplus::n.
|
inlinestatic |
Definition at line 190 of file BLAS.hpp.
References qmcplusplus::lda, qmcplusplus::Units::distance::m, qmcplusplus::n, and sgemv().
|
inlinestatic |
Definition at line 205 of file BLAS.hpp.
References qmcplusplus::lda, qmcplusplus::Units::distance::m, qmcplusplus::n, and zgemv().
|
inlinestatic |
Definition at line 220 of file BLAS.hpp.
References cgemv(), qmcplusplus::lda, qmcplusplus::Units::distance::m, and qmcplusplus::n.
|
inlinestatic |
Definition at line 147 of file BLAS.hpp.
References dgemv(), done, dzero, INCX, INCY, qmcplusplus::Units::distance::m, qmcplusplus::n, and TRANS.
Referenced by qmcplusplus::test_gemv(), qmcplusplus::test_gemv_batched(), and qmcplusplus::test_one_gemv().
|
inlinestatic |
Definition at line 152 of file BLAS.hpp.
References done, dzero, INCX, INCY, qmcplusplus::Units::distance::m, qmcplusplus::n, sgemv(), and TRANS.
|
inlinestatic |
Definition at line 157 of file BLAS.hpp.
References done, dzero, INCX, INCY, qmcplusplus::Units::distance::m, qmcplusplus::n, TRANS, and zgemv().
|
inlinestatic |
Definition at line 166 of file BLAS.hpp.
References cgemv(), done, dzero, INCX, INCY, qmcplusplus::Units::distance::m, qmcplusplus::n, and TRANS.
|
inlinestatic |
Definition at line 437 of file BLAS.hpp.
References dger(), qmcplusplus::lda, qmcplusplus::Units::distance::m, and qmcplusplus::n.
Referenced by DelayedUpdate< T, T_FP >::acceptRow(), DelayedUpdateSYCL< T, T_FP >::acceptRow(), DelayedUpdateCUDA< T, T_FP >::acceptRow(), qmcplusplus::det_col_update(), qmcplusplus::det_row_update(), qmcplusplus::test_ger(), qmcplusplus::test_ger_batched(), qmcplusplus::test_one_ger(), and DelayedUpdate< T, T_FP >::updateInvMat().
|
inlinestatic |
Definition at line 450 of file BLAS.hpp.
References qmcplusplus::lda, qmcplusplus::Units::distance::m, qmcplusplus::n, and sger().
|
inlinestatic |
Definition at line 463 of file BLAS.hpp.
References qmcplusplus::lda, qmcplusplus::Units::distance::m, qmcplusplus::n, and zgeru().
|
inlinestatic |
Definition at line 476 of file BLAS.hpp.
References cgeru(), qmcplusplus::lda, qmcplusplus::Units::distance::m, and qmcplusplus::n.
|
inlinestatic |
Definition at line 91 of file BLAS.hpp.
References qmcplusplus::n, and snrm2().
|
inlinestatic |
Definition at line 93 of file BLAS.hpp.
References qmcplusplus::n, and scnrm2().
|
inlinestatic |
Definition at line 95 of file BLAS.hpp.
References dnrm2(), and qmcplusplus::n.
|
inlinestatic |
Definition at line 97 of file BLAS.hpp.
References dznrm2(), and qmcplusplus::n.
|
inlinestatic |
Definition at line 99 of file BLAS.hpp.
References qmcplusplus::n, and sscal().
|
inlinestatic |
Definition at line 101 of file BLAS.hpp.
References cscal(), and qmcplusplus::n.
|
inlinestatic |
Definition at line 106 of file BLAS.hpp.
References dscal(), and qmcplusplus::n.
|
inlinestatic |
Definition at line 108 of file BLAS.hpp.
References qmcplusplus::n, and zscal().
|
inlinestatic |
Definition at line 113 of file BLAS.hpp.
References qmcplusplus::n, and zdscal().
|
inlinestatic |
Definition at line 115 of file BLAS.hpp.
References csscal(), and qmcplusplus::n.
constexpr std::complex<float> cone = 1.0e0 |
Definition at line 50 of file BLAS.hpp.
Referenced by DelayedUpdate< T, T_FP >::acceptRow(), DelayedUpdateSYCL< T, T_FP >::acceptRow(), DelayedUpdateCUDA< T, T_FP >::acceptRow(), qmcplusplus::bessel_steed_array_cpu(), NonLocalECPComponent::calculateProjector(), JeeIOrbitalSoA< FT >::computeU3_engine(), qmcplusplus::det_col_update(), MPC::evalSR(), PolynomialFunctor3D::evaluate(), SoaSphericalTensor< ST >::evaluate_bare(), AtomicOrbitals< ST >::evaluate_vgl(), PolynomialFunctor3D::evaluateDerivatives(), JeeIOrbitalSoA< FT >::evaluateDerivatives(), J1Spin< FT >::evaluateDerivativesWF(), TwoBodyJastrow< FT >::evaluateDerivativesWF(), J1OrbitalSoA< FT >::evaluateDerivativesWF(), JeeIOrbitalSoA< FT >::evaluateDerivativesWF(), NonLocalECPComponent::evaluateOneBodyOpMatrixContribution(), NonLocalECPComponent::evaluateOneBodyOpMatrixdRContribution(), NonLocalECPComponent::evaluateOneWithForces(), PolynomialFunctor3D::evaluateV(), SoaAtomicBasisSet< ROT, SH >::evaluateVGH(), SoaAtomicBasisSet< ROT, SH >::evaluateVGHGH(), SoaAtomicBasisSet< ROT, SH >::evaluateVGL(), PolynomialFunctor3D::evaluateVGL(), gemv(), DelayedUpdate< T, T_FP >::getInvRow(), DelayedUpdateSYCL< T, T_FP >::getInvRow(), DelayedUpdateCUDA< T, T_FP >::getInvRow(), HybridRepCenterOrbitals< SPLINEBASE::DataType >::interpolate_buffer_v(), HybridRepCenterOrbitals< SPLINEBASE::DataType >::interpolate_buffer_vgl(), SoaAtomicBasisSet< ROT, SH >::mw_evaluateVGL(), MultiDiracDeterminant::mw_InverseUpdateByColumn(), HybridRepCenterOrbitals< SPLINEBASE::DataType >::selectRegionAndComputeSmoothing(), QMCUpdateBase::setMultiplicity(), qmcplusplus::smoothing(), SoaSphericalTensor< ST >::SoaSphericalTensor(), DelayedUpdate< T, T_FP >::updateInvMat(), DelayedUpdateSYCL< T, T_FP >::updateInvMat(), and DelayedUpdateBatched< PL, VALUE >::updateRow().
constexpr std::complex<float> czero = 0.0e0 |
Definition at line 51 of file BLAS.hpp.
Referenced by DelayedUpdate< T, T_FP >::acceptRow(), DelayedUpdateSYCL< T, T_FP >::acceptRow(), DelayedUpdateCUDA< T, T_FP >::acceptRow(), SoaCartesianTensor< T >::batched_evaluateVGL(), NonLocalECPComponent::calculateProjector(), TwoBodyJastrow< FT >::computeU3(), J1OrbitalSoA< FT >::computeU3(), J1Spin< FT >::computeU3(), JeeIOrbitalSoA< FT >::computeU3(), JeeIOrbitalSoA< FT >::computeU3_engine(), CoulombPBCAB::evalSR(), CoulombPBCAB::evalSRwithForces(), PolynomialFunctor3D::evaluate(), SoaSphericalTensor< ST >::evaluate_bare(), AtomicOrbitals< ST >::evaluate_v(), AtomicOrbitals< ST >::evaluate_vgl(), MultiSlaterDetTableMethod::evaluate_vgl_impl(), CoulombPotential< T >::evaluateAB(), PolynomialFunctor3D::evaluateDerivatives(), JeeIOrbitalSoA< FT >::evaluateDerivatives(), JeeIOrbitalSoA< FT >::evaluateDerivativesWF(), JeeIOrbitalSoA< FT >::evaluateDerivRatios(), NonLocalECPComponent::evaluateOneBodyOpMatrixContribution(), NonLocalECPComponent::evaluateOneBodyOpMatrixdRContribution(), NonLocalECPComponent::evaluateOneWithForces(), PolynomialFunctor3D::evaluateV(), AtomicOrbitals< ST >::evaluateValues(), SoaCartesianTensor< T >::evaluateVGH(), SoaCartesianTensor< T >::evaluateVGHGH(), SoaCartesianTensor< T >::evaluateVGL(), PolynomialFunctor3D::evaluateVGL(), SoaSphericalTensor< ST >::evaluateVGL_impl(), gemv(), DelayedUpdate< T, T_FP >::getInvRow(), DelayedUpdateSYCL< T, T_FP >::getInvRow(), DelayedUpdateCUDA< T, T_FP >::getInvRow(), SlaterDet::mw_accept_rejectMove(), TrialWaveFunction::mw_evaluateDeltaLog(), TrialWaveFunction::mw_evaluateDeltaLogSetup(), SlaterDet::mw_evaluateGL(), TrialWaveFunction::mw_evaluateGL(), SlaterDet::mw_evaluateLog(), TrialWaveFunction::mw_evaluateLog(), SoaSphericalTensor< ST >::SoaSphericalTensor(), DelayedUpdate< T, T_FP >::updateInvMat(), DelayedUpdateSYCL< T, T_FP >::updateInvMat(), and DelayedUpdateBatched< PL, VALUE >::updateRow().
constexpr double done = 1.0e0 |
Definition at line 48 of file BLAS.hpp.
Referenced by BackflowBuilder::addRPA(), BackflowBuilder::addTwoBody(), axpy(), NonlinearFitClass< M, ModelType >::Fit(), gemv(), gemv_trans(), QPParser::getGeometry(), GamesAsciiParser::getGeometry(), Backflow_eI< FT >::makeClone(), Backflow_ee< FT >::makeClone(), GaussianFCHKParser::parse(), SimpleGrid::ReverseMap(), and GeneralGrid::ReverseMap().
constexpr double dzero = 0.0e0 |
Definition at line 49 of file BLAS.hpp.
Referenced by gemv(), and gemv_trans().
constexpr int INCX = 1 |
Definition at line 40 of file BLAS.hpp.
Referenced by axpy(), gemv(), and gemv_trans().
constexpr int INCY = 1 |
Definition at line 41 of file BLAS.hpp.
Referenced by axpy(), gemv(), and gemv_trans().
constexpr char TRANS = 'T' |
Definition at line 43 of file BLAS.hpp.
Referenced by gemv_trans().
constexpr char UPLO = 'L' |
Definition at line 42 of file BLAS.hpp.
Referenced by RotatedSPOs::exponentiate_antisym_matrix(), LAPACK::hevr(), and LAPACK::potrf().
constexpr std::complex<double> zone = 1.0e0 |
Definition at line 52 of file BLAS.hpp.
Referenced by gemv(), qmcplusplus::MatrixOperators::product(), and qmcplusplus::Product_ABt().