21#ifndef VOTCA_XTP_EIGEN_H
22#define VOTCA_XTP_EIGEN_H
36typedef Eigen::Array<votca::Index, Eigen::Dynamic, 1>
ArrayXl;
43 bool mkl_overload =
false;
44#ifdef EIGEN_USE_MKL_ALL
47 bool mkl_found =
false;
51 if (mkl_overload && mkl_found) {
93 AxA(
const Eigen::Vector3d& a) {
94 data_.segment<3>(0) = a.x() * a;
95 data_.segment<2>(3) = a.y() * a.segment<2>(1);
96 data_[5] = a.z() * a.z();
98 inline const double&
xx()
const {
return data_[0]; }
99 inline const double&
xy()
const {
return data_[1]; }
100 inline const double&
xz()
const {
return data_[2]; }
101 inline const double&
yy()
const {
return data_[3]; }
102 inline const double&
yz()
const {
return data_[4]; }
103 inline const double&
zz()
const {
return data_[5]; }
112#pragma omp declare reduction (+:Mat_p_Energy: omp_out=omp_out+omp_in)\
113 initializer(omp_priv=Mat_p_Energy(omp_orig.rows(),omp_orig.cols()))
115#pragma omp declare reduction (+: Eigen::VectorXd: omp_out=omp_out+omp_in)\
116 initializer(omp_priv=Eigen::VectorXd::Zero(omp_orig.size()))
118#pragma omp declare reduction (+: Eigen::MatrixXd: omp_out=omp_out+omp_in)\
119 initializer(omp_priv=Eigen::MatrixXd::Zero(omp_orig.rows(),omp_orig.cols()))
121#pragma omp declare reduction (+: Eigen::Matrix3d: omp_out=omp_out+omp_in)\
122 initializer(omp_priv=Eigen::Matrix3d::Zero())
124#pragma omp declare reduction (+: Eigen::Vector3d: omp_out=omp_out+omp_in)\
125 initializer(omp_priv=Eigen::Vector3d::Zero())
131 nthreads =
Index(omp_get_max_threads());
138 return omp_in_parallel();
146 thread_id =
Index(omp_get_thread_num());
154 omp_set_num_threads(
int(threads));
const double & zz() const
const double & yz() const
const double & yy() const
const double & xy() const
AxA(const Eigen::Vector3d &a)
const double & xz() const
const double & xx() const
Eigen::Matrix< double, 6, 1 > data_
Eigen::MatrixXd & matrix()
const Eigen::MatrixXd & matrix() const
Mat_p_Energy operator+(const Mat_p_Energy &other) const
Mat_p_Energy(Index rows, Index cols)
Mat_p_Energy(double e, const Eigen::MatrixXd &mat)
Mat_p_Energy(double e, Eigen::MatrixXd &&mat)
void setMaxThreads(Index)
bool InsideActiveParallelRegion()
bool XTP_HAS_MKL_OVERLOAD()
base class for all analysis tools
Eigen::Array< votca::Index, Eigen::Dynamic, 1 > ArrayXl
Eigen::Matrix< double, 9, 9 > Matrix9d
Eigen::Matrix< double, 9, 1 > Vector9d