Skip to content

Commit

Permalink
Add helper functions to evaluate exchange matrix with occupations giv…
Browse files Browse the repository at this point in the history
…en as armadillo vectors
  • Loading branch information
susilehtola committed Feb 2, 2024
1 parent 6373505 commit af9e19b
Show file tree
Hide file tree
Showing 4 changed files with 25 additions and 0 deletions.
8 changes: 8 additions & 0 deletions src/density_fitting.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -1095,6 +1095,10 @@ arma::mat DensityFit::calcK(const arma::mat & Corig, const std::vector<double> &
return K;
}

arma::mat DensityFit::calcK(const arma::mat & Corig, const arma::vec & occo, size_t fitmem) const {
return calcK(Corig, arma::conv_to<std::vector<double>>::from(occo), fitmem);
}

arma::cx_mat DensityFit::calcK(const arma::cx_mat & Corig, const std::vector<double> & occo, size_t fitmem) const {
if(Corig.n_rows != Nbf) {
std::ostringstream oss;
Expand Down Expand Up @@ -1136,6 +1140,10 @@ arma::cx_mat DensityFit::calcK(const arma::cx_mat & Corig, const std::vector<dou
return K;
}

arma::cx_mat DensityFit::calcK(const arma::cx_mat & Corig, const arma::vec & occo, size_t fitmem) const {
return calcK(Corig, arma::conv_to<std::vector<double>>::from(occo), fitmem);
}

size_t DensityFit::get_Norb() const {
return Nbf;
}
Expand Down
5 changes: 5 additions & 0 deletions src/density_fitting.h
Original file line number Diff line number Diff line change
Expand Up @@ -155,6 +155,11 @@ class DensityFit {
/// Get exchange matrix from orbitals with occupation numbers occs
arma::cx_mat calcK(const arma::cx_mat & C, const std::vector<double> & occs, size_t fitmem) const;

/// Get exchange matrix from orbitals with occupation numbers occs
arma::mat calcK(const arma::mat & C, const arma::vec & occs, size_t fitmem) const;
/// Get exchange matrix from orbitals with occupation numbers occs
arma::cx_mat calcK(const arma::cx_mat & C, const arma::vec & occs, size_t fitmem) const;

/// Get the number of orbital functions
size_t get_Norb() const;
/// Get the number of auxiliary functions
Expand Down
8 changes: 8 additions & 0 deletions src/erichol.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -768,6 +768,10 @@ arma::mat ERIchol::calcK(const arma::mat & C, const std::vector<double> & occs)
return K;
}

arma::mat ERIchol::calcK(const arma::mat & C, const arma::vec & occs) const {
return calcK(C, arma::conv_to<std::vector<double>>::from(occs));
}

arma::cx_mat ERIchol::calcK(const arma::cx_mat & C, const std::vector<double> & occs) const {
arma::cx_mat K(C.n_rows,C.n_rows);
K.zeros();
Expand All @@ -788,6 +792,10 @@ arma::cx_mat ERIchol::calcK(const arma::cx_mat & C, const std::vector<double> &
return K;
}

arma::cx_mat ERIchol::calcK(const arma::cx_mat & C, const arma::vec & occs) const {
return calcK(C, arma::conv_to<std::vector<double>>::from(occs));
}

void ERIchol::B_matrix(arma::mat & Br) const {
Br.zeros(Nbf*Nbf,B.n_cols);
for(size_t P=0;P<B.n_cols;P++)
Expand Down
4 changes: 4 additions & 0 deletions src/erichol.h
Original file line number Diff line number Diff line change
Expand Up @@ -96,11 +96,15 @@ class ERIchol {
arma::mat calcK(const arma::vec & C) const;
/// Form exchange matrix
arma::mat calcK(const arma::mat & C, const std::vector<double> & occs) const;
/// Form exchange matrix
arma::mat calcK(const arma::mat & C, const arma::vec & occs) const;

/// Form exchange matrix
arma::cx_mat calcK(const arma::cx_vec & C) const;
/// Form exchange matrix
arma::cx_mat calcK(const arma::cx_mat & C, const std::vector<double> & occs) const;
/// Form exchange matrix
arma::cx_mat calcK(const arma::cx_mat & C, const arma::vec & occs) const;

/// Get full B matrix
void B_matrix(arma::mat & B) const;
Expand Down

0 comments on commit af9e19b

Please sign in to comment.