Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Some changes #89

Merged
merged 4 commits into from
Jan 15, 2024
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
31 changes: 19 additions & 12 deletions src/basis.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -2318,7 +2318,23 @@ arma::mat BasisSet::kinetic() const {
}

arma::mat BasisSet::nuclear() const {
// Form nuclear attraction matrix
std::vector<std::tuple<int,double,double,double>> nuclear_data;
for(size_t inuc=0;inuc<nuclei.size();inuc++) {
if(nuclei[inuc].bsse)
continue;
// Nuclear charge
int Z=nuclei[inuc].Z;

// Coordinates of nucleus
double cx=nuclei[inuc].r.x;
double cy=nuclei[inuc].r.y;
double cz=nuclei[inuc].r.z;
nuclear_data.push_back(std::make_tuple(Z,cx,cy,cz));
}
return nuclear(nuclear_data);
}

arma::mat BasisSet::nuclear(const std::vector<std::tuple<int,double,double,double>> & nuclear_data) const {

// Size of basis set
size_t N=get_Nbf();
Expand All @@ -2332,18 +2348,9 @@ arma::mat BasisSet::nuclear() const {
#pragma omp parallel for schedule(dynamic)
#endif
for(size_t ip=0;ip<shellpairs.size();ip++)
for(size_t inuc=0;inuc<nuclei.size();inuc++) {
// If BSSE nucleus, do nothing
if(nuclei[inuc].bsse)
continue;

// Nuclear charge
int Z=nuclei[inuc].Z;
for(size_t inuc=0;inuc<nuclear_data.size();inuc++) {

// Coordinates of nucleus
double cx=nuclei[inuc].r.x;
double cy=nuclei[inuc].r.y;
double cz=nuclei[inuc].r.z;
auto [Z, cx, cy, cz] = nuclear_data[inuc];

// Shells in pair
size_t i=shellpairs[ip].is;
Expand Down
2 changes: 2 additions & 0 deletions src/basis.h
Original file line number Diff line number Diff line change
Expand Up @@ -447,6 +447,8 @@ class BasisSet {
arma::mat kinetic() const;
/// Calculate nuclear repulsion matrix
arma::mat nuclear() const;
/// Calculate nuclear repulsion matrix with given nuclei
arma::mat nuclear(const std::vector<std::tuple<int,double,double,double>> & nuclei) const;
/// Calculate electric potential matrix
arma::mat potential(coords_t r) const;
/// Calculate superposition of atomic potentials
Expand Down
35 changes: 12 additions & 23 deletions src/density_fitting.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -52,7 +52,7 @@ bool DensityFit::Bmat_enabled() const {
return Bmat;
}

size_t DensityFit::fill(const BasisSet & orbbas, const BasisSet & auxbas, bool dir, double erithr, double linthr, bool bmat) {
size_t DensityFit::fill(const BasisSet & orbbas, const BasisSet & auxbas, bool dir, double erithr, double linthr, double cholthr, bool bmat) {
// Construct density fitting basis

// Store amount of functions
Expand Down Expand Up @@ -125,28 +125,8 @@ size_t DensityFit::fill(const BasisSet & orbbas, const BasisSet & auxbas, bool d
}

if(Bmat) {
// Form ab^-1 and ab^-1/2
arma::mat abvec;
arma::vec abval;
eig_sym_ordered(abval,abvec,ab);

// Count linearly independent vectors
size_t Nind=0;
for(size_t i=0;i<abval.n_elem;i++)
if(abval(i)>=linthr)
Nind++;

// and drop the linearly dependent ones
abval=abval.subvec(abval.n_elem-Nind,abval.n_elem-1);
abvec=abvec.cols(abvec.n_cols-Nind,abvec.n_cols-1);

// Form matrices
ab_inv.zeros(abvec.n_rows,abvec.n_rows);
ab_invh.zeros(abvec.n_rows,abvec.n_rows);
for(size_t i=0;i<abval.n_elem;i++) {
ab_inv+=abvec.col(i)*arma::trans(abvec.col(i))/abval(i);
ab_invh+=abvec.col(i)*arma::trans(abvec.col(i))/sqrt(abval(i));
}
ab_invh = PartialCholeskyOrth(ab, cholthr, linthr);
ab_inv = ab_invh * ab_invh.t();
} else {
// Just RI-J(K), so use faster method from Eichkorn et al to form ab_inv only
ab_inv=arma::inv(ab + DELTA*arma::eye(ab.n_rows,ab.n_cols));
Expand Down Expand Up @@ -803,6 +783,15 @@ arma::mat DensityFit::calcJ(const arma::mat & P) const {

// Get the expansion coefficients
arma::vec c=compute_expansion(P);
return digestJ(c);
}

arma::mat DensityFit::calcJ(const arma::vec & c) const {
if(c.n_elem != Naux) {
std::ostringstream oss;
oss << "Error in DensityFit: Naux = " << Naux << ", c.n_elem = " << c.n_elem << "!\n";
throw std::logic_error(oss.str());
}

arma::mat J(Nbf,Nbf);
J.zeros();
Expand Down
4 changes: 3 additions & 1 deletion src/density_fitting.h
Original file line number Diff line number Diff line change
Expand Up @@ -135,7 +135,7 @@ class DensityFit {
* HF routine should be more tolerant of linear dependencies in the basis.
* Returns amount of significant orbital shell pairs.
*/
size_t fill(const BasisSet & orbbas, const BasisSet & auxbas, bool direct, double erithr, double linthr, bool bmat=false);
size_t fill(const BasisSet & orbbas, const BasisSet & auxbas, bool direct, double erithr, double linthr, double cholthr, bool bmat=false);

/// Compute estimate of necessary memory
size_t memory_estimate(const BasisSet & orbbas, const BasisSet & auxbas, double erithr, bool direct) const;
Expand All @@ -149,6 +149,8 @@ class DensityFit {
arma::mat calcJ(const arma::mat & P) const;
/// Get Coulomb matrix from P
std::vector<arma::mat> calcJ(const std::vector<arma::mat> & P) const;
/// Digest J matrix from computed expansion
arma::mat digest(const arma::vec & gamma) const;

/// Calculate force from P
arma::vec forceJ(const arma::mat & P);
Expand Down
3 changes: 3 additions & 0 deletions src/linalg.h
Original file line number Diff line number Diff line change
Expand Up @@ -104,6 +104,9 @@ arma::mat CanonicalOrth(const arma::mat & S, double cutoff);
/// Same, but use computed decomposition
arma::mat CanonicalOrth(const arma::mat & Svec, const arma::vec & Sval, double cutoff);

/// Pivoted Cholesky orthogonalization
arma::mat PartialCholeskyOrth(const arma::mat & S, double cholcut, double scut);

/// Automatic orthonormalization.
arma::mat BasOrth(const arma::mat & S, bool verbose);
/// Orthogonalize basis
Expand Down
6 changes: 3 additions & 3 deletions src/localization.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -1203,11 +1203,11 @@ void Pipek::cost_func_der(const arma::cx_mat & Wv, double & Dinv, arma::cx_mat &
Edmiston::Edmiston(const BasisSet & basis, const BasisSet & fitbas, const arma::mat & Cv, bool delocalize) : UnitaryFunction(4,!delocalize) {
// Store orbitals
C=Cv;
// Initialize fitting integrals. Direct computation, linear dependence threshold 1e-8, use Hartree-Fock routine since it has better tolerance for linear dependencies
// Initialize fitting integrals. Direct computation, linear dependence threshold 1e-8, Cholesky threshold 1e-9, use Hartree-Fock routine since it has better tolerance for linear dependencies
if(!fitbas.get_Nbf())
dfit.fill(basis,basis.density_fitting(),true,1e-8,false);
dfit.fill(basis,basis.density_fitting(),true,1e-8,1e-9,false);
else
dfit.fill(basis,fitbas,true,1e-8,false);
dfit.fill(basis,fitbas,true,1e-8,1e-9,false);

use_chol=false;
}
Expand Down
6 changes: 4 additions & 2 deletions src/scf-base.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -144,6 +144,8 @@ SCF::SCF(const BasisSet & basis, Checkpoint & chkpt) {
fitmem=1000000*settings.get_int("FittingMemory");
// Linear dependence threshold
fitthr=settings.get_double("FittingThreshold");
// Linear dependence threshold
fitcholthr=settings.get_double("FittingCholeskyThreshold");

// Use Cholesky?
cholesky=settings.get_bool("Cholesky");
Expand Down Expand Up @@ -347,7 +349,7 @@ SCF::SCF(const BasisSet & basis, Checkpoint & chkpt) {
}

// Calculate the fitting integrals, don't run in B-matrix mode
size_t Npairs=dfit.fill(*basisp,dfitbas,direct,intthr,fitthr,false);
size_t Npairs=dfit.fill(*basisp,dfitbas,direct,intthr,fitthr,fitcholthr,false);

if(verbose) {
printf("done (%s)\n",t.elapsed().c_str());
Expand Down Expand Up @@ -505,7 +507,7 @@ void SCF::fill_rs(double omega) {
}

t.set();
size_t Npairs=dfit_rs.fill(*basisp,dfitbas,direct,intthr,fitthr,dfit.Bmat_enabled());
size_t Npairs=dfit_rs.fill(*basisp,dfitbas,direct,intthr,fitthr,fitcholthr,dfit.Bmat_enabled());
if(verbose) {
printf("done (%s)\n",t.elapsed().c_str());
printf("%i shell pairs out of %i are significant.\n",(int) Npairs, (int) basisp->get_unique_shellpairs().size());
Expand Down
2 changes: 2 additions & 0 deletions src/scf.h
Original file line number Diff line number Diff line change
Expand Up @@ -287,6 +287,8 @@ class SCF {
size_t fitmem;
/// Threshold for density fitting
double fitthr;
/// Pivoted Cholesky threshold for density fitting
double fitcholthr;

/// Cholesky calculation?
bool cholesky;
Expand Down
1 change: 1 addition & 0 deletions src/settings.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -142,6 +142,7 @@ void Settings::add_scf_settings() {
add_int("FittingMemory", "Amount of memory in MB to use for exchange fitting",1000);
// Threshold for screening eigenvectors
add_double("FittingThreshold", "Linear dependence threshold for Coulomb integrals in density fitting",1e-8);
add_double("FittingCholeskyThreshold", "Linear dependence threshold for pivoted Cholesky of Coulomb integrals in density fitting",1e-9);

// SAP basis
add_string("SAPBasis", "Tabulated atomic effective potential \"basis set\"","helfem_large.gbs");
Expand Down
Loading