Skip to content
This repository has been archived by the owner on Sep 28, 2021. It is now read-only.

using eigen middleCol and middleRow functions #732

Merged
merged 6 commits into from
Aug 18, 2021
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 1 addition & 0 deletions CHANGELOG.rst
Original file line number Diff line number Diff line change
Expand Up @@ -48,6 +48,7 @@ Version 2022-dev
- Fixed exciton options checking (#726)
- added basis gpu runner and test to suite (#725)
- turned sigma choice into a factory (#731)
- use middleRows/Cols and refactor numerical integration (#732)

Version 2021.2 (released XX.07.21)
==================================
Expand Down
16 changes: 11 additions & 5 deletions include/votca/xtp/aoshell.h
Original file line number Diff line number Diff line change
Expand Up @@ -121,11 +121,17 @@ class AOShell {

double getMinDecay() const { return mindecay_; }

void EvalAOspace(Eigen::VectorBlock<Eigen::VectorXd>& AOvalues,
const Eigen::Vector3d& grid_pos) const;
void EvalAOspace(Eigen::VectorBlock<Eigen::VectorXd>& AOvalues,
Eigen::Block<Eigen::MatrixX3d>& AODervalues,
const Eigen::Vector3d& grid_pos) const;
struct AOValues {

AOValues(Index size) {
values = Eigen::VectorXd::Zero(size);
derivatives = Eigen::MatrixX3d::Zero(size, 3);
}
Eigen::VectorXd values;
Eigen::MatrixX3d derivatives;
};

AOValues EvalAOspace(const Eigen::Vector3d& grid_pos) const;

// iterator over pairs (decay constant; contraction coefficient)
using GaussianIterator = std::vector<AOGaussianPrimitive>::const_iterator;
Expand Down
15 changes: 15 additions & 0 deletions include/votca/xtp/cudamatrix.h
Original file line number Diff line number Diff line change
Expand Up @@ -125,6 +125,21 @@ class CudaMatrix {
return CudaMatrixBlock<CudaMatrix>(*this, rowoffset, coloffset, rows, cols);
}

CudaMatrixBlock<CudaMatrix> row(Index row) const {
return CudaMatrixBlock<CudaMatrix>(*this, row, 0, 1, cols());
}

CudaMatrixBlock<CudaMatrix> col(Index col) const {
return CudaMatrixBlock<CudaMatrix>(*this, 0, col, rows(), 1);
}

CudaMatrixBlock<CudaMatrix> middleRows(Index rowoffset, Index rows) const {
return CudaMatrixBlock<CudaMatrix>(*this, rowoffset, 0, rows, cols());
}
CudaMatrixBlock<CudaMatrix> middleCols(Index coloffset, Index cols) const {
return CudaMatrixBlock<CudaMatrix>(*this, 0, coloffset, rows(), cols);
}

static constexpr bool transposed() { return false; }

template <class T>
Expand Down
4 changes: 1 addition & 3 deletions include/votca/xtp/gridbox.h
Original file line number Diff line number Diff line change
Expand Up @@ -36,9 +36,7 @@ class GridBox {

public:
void FindSignificantShells(const AOBasis& basis);
Eigen::VectorXd CalcAOValue_and_Grad(Eigen::MatrixX3d& ao_grad,
const Eigen::Vector3d& point) const;
Eigen::VectorXd CalcAOValues(const Eigen::Vector3d& pos) const;
AOShell::AOValues CalcAOValues(const Eigen::Vector3d& point) const;

const std::vector<Eigen::Vector3d>& getGridPoints() const { return grid_pos; }

Expand Down
18 changes: 5 additions & 13 deletions src/libxtp/aoshell.cc
Original file line number Diff line number Diff line change
Expand Up @@ -88,13 +88,14 @@ void AOShell::normalizeContraction() {
return;
}

void AOShell::EvalAOspace(Eigen::VectorBlock<Eigen::VectorXd>& AOvalues,
Eigen::Block<Eigen::MatrixX3d>& gradAOvalues,
const Eigen::Vector3d& grid_pos) const {
AOShell::AOValues AOShell::EvalAOspace(const Eigen::Vector3d& grid_pos) const {

// need position of shell
const Eigen::Vector3d center = (grid_pos - pos_);
const double distsq = center.squaredNorm();
AOShell::AOValues AO(getNumFunc());
Eigen::VectorXd& AOvalues = AO.values;
Eigen::MatrixX3d& gradAOvalues = AO.derivatives;

// iterate over Gaussians in this shell
for (const AOGaussianPrimitive& gaussian : gaussians_) {
Expand Down Expand Up @@ -299,16 +300,7 @@ void AOShell::EvalAOspace(Eigen::VectorBlock<Eigen::VectorXd>& AOvalues,
break;
}
} // contractions
return;
} // namespace xtp

void AOShell::EvalAOspace(Eigen::VectorBlock<Eigen::VectorXd>& AOvalues,
const Eigen::Vector3d& grid_pos) const {

Eigen::MatrixX3d temp = Eigen::MatrixX3d::Zero(AOvalues.size(), 3);
Eigen::Block<Eigen::MatrixX3d> temp2 =
temp.block(0, 0, temp.rows(), temp.cols());
EvalAOspace(AOvalues, temp2, grid_pos);
return AO;
}

std::ostream& operator<<(std::ostream& out, const AOShell& shell) {
Expand Down
18 changes: 8 additions & 10 deletions src/libxtp/bsecoupling.cc
Original file line number Diff line number Diff line change
Expand Up @@ -139,14 +139,13 @@ Eigen::MatrixXd BSECoupling::SetupCTStates(Index bseA_vtotal, Index bseB_vtotal,

Index noAB = occA_ * unoccB_;
Index noBA = unoccA_ * occB_;
Index bseAB_total = bseAB_vtotal + bseAB_ctotal;
Index bseAB_size = bseAB_vtotal * bseAB_ctotal;
Eigen::MatrixXd CTstates = Eigen::MatrixXd::Zero(bseAB_size, noAB + noBA);

auto A_occ = A_AB.block(0, bseA_vtotal - occA_, bseAB_total, occA_);
auto A_unocc = A_AB.block(0, bseA_vtotal, bseAB_total, unoccA_);
auto B_occ = B_AB.block(0, bseB_vtotal - occB_, bseAB_total, occB_);
auto B_unocc = B_AB.block(0, bseB_vtotal, bseAB_total, unoccB_);
auto A_occ = A_AB.middleCols(bseA_vtotal - occA_, occA_);
auto A_unocc = A_AB.middleCols(bseA_vtotal, unoccA_);
auto B_occ = B_AB.middleCols(bseB_vtotal - occB_, occB_);
auto B_unocc = B_AB.middleCols(bseB_vtotal, unoccB_);

const Eigen::MatrixXd A_occ_occ = A_occ.topRows(bseAB_vtotal);
const Eigen::MatrixXd B_unocc_unocc = B_unocc.bottomRows(bseAB_ctotal);
Expand Down Expand Up @@ -348,7 +347,6 @@ void BSECoupling::CalculateCouplings(const Orbitals& orbitalsA,
Index bseAB_cmin = orbitalsAB.getBSEcmin();
Index bseAB_vmax = orbitalsAB.getBSEvmax();
Index bseAB_vmin = orbitalsAB.getBSEvmin();
Index basisAB = orbitalsAB.getBasisSetSize();
Index bseAB_vtotal = bseAB_vmax - bseAB_vmin + 1;
Index bseAB_ctotal = bseAB_cmax - bseAB_cmin + 1;
Index bseAB_total = bseAB_vtotal + bseAB_ctotal;
Expand All @@ -366,11 +364,11 @@ void BSECoupling::CalculateCouplings(const Orbitals& orbitalsA,
<< " to " << bseAB_cmax << " total: " << bseAB_total << flush;

Eigen::MatrixXd MOsA =
orbitalsA.MOs().eigenvectors().block(0, bseA_vmin, basisA, bseA_total);
orbitalsA.MOs().eigenvectors().middleCols(bseA_vmin, bseA_total);
Eigen::MatrixXd MOsB =
orbitalsB.MOs().eigenvectors().block(0, bseB_vmin, basisB, bseB_total);
Eigen::MatrixXd MOsAB = orbitalsAB.MOs().eigenvectors().block(
0, bseAB_vmin, basisAB, bseAB_total);
orbitalsB.MOs().eigenvectors().middleCols(bseB_vmin, bseB_total);
Eigen::MatrixXd MOsAB =
orbitalsAB.MOs().eigenvectors().middleCols(bseAB_vmin, bseAB_total);

XTP_LOG(Log::info, *pLog_)
<< TimeStamp() << " Calculating overlap matrix for basisset: "
Expand Down
8 changes: 4 additions & 4 deletions src/libxtp/davidsonsolver.cc
Original file line number Diff line number Diff line change
Expand Up @@ -455,8 +455,8 @@ void DavidsonSolver::gramschmidt(Eigen::MatrixXd &Q, Index nstart) const {
// orthogonalize vectors to each other
for (Index j = nstart + 1; j < Q.cols(); ++j) {
Index range = j - nstart;
Q.col(j) -= Q.block(0, nstart, Q.rows(), range) *
(Q.block(0, nstart, Q.rows(), range).transpose() * Q.col(j));
Q.col(j) -= Q.middleCols(nstart, range) *
(Q.middleCols(nstart, range).transpose() * Q.col(j));
Q.col(j).normalize();
}
// repeat again two is enough GS
Expand All @@ -470,8 +470,8 @@ void DavidsonSolver::gramschmidt(Eigen::MatrixXd &Q, Index nstart) const {

for (Index j = nstart + 1; j < Q.cols(); ++j) {
Index range = j - nstart;
Q.col(j) -= Q.block(0, nstart, Q.rows(), range) *
(Q.block(0, nstart, Q.rows(), range).transpose() * Q.col(j));
Q.col(j) -= Q.middleCols(nstart, range) *
(Q.middleCols(nstart, range).transpose() * Q.col(j));
if (Q.col(j).norm() <= 1E-12 * norms(range)) {
// info_ = Eigen::ComputationInfo::NumericalIssue;
throw std::runtime_error("Linear dependencies in Gram-Schmidt.");
Expand Down
8 changes: 4 additions & 4 deletions src/libxtp/dftcoupling.cc
Original file line number Diff line number Diff line change
Expand Up @@ -173,10 +173,10 @@ void DFTcoupling::CalculateCouplings(const Orbitals& orbitalsA,
}

// constructing merged orbitals
auto MOsA = orbitalsA.MOs().eigenvectors().block(0, Range_orbA.first, basisA,
Range_orbA.second);
auto MOsB = orbitalsB.MOs().eigenvectors().block(0, Range_orbB.first, basisB,
Range_orbB.second);
auto MOsA = orbitalsA.MOs().eigenvectors().middleCols(Range_orbA.first,
Range_orbA.second);
auto MOsB = orbitalsB.MOs().eigenvectors().middleCols(Range_orbB.first,
Range_orbB.second);

XTP_LOG(Log::info, *pLog_) << "Calculating overlap matrix for basisset: "
<< orbitalsAB.getDFTbasisName() << flush;
Expand Down
3 changes: 1 addition & 2 deletions src/libxtp/dftengine/dftengine.cc
Original file line number Diff line number Diff line change
Expand Up @@ -144,8 +144,7 @@ std::array<Eigen::MatrixXd, 2> DFTEngine::CalcERIs_EXX(
if (conv_accelerator_.getUseMixing() || MOCoeff.rows() == 0) {
return ERIs_.CalculateERIs_EXX_3c(Eigen::MatrixXd::Zero(0, 0), Dmat);
} else {
Eigen::MatrixXd occblock =
MOCoeff.block(0, 0, MOCoeff.rows(), numofelectrons_ / 2);
Eigen::MatrixXd occblock = MOCoeff.leftCols(numofelectrons_ / 2);
return ERIs_.CalculateERIs_EXX_3c(occblock, Dmat);
}
} else {
Expand Down
26 changes: 7 additions & 19 deletions src/libxtp/grids/gridbox.cc
Original file line number Diff line number Diff line change
Expand Up @@ -41,27 +41,15 @@ void GridBox::FindSignificantShells(const AOBasis& basis) {
}
}

Eigen::VectorXd GridBox::CalcAOValue_and_Grad(
Eigen::MatrixX3d& ao_grad, const Eigen::Vector3d& point) const {
Eigen::VectorXd ao = Eigen::VectorXd::Zero(Matrixsize());
AOShell::AOValues GridBox::CalcAOValues(const Eigen::Vector3d& point) const {
AOShell::AOValues result(Matrixsize());
for (Index j = 0; j < Shellsize(); ++j) {
Eigen::Block<Eigen::MatrixX3d> grad_block =
ao_grad.block(aoranges[j].start, 0, aoranges[j].size, 3);
Eigen::VectorBlock<Eigen::VectorXd> ao_block =
ao.segment(aoranges[j].start, aoranges[j].size);
significant_shells[j]->EvalAOspace(ao_block, grad_block, point);
const AOShell::AOValues val = significant_shells[j]->EvalAOspace(point);
result.derivatives.middleRows(aoranges[j].start, aoranges[j].size) =
val.derivatives;
result.values.segment(aoranges[j].start, aoranges[j].size) = val.values;
}
return ao;
}

Eigen::VectorXd GridBox::CalcAOValues(const Eigen::Vector3d& pos) const {
Eigen::VectorXd ao = Eigen::VectorXd::Zero(Matrixsize());
for (Index j = 0; j < Shellsize(); ++j) {
Eigen::VectorBlock<Eigen::VectorXd> ao_block =
ao.segment(aoranges[j].start, aoranges[j].size);
significant_shells[j]->EvalAOspace(ao_block, pos);
}
return ao;
return result;
}

void GridBox::AddtoBigMatrix(Eigen::MatrixXd& bigmatrix,
Expand Down
15 changes: 7 additions & 8 deletions src/libxtp/gwbse/bse_operator.cc
Original file line number Diff line number Diff line change
Expand Up @@ -68,23 +68,23 @@ Eigen::MatrixXd BSE_OPERATOR<cqp, cx, cd, cd2>::matmul(
// holds a reference to it
Eigen::MatrixXd Temp;
if (cd != 0) {
Temp = -cd * (Mmn_[c1 + cmin].block(cmin, 0, bse_ctotal_, auxsize));
Temp = -cd * (Mmn_[c1 + cmin].middleRows(cmin, bse_ctotal_));
transform.PrepareMatrix1(Temp, threadid);
} else if (cd2 != 0) {
Temp = -cd2 * (Mmn_[c1 + cmin].block(vmin, 0, bse_vtotal_, auxsize));
Temp = -cd2 * (Mmn_[c1 + cmin].middleRows(vmin, bse_vtotal_));
transform.PrepareMatrix1(Temp, threadid);
}

for (Index v1 = 0; v1 < bse_vtotal_; v1++) {
transform.SetTempZero(threadid);
if (cd != 0) {
transform.PrepareMatrix2(
Mmn_[v1 + vmin].block(vmin, 0, bse_vtotal_, auxsize), cd2 != 0,
Mmn_[v1 + vmin].middleRows(vmin, bse_vtotal_), cd2 != 0,
threadid);
}
if (cd2 != 0) {
transform.PrepareMatrix2(
Mmn_[v1 + vmin].block(cmin, 0, bse_ctotal_, auxsize), cd2 != 0,
Mmn_[v1 + vmin].middleRows(cmin, bse_ctotal_), cd2 != 0,
threadid);
}
if (cqp != 0) {
Expand All @@ -104,13 +104,12 @@ Eigen::MatrixXd BSE_OPERATOR<cqp, cx, cd, cd2>::matmul(
#pragma omp for schedule(dynamic)
for (Index v1 = 0; v1 < bse_vtotal_; v1++) {
Index va = v1 + vmin;
Eigen::MatrixXd Mmn1 =
cx * Mmn_[va].block(cmin, 0, bse_ctotal_, auxsize);
Eigen::MatrixXd Mmn1 = cx * Mmn_[va].middleRows(cmin, bse_ctotal_);
transform.PushMatrix1(Mmn1, threadid);
for (Index v2 = v1; v2 < bse_vtotal_; v2++) {
Index vb = v2 + vmin;
transform.MultiplyBlocks(
Mmn_[vb].block(cmin, 0, bse_ctotal_, auxsize), v1, v2, threadid);
transform.MultiplyBlocks(Mmn_[vb].middleRows(cmin, bse_ctotal_), v1,
v2, threadid);
}
}
}
Expand Down
3 changes: 1 addition & 2 deletions src/libxtp/gwbse/gwbse.cc
Original file line number Diff line number Diff line change
Expand Up @@ -532,9 +532,8 @@ Eigen::MatrixXd GWBSE::CalculateVXC(const AOBasis& dftbasis) {
<< TimeStamp() << " Set hybrid exchange factor: " << orbitals_.getScaHFX()
<< flush;
Index qptotal = gwopt_.qpmax - gwopt_.qpmin + 1;
Index basissize = Index(orbitals_.MOs().eigenvectors().rows());
Eigen::MatrixXd mos =
orbitals_.MOs().eigenvectors().block(0, gwopt_.qpmin, basissize, qptotal);
orbitals_.MOs().eigenvectors().middleCols(gwopt_.qpmin, qptotal);

Eigen::MatrixXd vxc = mos.transpose() * e_vxc_ao.matrix() * mos;
XTP_LOG(Log::error, *pLog_)
Expand Down
5 changes: 2 additions & 3 deletions src/libxtp/gwbse/rpa.cc
Original file line number Diff line number Diff line change
Expand Up @@ -222,19 +222,18 @@ Eigen::MatrixXd RPA::Calculate_H2p_ApB() const {
const Index n_occ = lumo - rpamin_;
const Index n_unocc = rpamax_ - lumo + 1;
const Index rpasize = n_occ * n_unocc;
const Index auxsize = Mmn_.auxsize();
vc2index vc = vc2index(0, 0, n_unocc);
Eigen::MatrixXd ApB = Eigen::MatrixXd::Zero(rpasize, rpasize);
#pragma omp parallel for schedule(guided)
for (Index v2 = 0; v2 < n_occ; v2++) {
Index i2 = vc.I(v2, 0);
const Eigen::MatrixXd Mmn_v2T =
Mmn_[v2].block(n_occ, 0, n_unocc, auxsize).transpose();
Mmn_[v2].middleRows(n_occ, n_unocc).transpose();
for (Index v1 = v2; v1 < n_occ; v1++) {
Index i1 = vc.I(v1, 0);
// Multiply with factor 2 to sum over both (identical) spin states
ApB.block(i1, i2, n_unocc, n_unocc) =
2 * 2 * Mmn_[v1].block(n_occ, 0, n_unocc, auxsize) * Mmn_v2T;
2 * 2 * Mmn_[v1].middleRows(n_occ, n_unocc) * Mmn_v2T;
}
}
ApB.diagonal() += Calculate_H2p_AmB();
Expand Down
4 changes: 2 additions & 2 deletions src/libxtp/numerical_integration/amplitude_integration.cc
Original file line number Diff line number Diff line change
Expand Up @@ -42,8 +42,8 @@ std::vector<std::vector<double> >
const std::vector<double>& weights = box.getGridWeights();
// iterate over gridpoints
for (Index p = 0; p < box.size(); p++) {
Eigen::VectorXd ao = box.CalcAOValues(points[p]);
result[i][p] = weights[p] * amplitude_here.dot(ao);
AOShell::AOValues ao = box.CalcAOValues(points[p]);
result[i][p] = weights[p] * amplitude_here.dot(ao.values);
}
}
return result;
Expand Down
8 changes: 4 additions & 4 deletions src/libxtp/numerical_integration/density_integration.cc
Original file line number Diff line number Diff line change
Expand Up @@ -84,8 +84,8 @@ double DensityIntegration<Grid>::IntegrateDensity(
const std::vector<double>& weights = box.getGridWeights();
// iterate over gridpoints
for (Index p = 0; p < box.size(); p++) {
Eigen::VectorXd ao = box.CalcAOValues(points[p]);
double rho = (ao.transpose() * DMAT_here * ao)(0, 0) * weights[p];
AOShell::AOValues ao = box.CalcAOValues(points[p]);
double rho = ao.values.dot(DMAT_here * ao.values) * weights[p];
densities_[i][p] = rho;
N += rho;
}
Expand Down Expand Up @@ -113,8 +113,8 @@ Gyrationtensor DensityIntegration<Grid>::IntegrateGyrationTensor(
const std::vector<double>& weights = box.getGridWeights();
// iterate over gridpoints
for (Index p = 0; p < box.size(); p++) {
Eigen::VectorXd ao = box.CalcAOValues(points[p]);
double rho = (ao.transpose() * DMAT_here * ao).value() * weights[p];
AOShell::AOValues ao = box.CalcAOValues(points[p]);
double rho = ao.values.dot(DMAT_here * ao.values) * weights[p];
densities_[i][p] = rho;
N += rho;
centroid += rho * points[p];
Expand Down
15 changes: 8 additions & 7 deletions src/libxtp/numerical_integration/vxc_potential.cc
Original file line number Diff line number Diff line change
Expand Up @@ -182,20 +182,21 @@ Mat_p_Energy Vxc_Potential<Grid>::IntegrateVXC(

// iterate over gridpoints
for (Index p = 0; p < box.size(); p++) {
Eigen::MatrixX3d ao_grad = Eigen::MatrixX3d::Zero(box.Matrixsize(), 3);
Eigen::VectorXd ao = box.CalcAOValue_and_Grad(ao_grad, points[p]);
const double rho = 0.5 * (ao.transpose() * DMAT_symm * ao).value();
AOShell::AOValues ao = box.CalcAOValues(points[p]);
const double rho = 0.5 * ao.values.dot(DMAT_symm * ao.values);
const double weight = weights[p];
if (rho * weight < 1.e-20) {
continue; // skip the rest, if density is very small
}
const Eigen::Vector3d rho_grad = ao.transpose() * DMAT_symm * ao_grad;
const Eigen::Vector3d rho_grad =
ao.values.transpose() * DMAT_symm * ao.derivatives;
const double sigma = (rho_grad.transpose() * rho_grad).value();
const Eigen::VectorXd grad = ao_grad * rho_grad;
const Eigen::VectorXd grad = ao.derivatives * rho_grad;
typename Vxc_Potential<Grid>::XC_entry xc = EvaluateXC(rho, sigma);
EXC_box += weight * rho * xc.f_xc;
auto addXC = weight * (0.5 * xc.df_drho * ao + 2.0 * xc.df_dsigma * grad);
Vxc_here.noalias() += addXC * ao.transpose();
auto addXC =
weight * (0.5 * xc.df_drho * ao.values + 2.0 * xc.df_dsigma * grad);
Vxc_here.noalias() += addXC * ao.values.transpose();
}
box.AddtoBigMatrix(vxc.matrix(), Vxc_here);
vxc.energy() += EXC_box;
Expand Down
Loading