MADNESS 0.10.1
SCF.h
Go to the documentation of this file.
1/*
2 This file is part of MADNESS.
3
4 Copyright (C) 2007,2010 Oak Ridge National Laboratory
5
6 This program is free software; you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation; either version 2 of the License, or
9 (at your option) any later version.
10
11 This program is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with this program; if not, write to the Free Software
18 Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19
20 For more information please contact:
21
22 Robert J. Harrison
23 Oak Ridge National Laboratory
24 One Bethel Valley Road
25 P.O. Box 2008, MS-6367
26
27 email: harrisonrj@ornl.gov
28 tel: 865-241-3937
29 fax: 865-572-0680
30
31
32 $Id$
33 */
34
35/// \file moldft.cc
36/// \brief Molecular HF and DFT code
37/// \defgroup moldft The molecular density functional and Hartree-Fock code
38
39
40#ifndef MADNESS_CHEM_SCF_H__INCLUDED
41#define MADNESS_CHEM_SCF_H__INCLUDED
42
43#include <memory>
44
46#include <madness/mra/mra.h>
47
59#include<madness/chem/pcm.h>
61
62#include <madness/tensor/tensor_json.hpp>
63#include <memory>
64
65namespace madness {
66
67typedef std::shared_ptr<WorldDCPmapInterface<Key<3> > > pmapT;
69typedef std::shared_ptr<FunctionFunctorInterface<double, 3> > functorT;
71typedef std::vector<functionT> vecfuncT;
72typedef std::pair<vecfuncT, vecfuncT> pairvecfuncT;
73typedef std::vector<pairvecfuncT> subspaceT;
78typedef std::shared_ptr<operatorT> poperatorT;
80typedef std::vector<complex_functionT> cvecfuncT;
81
82
83template<typename T, int NDIM>
84struct lbcost {
85 double leaf_value;
87
89
90 double operator()(const Key<NDIM>& key, const FunctionNode<T, NDIM>& node) const {
91 if (key.level() < 1) {
92 return 100.0 * (leaf_value + parent_value);
93 } else if (node.is_leaf()) {
94 return leaf_value;
95 } else {
96 return parent_value;
97 }
98 }
99};
100
101
102inline double mask1(double x) {
103 /* Iterated first beta function to switch smoothly
104 from 0->1 in [0,1]. n iterations produce 2*n-1
105 zero derivatives at the end points. Order of polyn
106 is 3^n.
107
108 Currently use one iteration so that first deriv.
109 is zero at interior boundary and is exactly representable
110 by low order multiwavelet without refinement */
111
112 x = (x * x * (3. - 2. * x));
113 return x;
114}
115
116static double mask3(const coordT& ruser) {
117 coordT rsim;
119 double x = rsim[0], y = rsim[1], z = rsim[2];
120 double lo = 0.0625, hi = 1.0 - lo, result = 1.0;
121 double rlo = 1.0 / lo;
122
123 if (x < lo)
124 result *= mask1(x * rlo);
125 else if (x > hi)
126 result *= mask1((1.0 - x) * rlo);
127 if (y < lo)
128 result *= mask1(y * rlo);
129 else if (y > hi)
130 result *= mask1((1.0 - y) * rlo);
131 if (z < lo)
132 result *= mask1(z * rlo);
133 else if (z > hi)
134 result *= mask1((1.0 - z) * rlo);
135
136 return result;
137}
138
139/// A MADNESS functor to compute either x, y, or z
140class DipoleFunctor : public FunctionFunctorInterface<double, 3> {
141private:
142 const int axis;
143public:
145
146 double operator()(const coordT& x) const {
147 return x[axis];
148 }
149};
150
151
152/// A MADNESS functor to compute the cartesian moment x^i * y^j * z^k (i, j, k integer and >= 0)
153class MomentFunctor : public FunctionFunctorInterface<double, 3> {
154private:
155 const int i, j, k;
156public:
157 MomentFunctor(int i, int j, int k) : i(i), j(j), k(k) {}
158
159 MomentFunctor(const std::vector<int>& x) : i(x[0]), j(x[1]), k(x[2]) {}
160
161 double operator()(const coordT& r) const {
162 double xi = 1.0, yj = 1.0, zk = 1.0;
163 for (int p = 0; p < i; ++p) xi *= r[0];
164 for (int p = 0; p < j; ++p) yj *= r[1];
165 for (int p = 0; p < k; ++p) zk *= r[2];
166 return xi * yj * zk;
167 }
168};
169
170 class scf_data {
171
172 std::map<std::string, std::vector<double>> e_data;
175 int iter;
176 public:
177
178 scf_data();
179
180 void to_json(json &j) const;
181
182 void print_data();
183
184 void add_data(std::map<std::string, double> values);
185
186 void add_gradient(const Tensor<double> &grad);
187 };
188
189
190class SCF {
191public:
192 std::shared_ptr<PotentialManager> potentialmanager;
193 std::shared_ptr<GTHPseudopotential<double> > gthpseudopotential;
200
202
203 /// alpha and beta molecular orbitals
205
206 /// sets of orbitals grouped by their orbital energies (for localization?)
207 /// only orbitals within the same set will be mixed to localize
208 std::vector<int> aset, bset;
209
210 /// MRA projection of the minimal basis set
212
213 std::vector<int> at_to_bf, at_nbf;
214
215 /// occupation numbers for alpha and beta orbitals
217
218 /// orbital energies for alpha and beta orbitals
221 std::vector<std::shared_ptr<real_derivative_3d> > gradop;
222 double vtol;
224 double converged_for_thresh=1.e10; ///< mos are converged for this threshold
225 double converged_for_dconv=1.e10; ///< mos are converged for this density
226 double converged_for_tconv=1.e10; ///< derivatives of mos are converged for this threshold
227
228 /// forwarding constructor
229 SCF(World& world, const commandlineparser& parser)
230 : SCF(world, CalculationParameters(world, parser), Molecule(world, parser)) {}
231
232 /// collective constructor for SCF uses contents of file \c filename and broadcasts to all nodes
233 SCF(World& world, const CalculationParameters& param, const Molecule& molecule);
234
235 void copy_data(World& world, const SCF& other);
236
237 static void help() {
238 print_header2("help page for MOLDFT ");
239 print("The moldft code computes Hartree-Fock and DFT energies and gradients, It is the fastest code in MADNESS");
240 print("and considered the reference implementation. No nuclear correlation factor can be used");
241 print("SCF orbitals are the basis for post-SCF calculations like");
242 print("excitation energies (cis), correlation energies (cc2), local potentials (oep), etc\n\n");
243 print("You can print all available calculation parameters by running\n");
244 print("moldft --print_parameters\n");
245 print("You can perform a simple calculation by running\n");
246 print("moldft --geometry=h2o.xyz\n");
247 print("provided you have an xyz file in your directory.");
248
249 }
250
251 static void print_parameters() {
253 print("default parameters for the moldft program are");
254 param.print("dft", "end");
255 print("\n\nthe molecular geometry must be specified in a separate block:");
257 }
258
259 void set_print_timings(const bool value);
260
261 template<std::size_t NDIM>
262 void set_protocol(World& world, double thresh) {
263 int k;
264 // Allow for imprecise conversion of threshold
265 if (thresh >= 0.9e-2)
266 k = 4;
267 else if (thresh >= 0.9e-4)
268 k = 6;
269 else if (thresh >= 0.9e-6)
270 k = 8;
271 else if (thresh >= 0.9e-8)
272 k = 10;
273 else
274 k = 12;
275
276 // k defaults to make sense with thresh, override by providing k in input file
277 if (param.k() == -1) {
279 // param.k=k;
280 } else {
282 }
283 // don't forget to adapt the molecular smoothing parameter!! NO ... it is independent
284 // molecule.set_eprec(std::min(thresh,molecule.get_eprec()));
288 // FunctionDefaults<NDIM>::set_truncate_mode(1);
294 double safety = 0.1;
298
299 // Update coefficients if using a different derivative
300 if (param.deriv() == "bspline") {
301 for (int i = 0; i < 3; ++i) (*gradop[i]).set_bspline1();
302 } else if (param.deriv() == "ble") {
303 for (int i = 0; i < 3; ++i) (*gradop[i]).set_ble1();
304 }
305
306 mask = functionT(factoryT(world).f(mask3).initial_level(4).norefine());
307 if (world.rank() == 0 and param.print_level() > 1) {
308 print("\nSolving NDIM=", NDIM, " with thresh", thresh, " k",
309 FunctionDefaults<NDIM>::get_k(), " conv", std::max(thresh, param.dconv()), "\n");
310 }
311 }
312
313 /// getter for the molecular orbitals, alpha spin
314 const vecfuncT& get_amo() const { return amo; }
315
316 /// getter for the molecular orbitals, beta spin
317 const vecfuncT& get_bmo() const { return bmo; }
318
319 /// getter for the occupation numbers, alpha spin
320 const tensorT& get_aocc() const { return aocc; }
321
322 /// getter for the occupation numbers, alpha spin
323 const tensorT& get_bocc() const { return bocc; }
324
325 bool is_spin_restricted() const { return param.get<bool>("spin_restricted"); }
326
327 void save_mos(World& world);
328
329 void load_mos(World& world);
330
331 bool restart_aos(World& world);
332
333 void do_plots(World& world);
334
335 void project(World& world);
336
337 void make_nuclear_potential(World& world);
338
340
342 const Molecule& molecule);
343
344 void reset_aobasis(const std::string& aobasisname) {
345 aobasis = AtomicBasisSet(); // reset
347 }
348
349 /// group orbitals into sets of similar orbital energies for localization
350
351 /// @param[in] eps orbital energies
352 /// @param[in] occ occupation numbers
353 /// @param[in] nmo number of MOs for the given spin
354 /// @return vector of length nmo with the set index for each MO
355 std::vector<int> group_orbital_sets(World& world, const tensorT& eps,
356 const tensorT& occ, const int nmo) const;
357
358 static void analyze_vectors(World& world, const vecfuncT& mo,
359 const vecfuncT& ao, double vtol,
360 const Molecule& molecule, const int print_level,
361 const AtomicBasisSet& aobasis, const tensorT& occ = tensorT(),
362 const tensorT& energy = tensorT(), const std::vector<int>& set = std::vector<int>());
363
364 distmatT kinetic_energy_matrix(World& world, const vecfuncT& v) const;
365
366 void get_initial_orbitals(World& world);
367
368 void initial_guess(World& world);
369
370 void initial_guess_from_nwchem(World& world);
371
372 void initial_load_bal(World& world);
373
374 functionT make_density(World& world, const tensorT& occ, const vecfuncT& v) const;
375
376 functionT make_density(World& world, const tensorT& occ, const cvecfuncT& v);
377
378 static std::vector<poperatorT> make_bsh_operators(World& world, const tensorT& evals,
380
381 // Used only for initial guess that is always spin-restricted LDA
382 static functionT make_lda_potential(World& world, const functionT& arho);
383
384
385 // functionT make_dft_potential(World & world, const vecfuncT& vf, int ispin, int what)
386 // {
387 // return multiop_values<double, xc_potential, 3>(xc_potential(xc, ispin, what), vf);
388 // }
389
390 double make_dft_energy(World& world, const vecfuncT& vf, int ispin) {
392 return vlda.trace();
393 }
394
395 vecfuncT apply_potential(World& world, const tensorT& occ,
396 const vecfuncT& amo,
397 const functionT& vlocal, double& exc, double& enl, int ispin);
398
399 tensorT derivatives(World& world, const functionT& rho) const;
400
401 /// compute the total dipole moment of the molecule
402
403 /// @param[in] rho the total (alpha + beta) density
404 /// @return the x,y,z components of the el. + nucl. dipole moment
405 tensorT dipole(World& world, const functionT& rho) const;
406
407 void vector_stats(const std::vector<double>& v, double& rms,
408 double& maxabsval) const;
409
410 vecfuncT compute_residual(World& world, tensorT& occ, tensorT& fock,
411 const vecfuncT& psi, vecfuncT& Vpsi, double& err);
412
414 const vecfuncT& Vpsi, const tensorT& occ,
415 double& ekinetic) const;
416
417 /// make the Coulomb potential given the total density
419 return apply(*coulop, rho);
420 }
421
422 /// Compute the two-electron integrals over the provided set of orbitals
423
424 /// Returned is a *replicated* tensor of \f$(ij|kl)\f$ with \f$i>=j\f$
425 /// and \f$k>=l\f$. The symmetry \f$(ij|kl)=(kl|ij)\f$ is enforced.
426 Tensor<double> twoint(World& world, const vecfuncT& psi) const;
427
428 /// compute the unitary transformation that diagonalizes the fock matrix
429
430 /// @param[in] world the world
431 /// @param[in] overlap the overlap matrix of the orbitals
432 /// @param[in,out] fock the fock matrix; diagonal upon exit
433 /// @param[out] evals the orbital energies
434 /// @param[in] occ the occupation numbers
435 /// @param[in] thresh_degenerate threshold for orbitals being degenerate
436 /// @return the unitary matrix U: U^T F U = evals
438 tensorT& fock, tensorT& evals, const tensorT& occ,
439 const double thresh_degenerate) const;
440
441
442 /// diagonalize the fock matrix, taking care of degenerate states
443
444 /// Vpsi is passed in to make sure orbitals and Vpsi are in phase
445 /// @param[in] world the world
446 /// @param[in,out] fock the fock matrix (diagonal upon exit)
447 /// @param[in,out] psi the orbitals
448 /// @param[in,out] Vpsi the orbital times the potential
449 /// @param[out] evals the orbital energies
450 /// @param[in] occ occupation numbers
451 /// @param[in] thresh threshold for rotation and truncation
452 /// @return the unitary matrix U: U^T F U = evals
455 const tensorT& occ, const double thresh) const;
456
457
460
461
462 void rotate_subspace(World& world, const tensorT& U, subspaceT& subspace,
463 int lo, int nfunc, double trantol) const;
464
465 void rotate_subspace(World& world, const distmatT& U, subspaceT& subspace,
466 int lo, int nfunc, double trantol) const;
467
468 void update_subspace(World& world,
472 double& bsh_residual, double& update_residual);
473
474 /// perform step restriction following the KAIN solver
475
476 /// undo the rotation from the KAIN solver if the rotation exceeds the
477 /// maxrotn parameter
478 /// @param[in] world the world
479 /// @param[in] mo vector of orbitals from previous iteration
480 /// @param[in,out] mo_new vector of orbitals from the KAIN solver
481 /// @param[in] spin "alpha" or "beta" for user information
482 /// @return max residual
483 double do_step_restriction(World& world, const vecfuncT& mo,
484 vecfuncT& mo_new, std::string spin) const;
485
486 /// orthonormalize the vectors
487
488 /// @param[in] world the world
489 /// @param[in,out] amo_new the vectors to be orthonormalized
490 void orthonormalize(World& world, vecfuncT& amo_new) const;
491
492 void orthonormalize(World& world, vecfuncT& amo_new, int nocc) const;
493
494 // For given protocol, solve the DFT/HF/response equations
495 void solve(World& world);
496
497 void output_calc_info_schema() const;
498
499 void output_scf_info_schema(const std::map<std::string, double> &vals,
500 const tensorT &dipole_T) const;
501
502};
503
504// Computes molecular energy as a function of the geometry
505// This is cludgy ... need better factorization of functionality
506// between calculation, main program and this ... or just merge it all.
510 mutable double coords_sum; // sum of square of coords at last solved geometry
511
512public:
515
516 std::string name() const { return "Molecularenerg"; }
517
518 bool selftest() { return true; }
519
520 bool provides_gradient() const { return true; }
521
522 double value(const Tensor<double>& x) {
523 double xsq = x.sumsq();
524 if (xsq == coords_sum) {
525 return calc.current_energy;
526 }
528 coords_sum = xsq;
529
530 // read converged wave function from disk if there is one
531 if (calc.param.no_compute()) {
535 return calc.current_energy;
536 }
537
538 // initialize the PCM solver for this geometry
539 if (calc.param.pcm_data() != "none") {
541 }
542
544
545 // AOs are needed for final analysis, and for localization
546 calc.reset_aobasis("sto-3g");
547 calc.ao.clear(); world.gop.fence();
549
550 // The below is missing convergence test logic, etc.
551
552 // Make the nuclear potential, initial orbitals, etc.
553 for (unsigned int proto = 0; proto < calc.param.protocol().size(); proto++) {
554
555 int nvalpha = calc.param.nmo_alpha() - calc.param.nalpha();
556 int nvbeta = calc.param.nmo_beta() - calc.param.nbeta();
558
559 //repeat with gradually decreasing nvirt, only for first protocol
560 if (proto == 0 && nvalpha > 0) {
561 nvalpha_start = nvalpha * calc.param.nv_factor();
562 } else {
563 nvalpha_start = nvalpha;
564 }
565
567
568 for (int nv = nvalpha_start; nv >= nvalpha; nv -= nvalpha) {
569
570 if (nv > 0 && world.rank() == 0) std::cout << "Running with " << nv << " virtual states" << std::endl;
571
573 // check whether this is sensible for spin restricted case
575 if (nvbeta == nvalpha) {
577 } else {
578 calc.param.set_user_defined_value("nmo_beta", calc.param.nbeta() + nv + nvbeta - nvalpha);
579 }
580 }
581
584
585 if (nv != nv_old) {
586 calc.amo.resize(calc.param.nmo_alpha());
587 calc.bmo.resize(calc.param.nmo_beta());
588
590 for (int i = 0; i < calc.param.nalpha(); ++i)
591 calc.aocc[i] = 1.0;
592
594 for (int i = 0; i < calc.param.nbeta(); ++i)
595 calc.bocc[i] = 1.0;
596
597 // might need to resize aset, bset, but for the moment this doesn't seem to be necessary
598
599 }
600
601 // project orbitals into higher k
602 if (proto > 0) calc.project(world);
603
604 // If the basis for the inital guess was not sto-3g
605 // switch to sto-3g since this is needed for analysis
606 // of the MOs and orbital localization
607 // Only do this if not starting from NWChem.
608 // analysis will be done on NWChem orbitals.
609
610 if (calc.param.aobasis() != "sto-3g") { // was also && calc.param.nwfile() == "none"
611 calc.reset_aobasis("sto-3g");
612 }
613 calc.ao.clear(); world.gop.fence();
616
617 if (calc.param.save())
619
620 nv_old = nv;
621 // exit loop over decreasing nvirt if nvirt=0
622 if (nv == 0) break;
623
624 }
625
626 }
627 return calc.current_energy;
628 }
629
631 value(x); // Ensures DFT equations are solved at this geometry
632
634 functionT brho = rho;
637 rho.gaxpy(1.0, brho, 1.0);
638
639 return calc.derivatives(world, rho);
640 }
641
642
644 value(molecule.get_all_coords().flat()); // Ensures DFT equations are solved at this geometry
645
647 functionT brho = rho;
650 rho.gaxpy(1.0, brho, 1.0);
651
654 }
655
657 nlohmann::json j = {};
661
663
664 nlohmann::json calc_precision={ };
666 calc_precision["dconv"]=calc.param.dconv();
667 calc_precision["econv"]=calc.param.econv();
670
671 auto mol_json=this->calc.molecule.to_json();
672
673 int_vals.push_back({"calcinfo_nmo", param.nmo_alpha() + param.nmo_beta()});
674 int_vals.push_back({"calcinfo_nalpha", param.nalpha()});
675 int_vals.push_back({"calcinfo_nbeta", param.nbeta()});
676 int_vals.push_back({"calcinfo_natom", calc.molecule.natom()});
677
678
679 to_json(j, int_vals);
680 double_vals.push_back({"return_energy", value(calc.molecule.get_all_coords().flat())});
682 double_tensor_vals.push_back({"scf_eigenvalues_a", calc.aeps});
683 if (param.nbeta() != 0 && !param.spin_restricted()) {
684 double_tensor_vals.push_back({"scf_eigenvalues_b", calc.beps});
685 }
686
688 param.to_json(j);
690
691 j["precision"]=calc_precision;
692 j["molecule"]=mol_json;
693
694 output_schema(param.prefix()+".calc_info", j);
695 }
696
697
698};
699}
700
701#endif /* SCF_H_ */
702
Operators for the molecular HF and DFT code.
A MADNESS functor to compute either x, y, or z.
Definition preal.cc:115
Contracted Gaussian basis.
Definition madness/chem/molecularbasis.h:469
void read_file(std::string filename)
read the atomic basis set from file
Definition molecularbasis.cc:119
DipoleFunctor(int axis)
Definition SCF.h:144
double operator()(const coordT &x) const
Definition SCF.h:146
const int axis
Definition solver.h:167
Manages data associated with a row/column/block distributed array.
Definition distributed_matrix.h:388
FunctionDefaults holds default paramaters as static class members.
Definition funcdefaults.h:100
static int get_k()
Returns the default wavelet order.
Definition funcdefaults.h:164
static void set_apply_randomize(bool value)
Sets the random load balancing for integral operators flag.
Definition funcdefaults.h:295
static void set_thresh(double value)
Sets the default threshold.
Definition funcdefaults.h:184
static void set_k(int value)
Sets the default wavelet order.
Definition funcdefaults.h:171
static const double & get_thresh()
Returns the default threshold.
Definition funcdefaults.h:177
static void set_autorefine(bool value)
Sets the default adaptive autorefinement flag.
Definition funcdefaults.h:261
static void set_project_randomize(bool value)
Sets the random load balancing for projection flag.
Definition funcdefaults.h:306
static void set_initial_level(int value)
Sets the default initial projection level.
Definition funcdefaults.h:201
static void set_cubic_cell(double lo, double hi)
Sets the user cell to be cubic with each dimension having range [lo,hi].
Definition funcdefaults.h:363
static void set_refine(bool value)
Sets the default adaptive refinement flag.
Definition funcdefaults.h:249
FunctionFactory implements the named-parameter idiom for Function.
Definition function_factory.h:86
Abstract base class interface required for functors used as input to Functions.
Definition function_interface.h:68
FunctionNode holds the coefficients, etc., at each node of the 2^NDIM-tree.
Definition funcimpl.h:127
A multiresolution adaptive numerical function.
Definition mra.h:139
Function< T, NDIM > & gaxpy(const T &alpha, const Function< Q, NDIM > &other, const R &beta, bool fence=true)
Inplace, general bi-linear operation in wavelet basis. No communication except for optional fence.
Definition mra.h:1055
Key is the index for a node of the 2^NDIM-tree.
Definition key.h:69
Level level() const
Definition key.h:168
Definition SCF.h:507
bool selftest()
Definition SCF.h:518
bool provides_gradient() const
Override this to return true if the derivative is implemented.
Definition SCF.h:520
double value(const Tensor< double > &x)
Should return the value of the objective function.
Definition SCF.h:522
World & world
Definition SCF.h:508
SCF & calc
Definition SCF.h:509
double coords_sum
Definition SCF.h:510
void output_calc_info_schema()
Definition SCF.h:656
madness::Tensor< double > gradient(const Tensor< double > &x)
Should return the derivative of the function.
Definition SCF.h:630
void energy_and_gradient(const Molecule &molecule, double &energy, Tensor< double > &gradient)
Definition SCF.h:643
std::string name() const
Definition SCF.h:516
MolecularEnergy(World &world, SCF &calc)
Definition SCF.h:513
Definition molecule.h:129
void set_all_coords(const madness::Tensor< double > &newcoords)
Definition molecule.cc:424
madness::Tensor< double > get_all_coords() const
Definition molecule.cc:402
size_t natom() const
Definition molecule.h:415
static void print_parameters()
Definition molecule.cc:115
json to_json() const
Definition molecule.cc:462
GeometryParameters parameters
Definition molecule.h:286
A MADNESS functor to compute the cartesian moment x^i * y^j * z^k (i, j, k integer and >= 0)
Definition SCF.h:153
const int j
Definition SCF.h:155
double operator()(const coordT &r) const
Definition SCF.h:161
MomentFunctor(int i, int j, int k)
Definition SCF.h:157
const int i
Definition SCF.h:155
MomentFunctor(const std::vector< int > &x)
Definition SCF.h:159
const int k
Definition SCF.h:155
interface class to the PCMSolver library
Definition pcm.h:52
void set_user_defined_value(const std::string &key, const T &value)
Definition QCCalculationParametersBase.h:522
T get(const std::string key) const
Definition QCCalculationParametersBase.h:305
void print(const std::string header="", const std::string footer="") const
print all parameters
Definition QCCalculationParametersBase.cc:22
class implementing properties of QC models
Definition QCPropertyInterface.h:11
Definition SCF.h:190
void copy_data(World &world, const SCF &other)
Definition SCF.cc:267
void do_plots(World &world)
Definition SCF.cc:503
tensorT derivatives(World &world, const functionT &rho) const
Definition SCF.cc:1432
void save_mos(World &world)
Definition SCF.cc:281
void output_scf_info_schema(const std::map< std::string, double > &vals, const tensorT &dipole_T) const
Definition SCF.cc:143
static vecfuncT project_ao_basis_only(World &world, const AtomicBasisSet &aobasis, const Molecule &molecule)
Definition SCF.cc:597
const tensorT & get_aocc() const
getter for the occupation numbers, alpha spin
Definition SCF.h:320
std::shared_ptr< GTHPseudopotential< double > > gthpseudopotential
Definition SCF.h:193
void make_nuclear_potential(World &world)
Definition SCF.cc:572
void get_initial_orbitals(World &world)
get the initial orbitals for a calculation
Definition SCF.cc:453
std::vector< int > aset
Definition SCF.h:208
static void print_parameters()
Definition SCF.h:251
AtomicBasisSet aobasis
Definition SCF.h:198
SCF(World &world, const commandlineparser &parser)
forwarding constructor
Definition SCF.h:229
vecfuncT ao
MRA projection of the minimal basis set.
Definition SCF.h:211
scf_data e_data
Definition SCF.h:201
void initial_guess(World &world)
Definition SCF.cc:985
vecfuncT apply_potential(World &world, const tensorT &occ, const vecfuncT &amo, const functionT &vlocal, double &exc, double &enl, int ispin)
Definition SCF.cc:1342
vecfuncT amo
alpha and beta molecular orbitals
Definition SCF.h:204
vecfuncT compute_residual(World &world, tensorT &occ, tensorT &fock, const vecfuncT &psi, vecfuncT &Vpsi, double &err)
Definition SCF.cc:1543
const tensorT & get_bocc() const
getter for the occupation numbers, alpha spin
Definition SCF.h:323
std::vector< int > at_to_bf
Definition SCF.h:213
std::vector< int > bset
Definition SCF.h:208
distmatT kinetic_energy_matrix(World &world, const vecfuncT &v) const
Definition SCF.cc:669
tensorT bocc
Definition SCF.h:216
double vtol
Definition SCF.h:222
tensorT dipole(World &world, const functionT &rho) const
compute the total dipole moment of the molecule
Definition SCF.cc:1505
PCM pcm
Definition SCF.h:197
bool is_spin_restricted() const
Definition SCF.h:325
poperatorT coulop
Definition SCF.h:220
double make_dft_energy(World &world, const vecfuncT &vf, int ispin)
Definition SCF.h:390
void output_calc_info_schema() const
Definition SCF.cc:158
static void help()
Definition SCF.h:237
void update_subspace(World &world, vecfuncT &Vpsia, vecfuncT &Vpsib, tensorT &focka, tensorT &fockb, subspaceT &subspace, tensorT &Q, double &bsh_residual, double &update_residual)
Definition SCF.cc:1898
Molecule molecule
Definition SCF.h:194
void load_mos(World &world)
Definition SCF.cc:321
tensorT beps
Definition SCF.h:219
tensorT diag_fock_matrix(World &world, tensorT &fock, vecfuncT &psi, vecfuncT &Vpsi, tensorT &evals, const tensorT &occ, const double thresh) const
diagonalize the fock matrix, taking care of degenerate states
Definition SCF.cc:1796
const vecfuncT & get_amo() const
getter for the molecular orbitals, alpha spin
Definition SCF.h:314
static void analyze_vectors(World &world, const vecfuncT &mo, const vecfuncT &ao, double vtol, const Molecule &molecule, const int print_level, const AtomicBasisSet &aobasis, const tensorT &occ=tensorT(), const tensorT &energy=tensorT(), const std::vector< int > &set=std::vector< int >())
Definition SCF.cc:611
XCfunctional xc
Definition SCF.h:196
void initial_load_bal(World &world)
Definition SCF.cc:1252
std::vector< int > at_nbf
Definition SCF.h:213
double do_step_restriction(World &world, const vecfuncT &mo, vecfuncT &mo_new, std::string spin) const
perform step restriction following the KAIN solver
Definition SCF.cc:2040
std::vector< std::shared_ptr< real_derivative_3d > > gradop
Definition SCF.h:221
void set_print_timings(const bool value)
Definition SCF.cc:263
double converged_for_dconv
mos are converged for this density
Definition SCF.h:225
tensorT get_fock_transformation(World &world, const tensorT &overlap, tensorT &fock, tensorT &evals, const tensorT &occ, const double thresh_degenerate) const
compute the unitary transformation that diagonalizes the fock matrix
Definition SCF.cc:1765
bool restart_aos(World &world)
Definition SCF.cc:699
tensorT aeps
orbital energies for alpha and beta orbitals
Definition SCF.h:219
functionT make_coulomb_potential(const functionT &rho) const
make the Coulomb potential given the total density
Definition SCF.h:418
tensorT make_fock_matrix(World &world, const vecfuncT &psi, const vecfuncT &Vpsi, const tensorT &occ, double &ekinetic) const
Definition SCF.cc:1677
double converged_for_tconv
derivatives of mos are converged for this threshold
Definition SCF.h:226
Tensor< double > twoint(World &world, const vecfuncT &psi) const
Compute the two-electron integrals over the provided set of orbitals.
Definition SCF.cc:1735
tensorT aocc
occupation numbers for alpha and beta orbitals
Definition SCF.h:216
void set_protocol(World &world, double thresh)
Definition SCF.h:262
static functionT make_lda_potential(World &world, const functionT &arho)
Definition SCF.cc:1334
functionT make_density(World &world, const tensorT &occ, const vecfuncT &v) const
Definition SCF.cc:1269
void rotate_subspace(World &world, const tensorT &U, subspaceT &subspace, int lo, int nfunc, double trantol) const
Definition SCF.cc:1864
void reset_aobasis(const std::string &aobasisname)
Definition SCF.h:344
void project(World &world)
Definition SCF.cc:550
double converged_for_thresh
mos are converged for this threshold
Definition SCF.h:224
CalculationParameters param
Definition SCF.h:195
double current_energy
Definition SCF.h:223
void loadbal(World &world, functionT &arho, functionT &brho, functionT &arho_old, functionT &brho_old, subspaceT &subspace)
Definition SCF.cc:1829
void vector_stats(const std::vector< double > &v, double &rms, double &maxabsval) const
Definition SCF.cc:1531
vecfuncT project_ao_basis(World &world, const AtomicBasisSet &aobasis)
Definition SCF.cc:589
std::vector< int > group_orbital_sets(World &world, const tensorT &eps, const tensorT &occ, const int nmo) const
group orbitals into sets of similar orbital energies for localization
Definition SCF.cc:1227
void initial_guess_from_nwchem(World &world)
Definition SCF.cc:752
functionT mask
Definition SCF.h:199
std::shared_ptr< PotentialManager > potentialmanager
Definition SCF.h:192
static std::vector< poperatorT > make_bsh_operators(World &world, const tensorT &evals, const CalculationParameters &param)
Definition SCF.cc:1310
const vecfuncT & get_bmo() const
getter for the molecular orbitals, beta spin
Definition SCF.h:317
vecfuncT bmo
Definition SCF.h:204
void solve(World &world)
Definition SCF.cc:2139
void orthonormalize(World &world, vecfuncT &amo_new) const
orthonormalize the vectors
Definition SCF.cc:2112
Convolutions in separated form (including Gaussian)
Definition operator.h:139
A tensor is a multidimensional array.
Definition tensor.h:317
Tensor< T > reshape(int ndimnew, const long *d)
Returns new view/tensor reshaping size/number of dimensions to conforming tensor.
Definition tensor.h:1384
T sumsq() const
Returns the sum of the squares of the elements.
Definition tensor.h:1669
Tensor< T > flat()
Returns new view/tensor rehshaping to flat (1-d) tensor.
Definition tensor.h:1555
A simple, fixed dimension vector.
Definition vector.h:64
void fence(bool debug=false)
Synchronizes all processes in communicator AND globally ensures no pending AM or tasks.
Definition worldgop.cc:161
A parallel world class.
Definition world.h:132
ProcessID rank() const
Returns the process rank in this World (same as MPI_Comm_rank()).
Definition world.h:320
WorldGopInterface & gop
Global operations.
Definition world.h:207
Simplified interface to XC functionals.
Definition xcfunctional.h:48
Definition SCF.h:170
std::map< std::string, std::vector< double > > e_data
Definition SCF.h:172
void add_gradient(const Tensor< double > &grad)
Definition SCF.cc:229
int iter
Definition SCF.h:175
void to_json(json &j) const
Definition SCF.cc:212
void add_data(std::map< std::string, double > values)
Definition SCF.cc:189
json hessian
Definition SCF.h:174
json gradient
Definition SCF.h:173
scf_data()
Definition SCF.cc:199
void print_data()
Definition SCF.cc:225
Declaration of core potential related class.
double(* energy)()
Definition derivatives.cc:58
char * p(char *buf, const char *name, int k, int initial_level, double thresh, int order)
Definition derivatives.cc:72
static double lo
Definition dirac-hatom.cc:23
double psi(const Vector< double, 3 > &r)
Definition hatom_energy.cc:78
static const double v
Definition hatom_sf_dirac.cc:20
Main include file for MADNESS and defines Function interface.
Namespace for all elements and tools of MADNESS.
Definition DFParameters.h:10
void print_header2(const std::string &s)
medium section heading
Definition print.cc:54
Vector< double, 3 > coordT
Definition corepotential.cc:54
static SeparatedConvolution< double, 3 > * CoulombOperatorPtr(World &world, double lo, double eps, const array_of_bools< 3 > &lattice_sum=FunctionDefaults< 3 >::get_bc().is_periodic(), int k=FunctionDefaults< 3 >::get_k())
Factory function generating separated kernel for convolution with 1/r in 3D.
Definition operator.h:1818
nlohmann::json json
Definition QCCalculationParametersBase.h:29
static void user_to_sim(const Vector< double, NDIM > &xuser, Vector< double, NDIM > &xsim)
Convert user coords (cell[][]) to simulation coords ([0,1]^ndim)
Definition funcdefaults.h:444
Tensor< double > tensorT
Definition distpm.cc:21
Function< std::complex< double >, 3 > complex_functionT
Definition SCF.h:79
std::vector< pairvecfuncT > subspaceT
Definition SCF.h:73
double mask1(double x)
Definition SCF.h:102
DistributedMatrix< double > distmatT
Definition SCF.h:75
std::pair< vecfuncT, vecfuncT > pairvecfuncT
Definition SCF.h:72
void print(const T &t, const Ts &... ts)
Print items to std::cout (items separated by spaces) and terminate with a new line.
Definition print.h:226
FunctionFactory< double, 3 > factoryT
Definition corepotential.cc:57
response_space apply(World &world, std::vector< std::vector< std::shared_ptr< real_convolution_3d > > > &op, response_space &f)
Definition basic_operators.cc:43
std::shared_ptr< operatorT> poperatorT
Definition SCF.h:78
Function< double, 3 > functionT
Definition corepotential.cc:56
NDIM & f
Definition mra.h:2528
void to_json(nlohmann::json &j)
std::shared_ptr< FunctionFunctorInterface< double, 3 > > functorT
Definition corepotential.cc:55
std::vector< complex_functionT > cvecfuncT
Definition SCF.h:80
vector< functionT > vecfuncT
Definition corepotential.cc:58
static double mask3(const coordT &ruser)
Definition SCF.h:116
std::shared_ptr< WorldDCPmapInterface< Key< 3 > > > pmapT
Definition SCF.h:67
std::vector< Function< T, NDIM > > grad(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
shorthand gradient operator
Definition vmra.h:2034
static XNonlinearSolver< std::vector< Function< T, NDIM > >, T, vector_function_allocator< T, NDIM > > nonlinear_vector_solver(World &world, const long nvec)
Definition nonlinsol.h:371
SeparatedConvolution< double, 3 > operatorT
Definition SCF.h:77
static const size_t nfunc
Definition pcr.cc:63
Declaration of molecule-related classes and functions.
double Q(double a)
Definition relops.cc:20
static const double thresh
Definition rk.cc:45
static const long k
Definition rk.cc:44
const double xi
Exponent for delta function approx.
Definition siam_example.cc:60
Defines interfaces for optimization and non-linear equation solvers.
std::string prefix
Definition tdse.cc:71
Definition CalculationParameters.h:51
std::vector< double > protocol() const
Definition CalculationParameters.h:210
int nv_factor() const
Definition CalculationParameters.h:175
int k() const
Definition CalculationParameters.h:187
int nmo_alpha() const
Definition CalculationParameters.h:177
bool save() const
Definition CalculationParameters.h:211
double dconv() const
Definition CalculationParameters.h:150
double econv() const
Definition CalculationParameters.h:149
int print_level() const
Definition CalculationParameters.h:196
std::string deriv() const
Definition CalculationParameters.h:201
double L() const
Definition CalculationParameters.h:186
std::string aobasis() const
Definition CalculationParameters.h:208
int nalpha() const
Definition CalculationParameters.h:170
int nbeta() const
Definition CalculationParameters.h:171
int nmo_beta() const
Definition CalculationParameters.h:178
bool no_compute() const
Definition CalculationParameters.h:183
bool spin_restricted() const
Definition CalculationParameters.h:182
double lo() const
Definition CalculationParameters.h:185
std::string pcm_data() const
Definition CalculationParameters.h:203
Definition convolution1d.h:992
double eprec() const
Definition molecule.h:258
The interface to be provided by functions to be optimized.
Definition solvers.h:176
very simple command line parser
Definition commandlineparser.h:15
Definition SCF.h:84
double parent_value
Definition SCF.h:86
lbcost(double leaf_value=1.0, double parent_value=0.0)
Definition SCF.h:88
double leaf_value
Definition SCF.h:85
double operator()(const Key< NDIM > &key, const FunctionNode< T, NDIM > &node) const
Definition SCF.h:90
Class to compute the energy functional.
Definition xcfunctional.h:365
InputParameters param
Definition tdse.cc:203
constexpr std::size_t NDIM
Definition testgconv.cc:54
static Molecule molecule
Definition testperiodicdft.cc:39
static Subspace * subspace
Definition testperiodicdft.cc:41