33#include <unordered_map>
37#ifdef LIBINT2_HAVE_BTAS
40#error "libint2::lcao requires BTAS"
44#include <libint2/chemistry/sto3g_atomic_density.h>
45#include <libint2/diis.h>
46#include <libint2/util/intpart_iter.h>
54typedef btas::RangeNd<CblasRowMajor, std::array<long, 2>> Range2;
55typedef btas::RangeNd<CblasRowMajor, std::array<long, 3>> Range3;
56typedef btas::RangeNd<CblasRowMajor, std::array<long, 4>> Range4;
57typedef btas::Tensor<double, Range2> Tensor2d;
58typedef btas::Tensor<double, Range3> Tensor3d;
59typedef btas::Tensor<double, Range3> Tensor4d;
61typedef Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>
66typedef Eigen::DiagonalMatrix<double, Eigen::Dynamic, Eigen::Dynamic>
69using libint2::BasisSet;
71using libint2::Operator;
73using libint2::libint2::Atom;
75std::vector<libint2::Atom> read_geometry(
const std::string& filename);
76Matrix compute_soad(
const std::vector<libint2::Atom>& atoms);
78Matrix compute_shellblock_norm(
const BasisSet& obs,
const Matrix& A);
80template <Operator obtype>
81std::array<Matrix, libint2::operator_traits<obtype>::nopers> compute_1body_ints(
83 const std::vector<libint2::Atom>& atoms = std::vector<libint2::Atom>());
85#if LIBINT2_DERIV_ONEBODY_ORDER
86template <Operator obtype>
87std::vector<Matrix> compute_1body_ints_deriv(
88 unsigned deriv_order,
const BasisSet& obs,
89 const std::vector<libint2::Atom>& atoms);
92template <lib
int2::Operator Kernel = lib
int2::Operator::coulomb>
93Matrix compute_schwarz_ints(
94 const BasisSet& bs1,
const BasisSet& bs2 = BasisSet(),
95 bool use_2norm =
false,
96 typename libint2::operator_traits<Kernel>::oper_params_type params =
97 libint2::operator_traits<Kernel>::default_params());
98Matrix compute_do_ints(
const BasisSet& bs1,
const BasisSet& bs2 = BasisSet(),
99 bool use_2norm =
false
102using shellpair_list_t = std::unordered_map<size_t, std::vector<size_t>>;
103shellpair_list_t obs_shellpair_list;
109shellpair_list_t compute_shellpair_list(
const BasisSet& bs1,
110 const BasisSet& bs2 = BasisSet(),
111 double threshold = 1e-12);
113Matrix compute_2body_fock(
114 const BasisSet& obs,
const Matrix& D,
115 double precision = std::numeric_limits<
117 const Matrix& Schwarz = Matrix()
121Matrix compute_2body_fock_general(
122 const BasisSet& obs,
const Matrix& D,
const BasisSet& D_bs,
123 bool D_is_sheldiagonal =
false,
124 double precision = std::numeric_limits<
128#if LIBINT2_DERIV_ERI_ORDER
129template <
unsigned deriv_order>
130std::vector<Matrix> compute_2body_fock_deriv(
131 const BasisSet& obs,
const std::vector<libint2::Atom>& atoms,
133 double precision = std::numeric_limits<
135 const Matrix& Schwarz = Matrix()
145std::tuple<Matrix, Matrix, double> conditioning_orthogonalizer(
146 const Matrix& S,
double S_condition_number_threshold);
148#ifdef LIBINT2_HAVE_BTAS
149#define HAVE_DENSITY_FITTING 1
152 const BasisSet& dfbs;
153 DFFockEngine(
const BasisSet& _obs,
const BasisSet& _dfbs)
154 : obs(_obs), dfbs(_dfbs) {}
156 typedef btas::RangeNd<CblasRowMajor, std::array<long, 3>> Range3d;
157 typedef btas::Tensor<double, Range3d> Tensor3d;
161 Matrix compute_2body_fock_dfC(
const Matrix& Cocc);
169template <
typename Lambda>
174 auto thread_id = omp_get_thread_num();
178 std::vector<std::thread> threads;
179 for (
int thread_id = 0; thread_id != libint2::nthreads; ++thread_id) {
180 if (thread_id != nthreads - 1)
181 threads.push_back(std::thread(lambda, thread_id));
185 for (
int thread_id = 0; thread_id < nthreads - 1; ++thread_id)
186 threads[thread_id].join();
191int main(
int argc,
char* argv[]) {
203 const auto filename = (argc > 1) ? argv[1] :
"h2o.xyz";
204 const auto basisname = (argc > 2) ? argv[2] :
"aug-cc-pVDZ";
205 bool do_density_fitting =
false;
206#ifdef HAVE_DENSITY_FITTING
207 do_density_fitting = (argc > 3);
208 const auto dfbasisname = do_density_fitting ? argv[3] :
"";
210 std::vector<libint2::Atom> atoms = read_geometry(filename);
214 using libint2::nthreads;
215 auto nthreads_cstr = getenv(
"LIBINT_NUM_THREADS");
217 if (nthreads_cstr && strcmp(nthreads_cstr,
"")) {
218 std::istringstream iss(nthreads_cstr);
220 if (nthreads > 1 << 16 || nthreads <= 0) nthreads = 1;
223 omp_set_num_threads(nthreads);
225 std::cout <<
"Will scale over " << nthreads
231 <<
" threads" << std::endl;
236 for (
auto i = 0; i < atoms.size(); ++i) nelectron += atoms[i].atomic_number;
237 const auto ndocc = nelectron / 2;
238 cout <<
"# of electrons = " << nelectron << endl;
242 for (
auto i = 0; i < atoms.size(); i++)
243 for (
auto j = i + 1; j < atoms.size(); j++) {
244 auto xij = atoms[i].x - atoms[j].x;
245 auto yij = atoms[i].y - atoms[j].y;
246 auto zij = atoms[i].z - atoms[j].z;
247 auto r2 = xij * xij + yij * yij + zij * zij;
249 enuc += atoms[i].atomic_number * atoms[j].atomic_number / r;
251 cout <<
"Nuclear repulsion energy = " << std::setprecision(15) << enuc
256 cout <<
"libint2::Atomic Cartesian coordinates (a.u.):" << endl;
257 for (
const auto& a : atoms)
258 std::cout << a.atomic_number <<
" " << a.x <<
" " << a.y <<
" " << a.z
261 BasisSet obs(basisname, atoms);
262 cout <<
"orbital basis set rank = " << obs.nbf() << endl;
264#ifdef HAVE_DENSITY_FITTING
266 if (do_density_fitting) {
267 dfbs = BasisSet(dfbasisname, atoms);
268 cout <<
"density-fitting basis set rank = " << dfbs.nbf() << endl;
281 obs_shellpair_list = compute_shellpair_list(obs);
283 for (
auto& sp : obs_shellpair_list) {
284 nsp += sp.second.size();
286 std::cout <<
"# of {all,non-negligible} shell-pairs = {"
287 << obs.size() * (obs.size() + 1) / 2 <<
"," << nsp <<
"}"
292 auto S = compute_1body_ints<Operator::overlap>(obs)[0];
293 auto T = compute_1body_ints<Operator::kinetic>(obs)[0];
294 auto V = compute_1body_ints<Operator::nuclear>(
302 double XtX_condition_number;
311 double S_condition_number_threshold =
312 1.0 / std::numeric_limits<double>::epsilon();
313 std::tie(X, Xinv, XtX_condition_number) =
314 conditioning_orthogonalizer(S, S_condition_number_threshold);
320 const auto tstart = std::chrono::high_resolution_clock::now();
322 auto D_minbs = compute_soad(atoms);
323 BasisSet minbs(
"STO-3G", atoms);
329 std::cout <<
"projecting SOAD into AO basis ... ";
331 F += compute_2body_fock_general(
332 obs, D_minbs, minbs,
true ,
333 std::numeric_limits<double>::epsilon()
340 Eigen::SelfAdjointEigenSolver<Matrix> eig_solver(X.transpose() * F * X);
341 auto C = X * eig_solver.eigenvectors();
344 C_occ = C.leftCols(ndocc);
345 D = C_occ * C_occ.transpose();
347 const auto tstop = std::chrono::high_resolution_clock::now();
348 const std::chrono::duration<double> time_elapsed = tstop - tstart;
349 std::cout <<
"done (" << time_elapsed.count() <<
" s)" << std::endl;
354 auto K = compute_schwarz_ints<>(obs, BasisSet(),
true);
357#ifdef HAVE_DENSITY_FITTING
358 std::unique_ptr<DFFockEngine> dffockengine(
359 do_density_fitting ?
new DFFockEngine(obs, dfbs) : nullptr);
366 const auto maxiter = 100;
367 const auto conv = 1e-12;
369 auto rms_error = 1.0;
370 auto ediff_rel = 0.0;
372 auto n2 = D.cols() * D.rows();
378 bool reset_incremental_fock_formation =
false;
379 bool incremental_Fbuild_started =
false;
380 double start_incremental_F_threshold = 1e-5;
381 double next_reset_threshold = 0.0;
382 size_t last_reset_iteration = 0;
384 if (do_density_fitting) start_incremental_F_threshold = 0.0;
387 const auto tstart = std::chrono::high_resolution_clock::now();
394 if (not incremental_Fbuild_started &&
395 rms_error < start_incremental_F_threshold) {
396 incremental_Fbuild_started =
true;
397 reset_incremental_fock_formation =
false;
398 last_reset_iteration = iter - 1;
399 next_reset_threshold = rms_error / 1e1;
400 std::cout <<
"== started incremental fock build" << std::endl;
402 if (reset_incremental_fock_formation || not incremental_Fbuild_started) {
406 if (reset_incremental_fock_formation && incremental_Fbuild_started) {
407 reset_incremental_fock_formation =
false;
408 last_reset_iteration = iter;
409 next_reset_threshold = rms_error / 1e1;
410 std::cout <<
"== reset incremental fock build" << std::endl;
414 if (not do_density_fitting) {
416 const auto precision_F = std::min(
417 std::min(1e-3 / XtX_condition_number, 1e-7),
418 std::max(rms_error / 1e4, std::numeric_limits<double>::epsilon()));
419 F += compute_2body_fock(obs, D_diff, precision_F, K);
421#if HAVE_DENSITY_FITTING
423 F = H + dffockengine->compute_2body_fock_dfC(C_occ);
433 ehf = D.cwiseProduct(H + F).sum();
434 ediff_rel = std::abs((ehf - ehf_last) / ehf);
437 Matrix FD_comm = F * D * S - S * D * F;
438 rms_error = FD_comm.norm() / n2;
439 if (rms_error < next_reset_threshold || iter - last_reset_iteration >= 8)
440 reset_incremental_fock_formation =
true;
446 diis.extrapolate(F_diis, FD_comm);
451 Eigen::SelfAdjointEigenSolver<Matrix> eig_solver(X.transpose() * F_diis *
453 evals = eig_solver.eigenvalues();
454 auto C = X * eig_solver.eigenvectors();
457 C_occ = C.leftCols(ndocc);
458 D = C_occ * C_occ.transpose();
461 const auto tstop = std::chrono::high_resolution_clock::now();
462 const std::chrono::duration<double> time_elapsed = tstop - tstart;
465 std::cout <<
"\n\nIter E(HF) D(E)/E "
466 "RMS([F,D])/nn Time(s)\n";
467 printf(
" %02d %20.12f %20.12e %20.12e %10.5lf\n", iter, ehf + enuc,
468 ediff_rel, rms_error, time_elapsed.count());
470 }
while (((ediff_rel > conv) || (rms_error > conv)) && (iter < maxiter));
472 printf(
"** Hartree-Fock energy = %20.12f\n", ehf + enuc);
474 auto Mu = compute_1body_ints<Operator::emultipole2>(obs);
475 std::array<double, 3> mu;
476 std::array<double, 6> qu;
477 for (
int xyz = 0; xyz != 3; ++xyz)
478 mu[xyz] = -2 * D.cwiseProduct(Mu[xyz + 1])
480 for (
int k = 0; k != 6; ++k)
481 qu[k] = -2 * D.cwiseProduct(Mu[k + 4])
483 std::cout <<
"** edipole = ";
484 std::copy(mu.begin(), mu.end(),
485 std::ostream_iterator<double>(std::cout,
" "));
486 std::cout << std::endl;
487 std::cout <<
"** equadrupole = ";
488 std::copy(qu.begin(), qu.end(),
489 std::ostream_iterator<double>(std::cout,
" "));
490 std::cout << std::endl;
493#if LIBINT2_DERIV_ONEBODY_ORDER
495 Matrix F1 = Matrix::Zero(atoms.size(), 3);
496 Matrix F_Pulay = Matrix::Zero(atoms.size(), 3);
500 auto T1 = compute_1body_ints_deriv<Operator::kinetic>(1, obs, atoms);
501 auto V1 = compute_1body_ints_deriv<Operator::nuclear>(1, obs, atoms);
502 for (
auto atom = 0, i = 0; atom != atoms.size(); ++atom) {
503 for (
auto xyz = 0; xyz != 3; ++xyz, ++i) {
504 auto force = 2 * (T1[i] + V1[i]).cwiseProduct(D).sum();
505 F1(atom, xyz) += force;
513 DiagonalMatrix evals_occ(evals.topRows(ndocc));
514 Matrix W = C_occ * evals_occ * C_occ.transpose();
515 auto S1 = compute_1body_ints_deriv<Operator::overlap>(1, obs, atoms);
516 for (
auto atom = 0, i = 0; atom != atoms.size(); ++atom) {
517 for (
auto xyz = 0; xyz != 3; ++xyz, ++i) {
518 auto force = 2 * S1[i].cwiseProduct(W).sum();
519 F_Pulay(atom, xyz) -= force;
523 std::cout <<
"** 1-body forces = ";
524 for (
int atom = 0; atom != atoms.size(); ++atom)
525 for (
int xyz = 0; xyz != 3; ++xyz) std::cout << F1(atom, xyz) <<
" ";
526 std::cout << std::endl;
527 std::cout <<
"** Pulay forces = ";
528 for (
int atom = 0; atom != atoms.size(); ++atom)
529 for (
int xyz = 0; xyz != 3; ++xyz)
530 std::cout << F_Pulay(atom, xyz) <<
" ";
531 std::cout << std::endl;
534#if LIBINT2_DERIV_ERI_ORDER
536 Matrix F2 = Matrix::Zero(atoms.size(), 3);
541 auto G1 = compute_2body_fock_deriv<1>(obs, atoms, D);
542 for (
auto atom = 0, i = 0; atom != atoms.size(); ++atom) {
543 for (
auto xyz = 0; xyz != 3; ++xyz, ++i) {
545 auto force = G1[i].cwiseProduct(D).sum();
546 F2(atom, xyz) += force;
550 std::cout <<
"** 2-body forces = ";
551 for (
int atom = 0; atom != atoms.size(); ++atom)
552 for (
int xyz = 0; xyz != 3; ++xyz) std::cout << F2(atom, xyz) <<
" ";
553 std::cout << std::endl;
558#if LIBINT2_DERIV_ONEBODY_ORDER && LIBINT2_DERIV_ERI_ORDER
560 Matrix FN = Matrix::Zero(atoms.size(), 3);
564 for (
auto a1 = 1; a1 != atoms.size(); ++a1) {
565 const auto& atom1 = atoms[a1];
566 for (
auto a2 = 0; a2 < a1; ++a2) {
567 const auto& atom2 = atoms[a2];
569 auto x12 = atom1.x - atom2.x;
570 auto y12 = atom1.y - atom2.y;
571 auto z12 = atom1.z - atom2.z;
572 auto r12_2 = x12 * x12 + y12 * y12 + z12 * z12;
573 auto r12 = sqrt(r12_2);
574 auto r12_3 = r12 * r12_2;
576 auto z1z2_over_r12_3 =
577 atom1.atomic_number * atom2.atomic_number / r12_3;
579 auto fx = -x12 * z1z2_over_r12_3;
580 auto fy = -y12 * z1z2_over_r12_3;
581 auto fz = -z12 * z1z2_over_r12_3;
591 std::cout <<
"** nuclear repulsion forces = ";
592 for (
int atom = 0; atom != atoms.size(); ++atom)
593 for (
int xyz = 0; xyz != 3; ++xyz) std::cout << FN(atom, xyz) <<
" ";
594 std::cout << std::endl;
596 auto F = F1 + F_Pulay + F2 + FN;
597 std::cout <<
"** Hartree-Fock forces = ";
598 for (
int atom = 0; atom != atoms.size(); ++atom)
599 for (
int xyz = 0; xyz != 3; ++xyz) std::cout << F(atom, xyz) <<
" ";
600 std::cout << std::endl;
605 const auto ncoords = 3 * atoms.size();
607 const auto nelem = ncoords * (ncoords + 1) / 2;
608#if LIBINT2_DERIV_ONEBODY_ORDER > 1
610 Matrix H1 = Matrix::Zero(ncoords, ncoords);
611 Matrix H_Pulay = Matrix::Zero(ncoords, ncoords);
615 auto T2 = compute_1body_ints_deriv<Operator::kinetic>(2, obs, atoms);
616 auto V2 = compute_1body_ints_deriv<Operator::nuclear>(2, obs, atoms);
617 for (
auto row = 0, i = 0; row != ncoords; ++row) {
618 for (
auto col = row; col != ncoords; ++col, ++i) {
619 auto hess = 2 * (T2[i] + V2[i]).cwiseProduct(D).sum();
620 H1(row, col) += hess;
628 DiagonalMatrix evals_occ(evals.topRows(ndocc));
629 Matrix W = C_occ * evals_occ * C_occ.transpose();
630 auto S2 = compute_1body_ints_deriv<Operator::overlap>(2, obs, atoms);
631 for (
auto row = 0, i = 0; row != ncoords; ++row) {
632 for (
auto col = row; col != ncoords; ++col, ++i) {
633 auto hess = 2 * S2[i].cwiseProduct(W).sum();
634 H_Pulay(row, col) -= hess;
638 std::cout <<
"** 1-body hessian = ";
639 for (
auto row = 0, i = 0; row != ncoords; ++row) {
640 for (
auto col = row; col != ncoords; ++col) {
641 std::cout << H1(row, col) <<
" ";
644 std::cout << std::endl;
646 std::cout <<
"** Pulay hessian = ";
647 for (
auto row = 0, i = 0; row != ncoords; ++row) {
648 for (
auto col = row; col != ncoords; ++col) {
649 std::cout << H_Pulay(row, col) <<
" ";
652 std::cout << std::endl;
655#if LIBINT2_DERIV_ERI_ORDER > 1
657 Matrix H2 = Matrix::Zero(ncoords, ncoords);
662 auto G2 = compute_2body_fock_deriv<2>(obs, atoms, D);
663 for (
auto row = 0, i = 0; row != ncoords; ++row) {
664 for (
auto col = row; col != ncoords; ++col, ++i) {
666 auto hess = G2[i].cwiseProduct(D).sum();
667 H2(row, col) += hess;
671 std::cout <<
"** 2-body hessian = ";
672 for (
auto row = 0, i = 0; row != ncoords; ++row) {
673 for (
auto col = row; col != ncoords; ++col) {
674 std::cout << H2(row, col) <<
" ";
677 std::cout << std::endl;
682#if LIBINT2_DERIV_ONEBODY_ORDER > 1 && LIBINT2_DERIV_ERI_ORDER > 1
685 Matrix HN = Matrix::Zero(ncoords, ncoords);
690 for (
auto a1 = 1; a1 != atoms.size(); ++a1) {
691 const auto& atom1 = atoms[a1];
692 for (
auto a2 = 0; a2 < a1; ++a2) {
693 const auto& atom2 = atoms[a2];
695 auto x12 = atom1.x - atom2.x;
696 auto y12 = atom1.y - atom2.y;
697 auto z12 = atom1.z - atom2.z;
698 auto x12_2 = x12 * x12;
699 auto y12_2 = y12 * y12;
700 auto z12_2 = z12 * z12;
701 auto r12_2 = x12 * x12 + y12 * y12 + z12 * z12;
702 auto r12 = sqrt(r12_2);
703 auto r12_5 = r12 * r12_2 * r12_2;
705 auto z1z2_over_r12_5 =
706 atom1.atomic_number * atom2.atomic_number / r12_5;
708 HN(3 * a1 + 0, 3 * a1 + 0) += z1z2_over_r12_5 * (3 * x12_2 - r12_2);
709 HN(3 * a1 + 1, 3 * a1 + 1) += z1z2_over_r12_5 * (3 * y12_2 - r12_2);
710 HN(3 * a1 + 2, 3 * a1 + 2) += z1z2_over_r12_5 * (3 * z12_2 - r12_2);
711 HN(3 * a1 + 0, 3 * a1 + 1) += z1z2_over_r12_5 * (3 * x12 * y12);
712 HN(3 * a1 + 0, 3 * a1 + 2) += z1z2_over_r12_5 * (3 * x12 * z12);
713 HN(3 * a1 + 1, 3 * a1 + 2) += z1z2_over_r12_5 * (3 * y12 * z12);
715 HN(3 * a2 + 0, 3 * a2 + 0) += z1z2_over_r12_5 * (3 * x12_2 - r12_2);
716 HN(3 * a2 + 1, 3 * a2 + 1) += z1z2_over_r12_5 * (3 * y12_2 - r12_2);
717 HN(3 * a2 + 2, 3 * a2 + 2) += z1z2_over_r12_5 * (3 * z12_2 - r12_2);
718 HN(3 * a2 + 0, 3 * a2 + 1) += z1z2_over_r12_5 * (3 * x12 * y12);
719 HN(3 * a2 + 0, 3 * a2 + 2) += z1z2_over_r12_5 * (3 * x12 * z12);
720 HN(3 * a2 + 1, 3 * a2 + 2) += z1z2_over_r12_5 * (3 * y12 * z12);
722 HN(3 * a2 + 0, 3 * a1 + 0) -= z1z2_over_r12_5 * (3 * x12_2 - r12_2);
723 HN(3 * a2 + 1, 3 * a1 + 1) -= z1z2_over_r12_5 * (3 * y12_2 - r12_2);
724 HN(3 * a2 + 2, 3 * a1 + 2) -= z1z2_over_r12_5 * (3 * z12_2 - r12_2);
725 HN(3 * a2 + 1, 3 * a1 + 0) -= z1z2_over_r12_5 * (3 * y12 * x12);
726 HN(3 * a2 + 2, 3 * a1 + 0) -= z1z2_over_r12_5 * (3 * z12 * x12);
727 HN(3 * a2 + 2, 3 * a1 + 1) -= z1z2_over_r12_5 * (3 * z12 * y12);
728 HN(3 * a2 + 0, 3 * a1 + 1) -= z1z2_over_r12_5 * (3 * x12 * y12);
729 HN(3 * a2 + 0, 3 * a1 + 2) -= z1z2_over_r12_5 * (3 * x12 * z12);
730 HN(3 * a2 + 1, 3 * a1 + 2) -= z1z2_over_r12_5 * (3 * y12 * z12);
734 std::cout <<
"** nuclear repulsion hessian = ";
735 for (
auto row = 0, i = 0; row != ncoords; ++row) {
736 for (
auto col = row; col != ncoords; ++col) {
737 std::cout << HN(row, col) <<
" ";
740 std::cout << std::endl;
742 auto H = H1 + H_Pulay + H2 + HN;
743 std::cout <<
"** Hartree-Fock hessian = ";
744 for (
auto row = 0, i = 0; row != ncoords; ++row) {
745 for (
auto col = row; col != ncoords; ++col) {
746 std::cout << H(row, col) <<
" ";
749 std::cout << std::endl;
758 catch (
const char* ex) {
759 cerr <<
"caught exception: " << ex << endl;
761 }
catch (std::string& ex) {
762 cerr <<
"caught exception: " << ex << endl;
764 }
catch (std::exception& ex) {
765 cerr << ex.what() << endl;
768 cerr <<
"caught unknown exception\n";
775std::vector<libint2::Atom> read_geometry(
const std::string& filename) {
776 std::cout <<
"Will read geometry from " << filename << std::endl;
777 std::ifstream is(filename);
779 char errmsg[256] =
"Could not open file ";
780 strncpy(errmsg + 20, filename.c_str(), 235);
782 throw std::ios_base::failure(errmsg);
789 std::ostringstream oss;
794 std::istringstream iss(oss.str());
798 if (filename.rfind(
".xyz") != std::string::npos)
801 throw std::invalid_argument{
"only .xyz files are accepted"};
807Matrix compute_soad(
const std::vector<libint2::Atom>& atoms) {
810 for (
const auto& atom : atoms) {
811 const auto Z = atom.atomic_number;
816 Matrix D = Matrix::Zero(nao, nao);
817 size_t ao_offset = 0;
818 for (
const auto& atom : atoms) {
819 const auto Z = atom.atomic_number;
821 for (
const auto& occ : occvec) {
822 D(ao_offset, ao_offset) = occ;
830Matrix compute_shellblock_norm(
const BasisSet& obs,
const Matrix& A) {
831 const auto nsh = obs.size();
832 Matrix Ash(nsh, nsh);
834 auto shell2bf = obs.shell2bf();
835 for (
size_t s1 = 0; s1 != nsh; ++s1) {
836 const auto& s1_first = shell2bf[s1];
837 const auto& s1_size = obs[s1].size();
838 for (
size_t s2 = 0; s2 != nsh; ++s2) {
839 const auto& s2_first = shell2bf[s2];
840 const auto& s2_size = obs[s2].size();
842 Ash(s1, s2) = A.block(s1_first, s2_first, s1_size, s2_size)
843 .lpNorm<Eigen::Infinity>();
850template <Operator obtype,
typename OperatorParams>
851std::array<Matrix, libint2::operator_traits<obtype>::nopers> compute_1body_ints(
852 const BasisSet& obs, OperatorParams params) {
853 const auto n = obs.nbf();
854 const auto nshells = obs.size();
855 using libint2::nthreads;
856 typedef std::array<Matrix, libint2::operator_traits<obtype>::nopers>
858 const unsigned int nopers = libint2::operator_traits<obtype>::nopers;
860 for (
auto& r : result) r = Matrix::Zero(n, n);
863 std::vector<libint2::Engine> engines(nthreads);
864 engines[0] = libint2::Engine(obtype, obs.max_nprim(), obs.max_l(), 0);
868 if (obtype == Operator::nuclear || obtype == Operator::erf_nuclear ||
869 obtype == Operator::erfc_nuclear) {
870 engines[0].set_params(params);
872 for (
size_t i = 1; i != nthreads; ++i) {
873 engines[i] = engines[0];
876 auto shell2bf = obs.shell2bf();
878 auto compute = [&](
int thread_id) {
879 const auto& buf = engines[thread_id].results();
884 for (
auto s1 = 0l, s12 = 0l; s1 != nshells; ++s1) {
885 auto bf1 = shell2bf[s1];
886 auto n1 = obs[s1].size();
887 auto s1_offset = s1 * (s1 + 1) / 2;
889 for (
auto s2 : obs_shellpair_list[s1]) {
890 auto s12 = s1_offset + s2;
891 if (s12 % nthreads != thread_id)
continue;
893 auto bf2 = shell2bf[s2];
894 auto n2 = obs[s2].size();
899 engines[thread_id].compute(obs[s1], obs[s2]);
901 for (
unsigned int op = 0; op != nopers; ++op) {
904 Eigen::Map<const Matrix> buf_mat(buf[op], n1, n2);
905 result[op].block(bf1, bf2, n1, n2) = buf_mat;
908 result[op].block(bf2, bf1, n2, n1) = buf_mat.transpose();
919#if LIBINT2_DERIV_ONEBODY_ORDER
920template <Operator obtype>
921std::vector<Matrix> compute_1body_ints_deriv(
922 unsigned deriv_order,
const BasisSet& obs,
923 const std::vector<libint2::Atom>& atoms) {
924 using libint2::nthreads;
925 const auto n = obs.nbf();
926 const auto nshells = obs.size();
927 constexpr auto nopers = libint2::operator_traits<obtype>::nopers;
928 const auto nresults =
929 nopers * libint2::num_geometrical_derivatives(atoms.size(), deriv_order);
930 typedef std::vector<Matrix> result_type;
931 result_type result(nresults);
932 for (
auto& r : result) r = Matrix::Zero(n, n);
935 std::vector<libint2::Engine> engines(nthreads);
937 libint2::Engine(obtype, obs.max_nprim(), obs.max_l(), deriv_order);
941 if (obtype == Operator::nuclear) {
942 std::vector<std::pair<double, std::array<double, 3>>> q;
943 for (
const auto& atom : atoms) {
944 q.push_back({
static_cast<double>(atom.atomic_number),
945 {{atom.x, atom.y, atom.z}}});
947 engines[0].set_params(q);
949 for (
size_t i = 1; i != nthreads; ++i) {
950 engines[i] = engines[0];
953 auto shell2bf = obs.shell2bf();
954 auto shell2atom = obs.shell2atom(atoms);
956 const auto natoms = atoms.size();
957 const auto two_times_ncoords = 6 * natoms;
958 const auto nderivcenters_shset =
959 2 + ((obtype == Operator::nuclear) ? natoms : 0);
961 auto compute = [&](
int thread_id) {
962 const auto& buf = engines[thread_id].results();
967 for (
auto s1 = 0l, s12 = 0l; s1 != nshells; ++s1) {
968 auto bf1 = shell2bf[s1];
969 auto n1 = obs[s1].size();
970 auto atom1 = shell2atom[s1];
973 auto s1_offset = s1 * (s1 + 1) / 2;
974 for (
auto s2 : obs_shellpair_list[s1]) {
975 auto s12 = s1_offset + s2;
976 if (s12 % nthreads != thread_id)
continue;
978 auto bf2 = shell2bf[s2];
979 auto n2 = obs[s2].size();
980 auto atom2 = shell2atom[s2];
985 engines[thread_id].compute(obs[s1], obs[s2]);
989 auto add_shellset_to_dest = [&](std::size_t op, std::size_t idx,
990 double scale = 1.0) {
993 Eigen::Map<const Matrix> buf_mat(buf[idx], n1, n2);
995 result[op].block(bf1, bf2, n1, n2) += buf_mat;
997 result[op].block(bf1, bf2, n1, n2) += scale * buf_mat;
1001 result[op].block(bf2, bf1, n2, n1) += buf_mat.transpose();
1003 result[op].block(bf2, bf1, n2, n1) += scale * buf_mat.transpose();
1007 switch (deriv_order) {
1009 for (std::size_t op = 0; op != nopers; ++op) {
1010 add_shellset_to_dest(op, op);
1028 std::size_t shellset_idx = 0;
1029 for (
auto c = 0; c != nderivcenters_shset; ++c) {
1030 auto atom = (c == 0) ? atom1 : ((c == 1) ? atom2 : c - 2);
1031 auto op_start = 3 * atom * nopers;
1032 auto op_fence = op_start + nopers;
1033 for (
auto xyz = 0; xyz != 3;
1034 ++xyz, op_start += nopers, op_fence += nopers) {
1035 for (
unsigned int op = op_start; op != op_fence;
1036 ++op, ++shellset_idx) {
1037 add_shellset_to_dest(op, shellset_idx);
1067 auto upper_triangle_index = [](
size_t n2,
long i,
long j) {
1068 return std::min(i, j) * (n2 - std::min(i, j) - 1) / 2 +
1073 std::size_t shellset_idx = 0;
1074 for (
auto c1 = 0; c1 != nderivcenters_shset; ++c1) {
1075 auto a1 = (c1 == 0) ? atom1 : ((c1 == 1) ? atom2 : c1 - 2);
1076 auto coord1 = 3 * a1;
1077 for (
auto xyz1 = 0; xyz1 != 3; ++xyz1, ++coord1) {
1078 for (
auto c2 = c1; c2 != nderivcenters_shset; ++c2) {
1079 auto a2 = (c2 == 0) ? atom1 : ((c2 == 1) ? atom2 : c2 - 2);
1080 auto xyz2_start = (c1 == c2) ? xyz1 : 0;
1081 auto coord2 = 3 * a2 + xyz2_start;
1082 for (
auto xyz2 = xyz2_start; xyz2 != 3; ++xyz2, ++coord2) {
1083 double scale = (coord1 == coord2 && c1 != c2) ? 2.0 : 1.0;
1085 const auto coord12 =
1086 upper_triangle_index(two_times_ncoords, coord1, coord2);
1087 auto op_start = coord12 * nopers;
1088 auto op_fence = op_start + nopers;
1089 for (
auto op = op_start; op != op_fence;
1090 ++op, ++shellset_idx) {
1091 add_shellset_to_dest(op, shellset_idx, scale);
1100 assert(
false &&
"not yet implemented");
1102 using ShellSetDerivIterator =
1104 std::vector<unsigned int>>;
1105 ShellSetDerivIterator shellset_diter(deriv_order,
1106 nderivcenters_shset);
1107 while (shellset_diter) {
1108 const auto& deriv = *shellset_diter;
1123template <lib
int2::Operator Kernel>
1124Matrix compute_schwarz_ints(
1125 const BasisSet& bs1,
const BasisSet& _bs2,
bool use_2norm,
1126 typename libint2::operator_traits<Kernel>::oper_params_type params) {
1127 const BasisSet& bs2 = (_bs2.empty() ? bs1 : _bs2);
1128 const auto nsh1 = bs1.size();
1129 const auto nsh2 = bs2.size();
1130 const auto bs1_equiv_bs2 = (&bs1 == &bs2);
1132 Matrix K = Matrix::Zero(nsh1, nsh2);
1135 using libint2::Engine;
1136 using libint2::nthreads;
1137 std::vector<Engine> engines(nthreads);
1141 engines[0] = Engine(Kernel, std::max(bs1.max_nprim(), bs2.max_nprim()),
1142 std::max(bs1.max_l(), bs2.max_l()), 0, epsilon, params);
1143 for (
size_t i = 1; i != nthreads; ++i) {
1144 engines[i] = engines[0];
1147 std::cout <<
"computing Schwarz bound prerequisites (kernel=" << (int)Kernel
1154 auto compute = [&](
int thread_id) {
1155 const auto& buf = engines[thread_id].results();
1158 for (
auto s1 = 0l, s12 = 0l; s1 != nsh1; ++s1) {
1159 auto n1 = bs1[s1].size();
1161 auto s2_max = bs1_equiv_bs2 ? s1 : nsh2 - 1;
1162 for (
auto s2 = 0; s2 <= s2_max; ++s2, ++s12) {
1163 if (s12 % nthreads != thread_id)
continue;
1165 auto n2 = bs2[s2].size();
1168 engines[thread_id].compute2<Kernel, BraKet::xx_xx, 0>(bs1[s1], bs2[s2],
1170 assert(buf[0] !=
nullptr &&
1171 "to compute Schwarz ints turn off primitive screening");
1174 Eigen::Map<const Matrix> buf_mat(buf[0], n12, n12);
1176 use_2norm ? buf_mat.norm() : buf_mat.lpNorm<Eigen::Infinity>();
1177 K(s1, s2) = std::sqrt(norm2);
1178 if (bs1_equiv_bs2) K(s2, s1) = K(s1, s2);
1186 std::cout <<
"done (" << timer.read(0) <<
" s)" << std::endl;
1191Matrix compute_do_ints(
const BasisSet& bs1,
const BasisSet& bs2,
1193 return compute_schwarz_ints<libint2::Operator::delta>(bs1, bs2, use_2norm);
1196shellpair_list_t compute_shellpair_list(
const BasisSet& bs1,
1197 const BasisSet& _bs2,
1198 const double threshold) {
1199 const BasisSet& bs2 = (_bs2.empty() ? bs1 : _bs2);
1200 const auto nsh1 = bs1.size();
1201 const auto nsh2 = bs2.size();
1202 const auto bs1_equiv_bs2 = (&bs1 == &bs2);
1204 using libint2::nthreads;
1207 using libint2::Engine;
1208 std::vector<Engine> engines;
1209 engines.reserve(nthreads);
1210 engines.emplace_back(Operator::overlap,
1211 std::max(bs1.max_nprim(), bs2.max_nprim()),
1212 std::max(bs1.max_l(), bs2.max_l()), 0);
1213 for (
size_t i = 1; i != nthreads; ++i) {
1214 engines.push_back(engines[0]);
1217 std::cout <<
"computing non-negligible shell-pair list ... ";
1223 shellpair_list_t result;
1227 auto compute = [&](
int thread_id) {
1228 auto& engine = engines[thread_id];
1229 const auto& buf = engine.results();
1232 for (
auto s1 = 0l, s12 = 0l; s1 != nsh1; ++s1) {
1234 if (result.find(s1) == result.end())
1235 result.insert(std::make_pair(s1, std::vector<size_t>()));
1238 auto n1 = bs1[s1].size();
1240 auto s2_max = bs1_equiv_bs2 ? s1 : nsh2 - 1;
1241 for (
auto s2 = 0; s2 <= s2_max; ++s2, ++s12) {
1242 if (s12 % nthreads != thread_id)
continue;
1244 auto on_same_center = (bs1[s1].O == bs2[s2].O);
1245 bool significant = on_same_center;
1246 if (not on_same_center) {
1247 auto n2 = bs2[s2].size();
1248 engines[thread_id].compute(bs1[s1], bs2[s2]);
1249 Eigen::Map<const Matrix> buf_mat(buf[0], n1, n2);
1250 auto norm = buf_mat.norm();
1251 significant = (norm >= threshold);
1256 result[s1].emplace_back(s2);
1267 auto sort = [&](
int thread_id) {
1268 for (
auto s1 = 0l; s1 != nsh1; ++s1) {
1269 if (s1 % nthreads == thread_id) {
1270 auto& list = result[s1];
1271 std::sort(list.begin(), list.end());
1279 std::cout <<
"done (" << timer.
read(0) <<
" s)" << std::endl;
1296std::tuple<Matrix, Matrix, size_t, double, double> gensqrtinv(
1297 const Matrix& S,
bool symmetric =
false,
1298 double max_condition_number = 1e8) {
1299 Eigen::SelfAdjointEigenSolver<Matrix> eig_solver(S);
1300 auto U = eig_solver.eigenvectors();
1301 auto s = eig_solver.eigenvalues();
1302 auto s_max = s.maxCoeff();
1303 auto condition_number = std::min(
1304 s_max / std::max(s.minCoeff(), std::numeric_limits<double>::min()),
1305 1.0 / std::numeric_limits<double>::epsilon());
1306 auto threshold = s_max / max_condition_number;
1309 for (
long i = n - 1; i >= 0; --i) {
1310 if (s(i) >= threshold) {
1316 auto sigma = s.bottomRows(n_cond);
1317 auto result_condition_number = sigma.maxCoeff() / sigma.minCoeff();
1318 auto sigma_sqrt = sigma.array().sqrt().matrix().asDiagonal();
1319 auto sigma_invsqrt = sigma.array().sqrt().inverse().matrix().asDiagonal();
1322 auto U_cond = U.block(0, n - n_cond, n, n_cond);
1323 Matrix X = U_cond * sigma_invsqrt;
1324 Matrix Xinv = U_cond * sigma_sqrt;
1327 X = X * U_cond.transpose();
1328 Xinv = Xinv * U_cond.transpose();
1330 return std::make_tuple(X, Xinv,
size_t(n_cond), condition_number,
1331 result_condition_number);
1334std::tuple<Matrix, Matrix, double> conditioning_orthogonalizer(
1335 const Matrix& S,
double S_condition_number_threshold) {
1337 double S_condition_number;
1338 double XtX_condition_number;
1341 assert(S.rows() == S.cols());
1343 std::tie(X, Xinv, obs_rank, S_condition_number, XtX_condition_number) =
1344 gensqrtinv(S,
false, S_condition_number_threshold);
1345 auto obs_nbf_omitted = (long)S.rows() - (long)obs_rank;
1346 std::cout <<
"overlap condition number = " << S_condition_number;
1347 if (obs_nbf_omitted > 0)
1348 std::cout <<
" (dropped " << obs_nbf_omitted <<
" "
1349 << (obs_nbf_omitted > 1 ?
"fns" :
"fn") <<
" to reduce to "
1350 << XtX_condition_number <<
")";
1351 std::cout << std::endl;
1353 if (obs_nbf_omitted > 0) {
1354 Matrix should_be_I = X.transpose() * S * X;
1355 Matrix I = Matrix::Identity(should_be_I.rows(), should_be_I.cols());
1356 std::cout <<
"||X^t * S * X - I||_2 = " << (should_be_I - I).norm()
1357 <<
" (should be 0)" << std::endl;
1360 return std::make_tuple(X, Xinv, XtX_condition_number);
1363Matrix compute_2body_2index_ints(
const BasisSet& bs) {
1364 using libint2::nthreads;
1365 const auto n = bs.nbf();
1366 const auto nshells = bs.size();
1367 Matrix result = Matrix::Zero(n, n);
1370 using libint2::Engine;
1371 std::vector<Engine> engines(nthreads);
1373 Engine(libint2::Operator::coulomb, bs.max_nprim(), bs.max_l(), 0);
1374 engines[0].set(BraKet::xs_xs);
1375 for (
size_t i = 1; i != nthreads; ++i) {
1376 engines[i] = engines[0];
1379 auto shell2bf = bs.shell2bf();
1382 auto compute = [&](
int thread_id) {
1383 auto& engine = engines[thread_id];
1384 const auto& buf = engine.results();
1389 for (
auto s1 = 0l, s12 = 0l; s1 != nshells; ++s1) {
1390 auto bf1 = shell2bf[s1];
1391 auto n1 = bs[s1].size();
1393 for (
auto s2 = 0; s2 <= s1; ++s2, ++s12) {
1394 if (s12 % nthreads != thread_id)
continue;
1396 auto bf2 = shell2bf[s2];
1397 auto n2 = bs[s2].size();
1400 engine.compute(bs[s1], bs[s2]);
1401 if (buf[0] ==
nullptr)
1406 Eigen::Map<const Matrix> buf_mat(buf[0], n1, n2);
1407 result.block(bf1, bf2, n1, n2) = buf_mat;
1410 result.block(bf2, bf1, n2, n1) = buf_mat.transpose();
1420Matrix compute_2body_fock(
const BasisSet& obs,
const Matrix& D,
1421 double precision,
const Matrix& Schwarz) {
1422 const auto n = obs.nbf();
1423 const auto nshells = obs.size();
1424 using libint2::nthreads;
1425 std::vector<Matrix> G(nthreads, Matrix::Zero(n, n));
1427 const auto do_schwarz_screen =
Schwarz.cols() != 0 &&
Schwarz.rows() != 0;
1428 Matrix D_shblk_norm =
1429 compute_shellblock_norm(obs, D);
1431 auto fock_precision = precision;
1434 auto max_nprim = obs.max_nprim();
1435 auto max_nprim4 = max_nprim * max_nprim * max_nprim * max_nprim;
1436 auto engine_precision = std::min(fock_precision / D_shblk_norm.maxCoeff(),
1437 std::numeric_limits<double>::epsilon()) /
1441 using libint2::Engine;
1442 std::vector<Engine> engines(nthreads);
1443 engines[0] = Engine(Operator::coulomb, obs.max_nprim(), obs.max_l(), 0);
1444 engines[0].set_precision(engine_precision);
1448 std::cout <<
"compute_2body_fock:precision = " << precision << std::endl;
1449 std::cout <<
"Engine::precision = " << engines[0].precision() << std::endl;
1450 for (
size_t i = 1; i != nthreads; ++i) {
1451 engines[i] = engines[0];
1453 std::atomic<size_t> num_ints_computed{0};
1455#if defined(REPORT_INTEGRAL_TIMINGS)
1456 std::vector<libint2::Timers<1>> timers(nthreads);
1459 auto shell2bf = obs.shell2bf();
1461 auto lambda = [&](
int thread_id) {
1462 auto& engine = engines[thread_id];
1463 auto& g = G[thread_id];
1464 const auto& buf = engine.results();
1466#if defined(REPORT_INTEGRAL_TIMINGS)
1467 auto& timer = timers[thread_id];
1473 for (
auto s1 = 0l, s1234 = 0l; s1 != nshells; ++s1) {
1474 auto bf1_first = shell2bf[s1];
1475 auto n1 = obs[s1].size();
1477 for (
const auto& s2 : obs_shellpair_list[s1]) {
1478 auto bf2_first = shell2bf[s2];
1479 auto n2 = obs[s2].size();
1481 const auto Dnorm12 = do_schwarz_screen ? D_shblk_norm(s1, s2) : 0.;
1483 for (
auto s3 = 0; s3 <= s1; ++s3) {
1484 auto bf3_first = shell2bf[s3];
1485 auto n3 = obs[s3].size();
1487 const auto Dnorm123 =
1489 ? std::max(D_shblk_norm(s1, s3),
1490 std::max(D_shblk_norm(s2, s3), Dnorm12))
1493 const auto s4_max = (s1 == s3) ? s2 : s3;
1494 for (
const auto& s4 : obs_shellpair_list[s3]) {
1499 if ((s1234++) % nthreads != thread_id)
continue;
1501 const auto Dnorm1234 =
1504 D_shblk_norm(s1, s4),
1505 std::max(D_shblk_norm(s2, s4),
1506 std::max(D_shblk_norm(s3, s4), Dnorm123)))
1509 if (do_schwarz_screen &&
1513 auto bf4_first = shell2bf[s4];
1514 auto n4 = obs[s4].size();
1516 num_ints_computed += n1 * n2 * n3 * n4;
1520 auto s12_deg = (s1 == s2) ? 1.0 : 2.0;
1521 auto s34_deg = (s3 == s4) ? 1.0 : 2.0;
1522 auto s12_34_deg = (s1 == s3) ? (s2 == s4 ? 1.0 : 2.0) : 2.0;
1523 auto s1234_deg = s12_deg * s34_deg * s12_34_deg;
1525#if defined(REPORT_INTEGRAL_TIMINGS)
1529 engine.compute2<Operator::coulomb, BraKet::xx_xx, 0>(
1530 obs[s1], obs[s2], obs[s3], obs[s4]);
1531 const auto* buf_1234 = buf[0];
1532 if (buf_1234 ==
nullptr)
1535#if defined(REPORT_INTEGRAL_TIMINGS)
1539 Eigen::Map<MatrixXd> buf_1234_map(buf_1234, n12, n34);
1554 for (
auto f1 = 0, f1234 = 0; f1 != n1; ++f1) {
1555 const auto bf1 = f1 + bf1_first;
1556 for (
auto f2 = 0; f2 != n2; ++f2) {
1557 const auto bf2 = f2 + bf2_first;
1558 for (
auto f3 = 0; f3 != n3; ++f3) {
1559 const auto bf3 = f3 + bf3_first;
1560 for (
auto f4 = 0; f4 != n4; ++f4, ++f1234) {
1561 const auto bf4 = f4 + bf4_first;
1563 const auto value = buf_1234[f1234];
1565 const auto value_scal_by_deg = value * s1234_deg;
1567 g(bf1, bf2) += D(bf3, bf4) * value_scal_by_deg;
1568 g(bf3, bf4) += D(bf1, bf2) * value_scal_by_deg;
1569 g(bf1, bf3) -= 0.25 * D(bf2, bf4) * value_scal_by_deg;
1570 g(bf2, bf4) -= 0.25 * D(bf1, bf3) * value_scal_by_deg;
1571 g(bf1, bf4) -= 0.25 * D(bf2, bf3) * value_scal_by_deg;
1572 g(bf2, bf3) -= 0.25 * D(bf1, bf4) * value_scal_by_deg;
1586 for (
size_t i = 1; i != nthreads; ++i) {
1590#if defined(REPORT_INTEGRAL_TIMINGS)
1591 double time_for_ints = 0.0;
1592 for (
auto& t : timers) {
1593 time_for_ints += t.read(0);
1595 std::cout <<
"time for integrals = " << time_for_ints << std::endl;
1596 for (
int t = 0; t != nthreads; ++t) engines[t].print_timers();
1599 Matrix GG = 0.5 * (G[0] + G[0].transpose());
1601 std::cout <<
"# of integrals = " << num_ints_computed << std::endl;
1607#if LIBINT2_DERIV_ERI_ORDER
1608template <
unsigned deriv_order>
1609std::vector<Matrix> compute_2body_fock_deriv(
1610 const BasisSet& obs,
const std::vector<libint2::Atom>& atoms,
1611 const Matrix& D,
double precision,
const Matrix& Schwarz) {
1612 const auto n = obs.nbf();
1613 const auto nshells = obs.size();
1614 const auto nderiv_shellset = libint2::num_geometrical_derivatives(
1616 const auto nderiv = libint2::num_geometrical_derivatives(
1617 atoms.size(), deriv_order);
1618 const auto ncoords_times_two = (atoms.size() * 3) * 2;
1619 using libint2::nthreads;
1620 std::vector<Matrix> G(nthreads * nderiv, Matrix::Zero(n, n));
1622 const auto do_schwarz_screen =
Schwarz.cols() != 0 &&
Schwarz.rows() != 0;
1623 Matrix D_shblk_norm =
1624 compute_shellblock_norm(obs, D);
1626 auto fock_precision = precision;
1629 auto max_nprim = obs.max_nprim();
1630 auto max_nprim4 = max_nprim * max_nprim * max_nprim * max_nprim;
1631 auto engine_precision = std::min(fock_precision / D_shblk_norm.maxCoeff(),
1632 std::numeric_limits<double>::epsilon()) /
1636 using libint2::Engine;
1637 std::vector<Engine> engines(nthreads);
1639 Engine(Operator::coulomb, obs.max_nprim(), obs.max_l(), deriv_order);
1640 engines[0].set_precision(engine_precision);
1644 std::cout <<
"compute_2body_fock:precision = " << precision << std::endl;
1645 std::cout <<
"Engine::precision = " << engines[0].precision() << std::endl;
1646 for (
size_t i = 1; i != nthreads; ++i) {
1647 engines[i] = engines[0];
1649 std::atomic<size_t> num_ints_computed{0};
1651#if defined(REPORT_INTEGRAL_TIMINGS)
1652 std::vector<libint2::Timers<1>> timers(nthreads);
1655 auto shell2bf = obs.shell2bf();
1656 auto shell2atom = obs.shell2atom(atoms);
1658 auto lambda = [&](
int thread_id) {
1659 auto& engine = engines[thread_id];
1660 const auto& buf = engine.results();
1662#if defined(REPORT_INTEGRAL_TIMINGS)
1663 auto& timer = timers[thread_id];
1668 size_t shell_atoms[4];
1671 for (
auto s1 = 0l, s1234 = 0l; s1 != nshells; ++s1) {
1672 auto bf1_first = shell2bf[s1];
1673 auto n1 = obs[s1].size();
1674 shell_atoms[0] = shell2atom[s1];
1676 for (
const auto& s2 : obs_shellpair_list[s1]) {
1677 auto bf2_first = shell2bf[s2];
1678 auto n2 = obs[s2].size();
1679 shell_atoms[1] = shell2atom[s2];
1681 const auto Dnorm12 = do_schwarz_screen ? D_shblk_norm(s1, s2) : 0.;
1683 for (
auto s3 = 0; s3 <= s1; ++s3) {
1684 auto bf3_first = shell2bf[s3];
1685 auto n3 = obs[s3].size();
1686 shell_atoms[2] = shell2atom[s3];
1688 const auto Dnorm123 =
1690 ? std::max(D_shblk_norm(s1, s3),
1691 std::max(D_shblk_norm(s2, s3), Dnorm12))
1694 const auto s4_max = (s1 == s3) ? s2 : s3;
1695 for (
const auto& s4 : obs_shellpair_list[s3]) {
1700 if ((s1234++) % nthreads != thread_id)
continue;
1702 const auto Dnorm1234 =
1705 D_shblk_norm(s1, s4),
1706 std::max(D_shblk_norm(s2, s4),
1707 std::max(D_shblk_norm(s3, s4), Dnorm123)))
1710 if (do_schwarz_screen &&
1714 auto bf4_first = shell2bf[s4];
1715 auto n4 = obs[s4].size();
1716 shell_atoms[3] = shell2atom[s4];
1718 const auto n1234 = n1 * n2 * n3 * n4;
1722 auto s12_deg = (s1 == s2) ? 1.0 : 2.0;
1723 auto s34_deg = (s3 == s4) ? 1.0 : 2.0;
1724 auto s12_34_deg = (s1 == s3) ? (s2 == s4 ? 1.0 : 2.0) : 2.0;
1725 auto s1234_deg = s12_deg * s34_deg * s12_34_deg;
1729 auto add_shellset_to_dest = [&](std::size_t op, std::size_t idx,
1730 int coord1,
int coord2,
1731 double scale = 1.0) {
1733 auto shset = buf[idx];
1734 const auto weight = scale * s1234_deg;
1736 for (
auto f1 = 0, f1234 = 0; f1 != n1; ++f1) {
1737 const auto bf1 = f1 + bf1_first;
1738 for (
auto f2 = 0; f2 != n2; ++f2) {
1739 const auto bf2 = f2 + bf2_first;
1740 for (
auto f3 = 0; f3 != n3; ++f3) {
1741 const auto bf3 = f3 + bf3_first;
1742 for (
auto f4 = 0; f4 != n4; ++f4, ++f1234) {
1743 const auto bf4 = f4 + bf4_first;
1745 const auto value = shset[f1234];
1746 const auto wvalue = value * weight;
1748 g(bf1, bf2) += D(bf3, bf4) * wvalue;
1749 g(bf3, bf4) += D(bf1, bf2) * wvalue;
1750 g(bf1, bf3) -= 0.25 * D(bf2, bf4) * wvalue;
1751 g(bf2, bf4) -= 0.25 * D(bf1, bf3) * wvalue;
1752 g(bf1, bf4) -= 0.25 * D(bf2, bf3) * wvalue;
1753 g(bf2, bf3) -= 0.25 * D(bf1, bf4) * wvalue;
1760#if defined(REPORT_INTEGRAL_TIMINGS)
1764 engine.compute2<Operator::coulomb, BraKet::xx_xx, deriv_order>(
1765 obs[s1], obs[s2], obs[s3], obs[s4]);
1766 if (buf[0] ==
nullptr)
1768 num_ints_computed += nderiv_shellset * n1234;
1770#if defined(REPORT_INTEGRAL_TIMINGS)
1774 switch (deriv_order) {
1776 int coord1 = 0, coord2 = 0;
1777 add_shellset_to_dest(thread_id, 0, coord1, coord2);
1781 for (
auto d = 0; d != 12; ++d) {
1782 const int a = d / 3;
1783 const int xyz = d % 3;
1785 auto coord = shell_atoms[a] * 3 + xyz;
1786 auto& g = G[thread_id * nderiv + coord];
1788 int coord1 = 0, coord2 = 0;
1790 add_shellset_to_dest(thread_id * nderiv + coord, d, coord1,
1800 auto upper_triangle_index = [](
size_t n2,
size_t i,
size_t j) {
1801 return std::min(i, j) * (n2 - std::min(i, j) - 1) / 2 +
1805 std::size_t shellset_idx = 0;
1806 for (
auto c1 = 0; c1 != 4; ++c1) {
1807 auto a1 = shell_atoms[c1];
1808 auto coord1 = 3 * a1;
1809 for (
auto xyz1 = 0; xyz1 != 3; ++xyz1, ++coord1) {
1810 for (
auto c2 = c1; c2 != 4; ++c2) {
1811 auto a2 = shell_atoms[c2];
1812 auto xyz2_start = (c1 == c2) ? xyz1 : 0;
1813 auto coord2 = 3 * a2 + xyz2_start;
1814 for (
auto xyz2 = xyz2_start; xyz2 != 3;
1817 (coord1 == coord2 && c1 != c2) ? 2.0 : 1.0;
1819 const auto coord12 = upper_triangle_index(
1820 ncoords_times_two, coord1, coord2);
1821 auto op = thread_id * nderiv + coord12;
1822 add_shellset_to_dest(op, shellset_idx, coord1, coord2,
1832 assert(deriv_order <= 2 &&
1833 "support for 3rd and higher derivatives of the Fock "
1834 "matrix not yet implemented");
1845 for (
size_t t = 1; t != nthreads; ++t) {
1846 for (
auto d = 0; d != nderiv; ++d) {
1847 G[d] += G[t * nderiv + d];
1851#if defined(REPORT_INTEGRAL_TIMINGS)
1852 double time_for_ints = 0.0;
1853 for (
auto& t : timers) {
1854 time_for_ints += t.read(0);
1856 std::cout <<
"time for integrals = " << time_for_ints << std::endl;
1857 for (
int t = 0; t != nthreads; ++t) engines[t].print_timers();
1860 std::vector<Matrix> GG(nderiv);
1861 for (
auto d = 0; d != nderiv; ++d) {
1862 GG[d] = 0.5 * (G[d] + G[d].transpose());
1865 std::cout <<
"# of integrals = " << num_ints_computed << std::endl;
1873Matrix compute_2body_fock_general(
const BasisSet& obs,
const Matrix& D,
1874 const BasisSet& D_bs,
bool D_is_shelldiagonal,
1876 const auto n = obs.nbf();
1877 const auto nshells = obs.size();
1878 const auto n_D = D_bs.nbf();
1879 assert(D.cols() == D.rows() && D.cols() == n_D);
1881 using libint2::nthreads;
1882 std::vector<Matrix> G(nthreads, Matrix::Zero(n, n));
1885 using libint2::Engine;
1886 std::vector<Engine> engines(nthreads);
1887 engines[0] = Engine(libint2::Operator::coulomb,
1888 std::max(obs.max_nprim(), D_bs.max_nprim()),
1889 std::max(obs.max_l(), D_bs.max_l()), 0);
1890 engines[0].set_precision(precision);
1894 for (
size_t i = 1; i != nthreads; ++i) {
1895 engines[i] = engines[0];
1897 auto shell2bf = obs.shell2bf();
1898 auto shell2bf_D = D_bs.shell2bf();
1900 auto lambda = [&](
int thread_id) {
1901 auto& engine = engines[thread_id];
1902 auto& g = G[thread_id];
1903 const auto& buf = engine.results();
1906 for (
auto s1 = 0l, s1234 = 0l; s1 != nshells; ++s1) {
1907 auto bf1_first = shell2bf[s1];
1908 auto n1 = obs[s1].size();
1910 for (
auto s2 = 0; s2 <= s1; ++s2) {
1911 auto bf2_first = shell2bf[s2];
1912 auto n2 = obs[s2].size();
1914 for (
auto s3 = 0; s3 < D_bs.size(); ++s3) {
1915 auto bf3_first = shell2bf_D[s3];
1916 auto n3 = D_bs[s3].size();
1918 auto s4_begin = D_is_shelldiagonal ? s3 : 0;
1919 auto s4_fence = D_is_shelldiagonal ? s3 + 1 : D_bs.size();
1921 for (
auto s4 = s4_begin; s4 != s4_fence; ++s4, ++s1234) {
1922 if (s1234 % nthreads != thread_id)
continue;
1924 auto bf4_first = shell2bf_D[s4];
1925 auto n4 = D_bs[s4].size();
1929 auto s12_deg = (s1 == s2) ? 1.0 : 2.0;
1932 auto s34_deg = (s3 == s4) ? 1.0 : 2.0;
1933 auto s1234_deg = s12_deg * s34_deg;
1935 engine.compute2<Operator::coulomb, BraKet::xx_xx, 0>(
1936 obs[s1], obs[s2], D_bs[s3], D_bs[s4]);
1937 const auto* buf_1234 = buf[0];
1938 if (buf_1234 !=
nullptr) {
1939 for (
auto f1 = 0, f1234 = 0; f1 != n1; ++f1) {
1940 const auto bf1 = f1 + bf1_first;
1941 for (
auto f2 = 0; f2 != n2; ++f2) {
1942 const auto bf2 = f2 + bf2_first;
1943 for (
auto f3 = 0; f3 != n3; ++f3) {
1944 const auto bf3 = f3 + bf3_first;
1945 for (
auto f4 = 0; f4 != n4; ++f4, ++f1234) {
1946 const auto bf4 = f4 + bf4_first;
1948 const auto value = buf_1234[f1234];
1949 const auto value_scal_by_deg = value * s1234_deg;
1950 g(bf1, bf2) += 2.0 * D(bf3, bf4) * value_scal_by_deg;
1958 engine.compute2<Operator::coulomb, BraKet::xx_xx, 0>(
1959 obs[s1], D_bs[s3], obs[s2], D_bs[s4]);
1960 const auto* buf_1324 = buf[0];
1961 if (buf_1324 ==
nullptr)
1964 for (
auto f1 = 0, f1324 = 0; f1 != n1; ++f1) {
1965 const auto bf1 = f1 + bf1_first;
1966 for (
auto f3 = 0; f3 != n3; ++f3) {
1967 const auto bf3 = f3 + bf3_first;
1968 for (
auto f2 = 0; f2 != n2; ++f2) {
1969 const auto bf2 = f2 + bf2_first;
1970 for (
auto f4 = 0; f4 != n4; ++f4, ++f1324) {
1971 const auto bf4 = f4 + bf4_first;
1973 const auto value = buf_1324[f1324];
1974 const auto value_scal_by_deg = value * s12_deg;
1975 g(bf1, bf2) -= D(bf3, bf4) * value_scal_by_deg;
1989 for (
size_t i = 1; i != nthreads; ++i) {
1994 return 0.5 * (G[0] + G[0].transpose());
1997#ifdef HAVE_DENSITY_FITTING
1999Matrix DFFockEngine::compute_2body_fock_dfC(
const Matrix& Cocc) {
2000 using libint2::nthreads;
2002 const auto n = obs.nbf();
2003 const auto ndf = dfbs.nbf();
2007 std::vector<libint2::Timers<5>> timers(nthreads);
2008 for (
auto& timer : timers) timer.set_now_overhead(25);
2010 typedef btas::RangeNd<CblasRowMajor, std::array<long, 1>> Range1d;
2011 typedef btas::RangeNd<CblasRowMajor, std::array<long, 2>> Range2d;
2012 typedef btas::Tensor<double, Range1d> Tensor1d;
2013 typedef btas::Tensor<double, Range2d> Tensor2d;
2017 if (xyK.size() == 0) {
2018 wall_timer.
start(0);
2020 const auto nshells = obs.size();
2021 const auto nshells_df = dfbs.size();
2027 std::vector<libint2::Engine> engines(nthreads);
2028 engines[0] = libint2::Engine(libint2::Operator::coulomb,
2029 std::max(obs.max_nprim(), dfbs.max_nprim()),
2030 std::max(obs.max_l(), dfbs.max_l()), 0);
2031 engines[0].set(BraKet::xs_xx);
2032 for (
size_t i = 1; i != nthreads; ++i) {
2033 engines[i] = engines[0];
2036 auto shell2bf = obs.shell2bf();
2037 auto shell2bf_df = dfbs.shell2bf();
2039 Tensor3d Zxy{ndf, n, n};
2041 auto lambda = [&](
int thread_id) {
2042 auto& engine = engines[thread_id];
2043 auto& timer = timers[thread_id];
2044 const auto& results = engine.results();
2047 for (
auto s1 = 0l, s123 = 0l; s1 != nshells_df; ++s1) {
2048 auto bf1_first = shell2bf_df[s1];
2049 auto n1 = dfbs[s1].size();
2051 for (
auto s2 = 0; s2 != nshells; ++s2) {
2052 auto bf2_first = shell2bf[s2];
2053 auto n2 = obs[s2].size();
2054 const auto n12 = n1 * n2;
2056 for (
auto s3 = 0; s3 != nshells; ++s3, ++s123) {
2057 if (s123 % nthreads != thread_id)
continue;
2059 auto bf3_first = shell2bf[s3];
2060 auto n3 = obs[s3].size();
2061 const auto n123 = n12 * n3;
2065 engine.compute2<Operator::coulomb, BraKet::xs_xx, 0>(
2066 dfbs[s1], unitshell, obs[s2], obs[s3]);
2067 const auto* buf = results[0];
2068 if (buf ==
nullptr)
continue;
2073 auto lower_bound = {bf1_first, bf2_first, bf3_first};
2074 auto upper_bound = {bf1_first + n1, bf2_first + n2, bf3_first + n3};
2075 auto view = btas::make_view(
2076 Zxy.range().slice(lower_bound, upper_bound), Zxy.storage());
2077 std::copy(buf, buf + n123, view.begin());
2089 double ints_time = 0;
2090 for (
const auto& timer : timers) ints_time += timer.read(0);
2091 std::cout <<
"time for Zxy integrals = " << ints_time
2092 <<
" (total from all threads)" << std::endl;
2093 double copy_time = 0;
2094 for (
const auto& timer : timers) copy_time += timer.read(1);
2095 std::cout <<
"time for copying into BTAS = " << copy_time
2096 <<
" (total from all threads)" << std::endl;
2097 std::cout <<
"wall time for Zxy integrals + copy = " << wall_timer.
read(0)
2102 Matrix V = compute_2body_2index_ints(dfbs);
2103 Eigen::LLT<Matrix> V_LLt(V);
2104 Matrix I = Matrix::Identity(ndf, ndf);
2105 auto L = V_LLt.matrixL();
2107 Matrix Linv = L.solve(I).transpose();
2116 Tensor2d K{ndf, ndf};
2117 std::copy(Linv.data(), Linv.data() + ndf * ndf, K.begin());
2119 xyK = Tensor3d{n, n, ndf};
2120 btas::contract(1.0, Zxy, {1, 2, 3}, K, {1, 4}, 0.0, xyK, {2, 3, 4});
2121 Zxy = Tensor3d{0, 0, 0};
2124 std::cout <<
"time for integrals metric tform = " << timers[0].read(2)
2131 const auto nocc = Cocc.cols();
2132 Tensor2d Co{n, nocc};
2133 std::copy(Cocc.data(), Cocc.data() + n * nocc, Co.begin());
2134 Tensor3d xiK{n, nocc, ndf};
2135 btas::contract(1.0, xyK, {1, 2, 3}, Co, {2, 4}, 0.0, xiK, {1, 4, 3});
2138 btas::contract(1.0, xiK, {1, 2, 3}, xiK, {4, 2, 3}, 0.0, G, {1, 4});
2141 std::cout <<
"time for exchange = " << timers[0].read(3) << std::endl;
2147 btas::contract(1.0, xiK, {1, 2, 3}, Co, {1, 2}, 0.0, Jtmp, {3});
2148 xiK = Tensor3d{0, 0, 0};
2149 btas::contract(2.0, xyK, {1, 2, 3}, Jtmp, {3}, -1.0, G, {1, 2});
2152 std::cout <<
"time for coulomb = " << timers[0].read(4) << std::endl;
2155 Matrix result(n, n);
2156 std::copy(G.cbegin(), G.cend(), result.data());
2162void api_basic_compile_test(
const BasisSet& obs) {
2164 Engine onebody_engine(
2170 auto shell2bf = obs.shell2bf();
2171 const auto& results = onebody_engine.results();
2172 for (
auto s1 = 0; s1 != obs.size(); ++s1) {
2173 for (
auto s2 = 0; s2 != obs.size(); ++s2) {
2174 std::cout <<
"compute shell set {" << s1 <<
"," << s2 <<
"} ... ";
2175 onebody_engine.compute(obs[s1], obs[s2]);
2176 const auto* ints_shellset = results[0];
2177 std::cout <<
"done" << std::endl;
2179 auto bf1 = shell2bf[s1];
2180 auto n1 = obs[s1].size();
2181 auto bf2 = shell2bf[s2];
2182 auto n2 = obs[s2].size();
2186 for (
auto f1 = 0; f1 != n1; ++f1)
2187 for (
auto f2 = 0; f2 != n2; ++f2)
2188 std::cout <<
" " << bf1 + f1 <<
" " << bf2 + f2 <<
" "
2189 << ints_shellset[f1 * n2 + f2] << std::endl;
2193 using libint2::Operator;
2195 std::vector<std::pair<double, double>> cgtg_params{
2196 {0.1, 0.2}, {0.3, 0.4}, {0.5, 0.6}};
2198 auto K = compute_schwarz_ints<Operator::cgtg>(obs, obs,
false, cgtg_params);
2199 std::cout <<
"cGTG Schwarz ints\n" << K << std::endl;
2202 auto K = compute_schwarz_ints<Operator::cgtg_x_coulomb>(obs, obs,
false,
2204 std::cout <<
"cGTG/r12 Schwarz ints\n" << K << std::endl;
2208 compute_schwarz_ints<Operator::delcgtg2>(obs, obs,
false, cgtg_params);
2209 std::cout <<
"||Del.cGTG||^2 Schwarz ints\n" << K << std::endl;
2213 Engine eri4_engine(Operator::coulomb, obs.max_nprim(), obs.max_l());
2214 Engine eri2_engine = eri4_engine;
2215 eri2_engine.set(BraKet::xs_xs);
2216 auto shell2bf = obs.shell2bf();
2217 const auto& results4 = eri4_engine.results();
2218 const auto& results2 = eri2_engine.results();
2219 for (
auto s1 = 0; s1 != obs.size(); ++s1) {
2220 for (
auto s2 = 0; s2 != obs.size(); ++s2) {
2222 eri2_engine.compute(obs[s1], obs[s2]);
2224 auto bf1 = shell2bf[s1];
2225 auto n1 = obs[s1].size();
2226 auto bf2 = shell2bf[s2];
2227 auto n2 = obs[s2].size();
2229 const auto* buf4 = results4[0];
2230 const auto* buf2 = results2[0];
2234 for (
auto f1 = 0, f12 = 0; f1 != n1; ++f1)
2235 for (
auto f2 = 0; f2 != n2; ++f2, ++f12)
2236 assert(std::abs(buf4[f12] - buf2[f12]) < 1e-12 &&
2237 "2-center ints test failed");
2242 Engine eri4_engine(Operator::coulomb, obs.max_nprim(), obs.max_l());
2243 Engine eri3_engine = eri4_engine;
2244 eri3_engine.set(BraKet::xs_xx);
2245 auto shell2bf = obs.shell2bf();
2246 const auto& results4 = eri4_engine.results();
2247 const auto& results3 = eri3_engine.results();
2248 for (
auto s1 = 0; s1 != obs.size(); ++s1) {
2249 for (
auto s2 = 0; s2 != obs.size(); ++s2) {
2250 for (
auto s3 = 0; s3 != obs.size(); ++s3) {
2251 eri4_engine.compute(obs[s1],
Shell::unit(), obs[s2], obs[s3]);
2252 eri3_engine.compute(obs[s1], obs[s2], obs[s3]);
2254 auto bf1 = shell2bf[s1];
2255 auto n1 = obs[s1].size();
2256 auto bf2 = shell2bf[s2];
2259 auto bf3 = shell2bf[s3];
2260 auto n3 = obs[s3].size();
2262 const auto* buf4 = results4[0];
2263 const auto* buf3 = results3[0];
2267 for (
auto f1 = 0, f123 = 0; f1 != n1; ++f1)
2268 for (
auto f2 = 0; f2 != n2; ++f2)
2269 for (
auto f3 = 0; f3 != n3; ++f3, ++f123)
2270 assert(std::abs(buf4[f123] - buf3[f123]) < 1e-12 &&
2271 "3-center ints test failed");
2277#if LIBINT2_DERIV_ERI_ORDER
2279 Engine eri4_engine(Operator::coulomb, obs.max_nprim(), obs.max_l(), 1);
2280 Engine eri2_engine = eri4_engine;
2281 eri2_engine.set(BraKet::xs_xs);
2282 auto shell2bf = obs.shell2bf();
2283 const auto& results4 = eri4_engine.results();
2284 const auto& results2 = eri2_engine.results();
2285 for (
auto s1 = 0; s1 != obs.size(); ++s1) {
2286 for (
auto s2 = 0; s2 != obs.size(); ++s2) {
2288 eri2_engine.compute(obs[s1], obs[s2]);
2290 auto bf1 = shell2bf[s1];
2291 auto n1 = obs[s1].size();
2292 auto bf2 = shell2bf[s2];
2293 auto n2 = obs[s2].size();
2296 for (
auto d = 0; d != 6; ++d) {
2297 const auto* buf4 = results4[d < 3 ? d : d + 3];
2298 const auto* buf2 = results2[d];
2302 for (
auto f1 = 0, f12 = 0; f1 != n1; ++f1)
2303 for (
auto f2 = 0; f2 != n2; ++f2, ++f12)
2304 assert(std::abs(buf4[f12] - buf2[f12]) < 1e-12 &&
2305 "deriv 2-center ints test failed");
2312#if LIBINT2_DERIV_ERI_ORDER > 1
2314 Engine eri4_engine(Operator::coulomb, obs.max_nprim(), obs.max_l(), 2);
2315 Engine eri2_engine = eri4_engine;
2316 eri2_engine.set(BraKet::xs_xs);
2317 auto shell2bf = obs.shell2bf();
2318 const auto& results4 = eri4_engine.results();
2319 const auto& results2 = eri2_engine.results();
2320 for (
auto s1 = 0; s1 != obs.size(); ++s1) {
2321 for (
auto s2 = 0; s2 != obs.size(); ++s2) {
2323 eri2_engine.compute(obs[s1], obs[s2]);
2325 auto bf1 = shell2bf[s1];
2326 auto n1 = obs[s1].size();
2327 auto bf2 = shell2bf[s2];
2328 auto n2 = obs[s2].size();
2331 for (
auto d1 = 0, d12 = 0; d1 != 6; ++d1) {
2332 const auto dd1 = d1 < 3 ? d1 : d1 + 3;
2333 for (
auto d2 = d1; d2 != 6; ++d2, ++d12) {
2334 const auto dd2 = d2 < 3 ? d2 : d2 + 3;
2335 const auto dd12 = dd1 * (24 - dd1 - 1) / 2 + dd2;
2336 const auto* buf4 = results4[dd12];
2337 const auto* buf2 = results2[d12];
2342 for (
auto f1 = 0, f12 = 0; f1 != n1; ++f1)
2343 for (
auto f2 = 0; f2 != n2; ++f2, ++f12)
2344 assert(std::abs(buf4[f12] - buf2[f12]) < 1e-12 &&
2345 "2nd deriv 2-center ints test failed");
DIIS (`‘direct inversion of iterative subspace’') extrapolation.
Definition diis.h:140
Timers aggregates N C++11 "timers"; used to high-resolution profile stages of integral computation.
Definition timer.h:39
double read(size_t t) const
reads value (in seconds) of timer t , converted to double
Definition timer.h:72
dur_t stop(size_t t)
stops timer t
Definition timer.h:65
void clear()
resets timers to zero
Definition timer.h:74
void set_now_overhead(size_t ns)
use this to report the overhead of now() call; if set, the reported timings will be adjusted for this...
Definition timer.h:58
void start(size_t t)
starts timer t
Definition timer.h:61
Defaults definitions for various parameters assumed by Libint.
Definition algebra.cc:24
void parallel_do(Lambda &lambda)
fires off nthreads instances of lambda in parallel
Definition 1body.h:170
@ Schwarz
Schwarz screening method using Frobenius norm:
void initialize(bool verbose=false)
initializes the libint library if not currently initialized, no-op otherwise
Definition initialize.h:76
std::vector< std::pair< double, std::array< double, 3 > > > make_point_charges(const std::vector< libint2::Atom > &atoms)
converts a vector of Atoms to a vector of point charges
Definition atom.h:241
BraKet
types of shell sets supported by Engine, in chemist notation (i.e.
Definition include/libint2/braket.h:36
std::vector< Atom > read_dotxyz(std::istream &is, const double bohr_to_angstrom=constants::codata_2018::bohr_to_angstrom)
reads the list of atoms from a file in the standard XYZ format supported by most chemistry software (...
Definition atom.h:209
void finalize()
finalizes the libint library.
Definition initialize.h:105
const std::vector< Real > & sto3g_ao_occupation_vector(size_t Z)
computes average orbital occupancies in the ground state of a neutral atoms
Definition sto3g_atomic_density.h:86
size_t sto3g_num_ao(size_t Z)
Definition sto3g_atomic_density.h:57
Iterates over all partitions of a non-negative integer into nonnegative integers in reverse lexicog...
Definition intpart_iter.h:74
generally-contracted Solid-Harmonic/Cartesion Gaussian Shell
Definition shell.h:720
static bool do_enforce_unit_normalization(defaultable_boolean flag=defaultable_boolean())
sets and/or reports whether the auto-renormalization to unity is set if called without arguments,...
Definition shell.h:895
static const Shell & unit()
Definition shell.h:906