| MADNESS 0.10.1
    | 
#include <nemo.h>


| Classes | |
| struct | NemoCalculationParameters | 
| class holding parameters for a nemo calculation beyond the standard dft parameters from moldft  More... | |
| Public Member Functions | |
| Nemo (World &world, const CalculationParameters ¶m, const NemoCalculationParameters &nemo_param, const Molecule &molecule) | |
| Nemo (World &world, const commandlineparser &parser) | |
| ctor | |
| virtual nlohmann::json | analyze () const | 
| compute dipole moment and gradient at the current geometry | |
| bool | check_converged (const Tensor< double > &x) const | 
| std::vector< vecfuncT > | compute_all_cphf () | 
| solve the CPHF equation for all displacements | |
| vecfuncT | compute_cphf_parallel_term (const size_t iatom, const int iaxis) const | 
| tensorT | compute_fock_matrix (const vecfuncT &nemo, const tensorT &occ) const | 
| compute the Fock matrix from scratch | |
| Tensor< double > | compute_IR_intensities (const Tensor< double > &normalmodes, const vecfuncT &dens_pt) const | 
| compute the IR intensities in the double harmonic approximation | |
| bool | do_ac () const | 
| bool | do_pcm () const | 
| bool | do_symmetry () const | 
| AC< 3 > | get_ac () const | 
| std::shared_ptr< SCF > | get_calc () const | 
| const CalculationParameters & | get_calc_param () const | 
| const NemoCalculationParameters & | get_nemo_param () const | 
| PCM | get_pcm () const | 
| projector_irrep | get_symmetry_projector () const | 
| return the symmetry_projector | |
| Tensor< double > | gradient (const Tensor< double > &x) | 
| compute the nuclear gradients | |
| Tensor< double > | hessian (const Tensor< double > &x) | 
| returns the molecular hessian matrix at structure x | |
| bool | is_dft () const | 
| real_function_3d | kinetic_energy_potential (const vecfuncT &nemo) const | 
| void | load_mos (World &w) | 
| real_function_3d | make_ddensity (const real_function_3d &rhonemo, const int axis) const | 
| make the derivative of the density | |
| real_function_3d | make_density (const Tensor< double > &occ, const vecfuncT &nemo) const | 
| make the density (alpha or beta) | |
| real_function_3d | make_density (const tensorT &occ, const vecfuncT &bra, const vecfuncT &ket, const bool refine=false) const | 
| make the density using different bra and ket vectors | |
| virtual std::shared_ptr< Fock< double, 3 > > | make_fock_operator () const | 
| construct the fock operator based on the calculation parameters (K or XC?) | |
| real_function_3d | make_laplacian_density (const real_function_3d &rhonemo) const | 
| the Laplacian of the density | |
| real_function_3d | make_sigma (const real_function_3d &rho1, const real_function_3d &rho2) const | 
| compute the reduced densities sigma (gamma) for GGA functionals | |
| Molecule & | molecule () | 
| return a reference to the molecule | |
| Molecule & | molecule () const | 
| return a reference to the molecule | |
| std::string | name () const | 
| bool | provides_gradient () const | 
| Override this to return true if the derivative is implemented. | |
| Tensor< double > | purify_hessian (const Tensor< double > &hessian) const | 
| purify and symmetrize the hessian | |
| bool | selftest () | 
| vecfuncT | solve_cphf (const size_t iatom, const int iaxis, const Tensor< double > fock, const vecfuncT &guess, const vecfuncT &rhsconst, const Tensor< double > incomplete_hessian, const vecfuncT ¶llel, const SCFProtocol &p, const std::string &xc_data) const | 
| solve the CPHF equations for the nuclear displacements | |
| virtual double | value () | 
| virtual double | value (const Tensor< double > &x) | 
| Should return the value of the objective function. | |
|  Public Member Functions inherited from madness::NemoBase | |
| NemoBase (World &w) | |
| virtual | ~NemoBase () | 
| bool | check_convergence (const std::vector< double > energies, const std::vector< double > oldenergies, const double bsh_norm, const double delta_density, const CalculationParameters ¶m, const double econv, const double dconv) const | 
| template<typename T , std::size_t NDIM> | |
| Function< typename Tensor< T >::scalar_type, NDIM > | compute_density (const std::vector< Function< T, NDIM > > nemo) const | 
| Tensor< double > | compute_gradient (const real_function_3d &rhonemo, const Molecule &molecule) const | 
| compute the nuclear gradients | |
| template<typename T , std::size_t NDIM> | |
| double | compute_kinetic_energy (const std::vector< Function< T, NDIM > > &nemo) const | 
| compute kinetic energy as square of the "analytical" expectation value | |
| template<typename T , std::size_t NDIM> | |
| double | compute_kinetic_energy1 (const std::vector< Function< T, NDIM > > &nemo) const | 
| compute kinetic energy as square of the "analytical" derivative of the orbitals | |
| template<typename T , std::size_t NDIM> | |
| double | compute_kinetic_energy1a (const std::vector< Function< T, NDIM > > &nemo) const | 
| compute kinetic energy as square of the "analytical" derivative of the orbitals | |
| template<typename T , std::size_t NDIM> | |
| double | compute_kinetic_energy2 (const std::vector< Function< T, NDIM > > &nemo) const | 
| compute kinetic energy as direct derivative of the orbitals (probably imprecise) | |
| void | construct_nuclear_correlation_factor (const Molecule &molecule, const std::shared_ptr< PotentialManager > pm, const std::pair< std::string, double > ncf_parameter) | 
| std::shared_ptr< NuclearCorrelationFactor > | get_ncf_ptr () const | 
| create an instance of the derived object based on the input parameters | |
| virtual void | invalidate_factors_and_potentials () | 
| virtual bool | need_recompute_factors_and_potentials (const double thresh) const | 
| template<typename T , std::size_t NDIM> | |
| void | orthonormalize (std::vector< Function< T, NDIM > > &nemo, const Function< double, NDIM > metric=Function< double, NDIM >(), const double trantol=FunctionDefaults< NDIM >::get_thresh() *0.01) const | 
| orthonormalize the vectors | |
|  Public Member Functions inherited from madness::OptimizationTargetInterface | |
| virtual | ~OptimizationTargetInterface () | 
| double | test_gradient (Tensor< double > &x, double value_precision, bool doprint=true) | 
| Numerical test of the derivative ... optionally prints to stdout, returns max abs error. | |
| virtual void | value_and_gradient (const Tensor< double > &x, double &value, Tensor< double > &gradient) | 
| Reimplement if more efficient to evaluate both value and gradient in one call. | |
| Static Public Member Functions | |
| static void | help () | 
| static void | print_parameters () | 
| static void | smoothen (real_function_3d &f) | 
| smooth a function by projecting it onto k-1 and then average with k | |
|  Static Public Member Functions inherited from madness::NemoBase | |
| template<typename T , std::size_t NDIM> | |
| static void | normalize (std::vector< Function< T, NDIM > > &nemo, const Function< double, NDIM > metric=Function< double, NDIM >()) | 
| normalize the nemos | |
| template<typename T > | |
| static Tensor< T > | Q2 (const Tensor< T > &s) | 
| Public Attributes | |
| const NemoCalculationParameters | nemo_param | 
|  Public Attributes inherited from madness::NemoBase | |
| std::shared_ptr< NuclearCorrelationFactor > | ncf | 
| the nuclear correlation factor | |
| real_function_3d | R | 
| the nuclear correlation factor | |
| real_function_3d | R_square | 
| the square of the nuclear correlation factor | |
| World & | world | 
| Protected Member Functions | |
| std::vector< double > | compute_energy_regularized (const vecfuncT &nemo, const vecfuncT &Jnemo, const vecfuncT &Knemo, const vecfuncT &Unemo) const | 
| given nemos, compute the HF energy using the regularized expressions for T and V | |
| void | compute_nemo_potentials (const vecfuncT &nemo, vecfuncT &Jnemo, vecfuncT &Knemo, vecfuncT &xcnemo, vecfuncT &pcmnemo, vecfuncT &Unemo) const | 
| compute the reconstructed orbitals, and all potentials applied on nemo | |
| real_function_3d | get_coulomb_potential (const vecfuncT &psi) const | 
| return the Coulomb potential | |
| template<typename T , size_t NDIM> | |
| void | load_function (std::vector< Function< T, NDIM > > &f, const std::string name) const | 
| load a function | |
| vecfuncT | localize (const vecfuncT &nemo, const double dconv, const bool randomize) const | 
| localize the nemo orbitals | |
| vecfuncT | make_cphf_constant_term (const size_t iatom, const int iaxis, const vecfuncT &R2nemo, const real_function_3d &rhonemo) const | 
| compute the constant term for the CPHF equations | |
| Tensor< double > | make_incomplete_hessian () const | 
| compute the incomplete hessian | |
| Tensor< double > | make_incomplete_hessian_response_part (const std::vector< vecfuncT > &xi) const | 
| compute the complementary incomplete hessian | |
| void | make_plots (const real_function_3d &f, const std::string &name="function") const | 
| template<typename T , size_t NDIM> | |
| void | save_function (const std::vector< Function< T, NDIM > > &f, const std::string name) const | 
| save a function | |
| void | set_protocol (const double thresh) | 
| adapt the thresholds consistently to a common value | |
| double | solve (const SCFProtocol &proto) | 
| solve the HF equations | |
| double | trantol () const | 
| return the threshold for vanishing elements in orbital rotations | |
| Protected Attributes | |
| AC< 3 > | ac | 
| asymptotic correction for DFT | |
| std::shared_ptr< SCF > | calc | 
| std::shared_ptr< real_convolution_3d > | poisson | 
| a poisson solver | |
| projector_irrep | symmetry_projector | 
| Private Types | |
| typedef std::shared_ptr< real_convolution_3d > | poperatorT | 
| Private Attributes | |
| double | coords_sum | 
| sum of square of coords at last solved geometry | |
| PCM | pcm | 
| polarizable continuum model | |
| Friends | |
| class | PNO | 
| class | TDHF | 
The Nemo class.
| 
 | private | 
| madness::Nemo::Nemo | ( | World & | world, | 
| const commandlineparser & | parser | ||
| ) | 
ctor
| [in] | world1 | the world | 
| [in] | calc | the SCF | 
References do_pcm(), get_calc_param(), madness::projector_irrep::get_verbosity(), molecule(), pcm, madness::projector_irrep::print_character_table(), madness::projector_irrep::set_ordering(), madness::projector_irrep::set_orthonormalize_irreps(), madness::projector_irrep::set_verbosity(), symmetry_projector, and madness::NemoBase::world.
| madness::Nemo::Nemo | ( | World & | world, | 
| const CalculationParameters & | param, | ||
| const NemoCalculationParameters & | nemo_param, | ||
| const Molecule & | molecule | ||
| ) | 
References do_pcm(), get_calc_param(), madness::projector_irrep::get_verbosity(), molecule(), pcm, madness::projector_irrep::print_character_table(), madness::projector_irrep::set_ordering(), madness::projector_irrep::set_orthonormalize_irreps(), madness::projector_irrep::set_verbosity(), symmetry_projector, and madness::NemoBase::world.
| 
 | virtual | 
compute dipole moment and gradient at the current geometry
Reimplemented in madness::OEP.
References calc, madness::NemoBase::compute_gradient(), madness::PropertyResults::dipole, madness::timer::end(), madness::PropertyResults::energy, madness::grad(), madness::PropertyResults::gradient, make_density(), molecule(), madness::NemoBase::R_square, madness::PropertyResults::to_json(), and madness::NemoBase::world.
Referenced by madness::OEP::solve().
| 
 | inline | 
References coords_sum, and madness::Tensor< T >::sumsq().
Referenced by value().
| std::vector< vecfuncT > madness::Nemo::compute_all_cphf | ( | ) | 
solve the CPHF equation for all displacements
this function computes the nemo response F^X
![\[
    F^X = F^\perp + F^\parallel
\]](form_293.png)
To reconstruct the unregularized orbital response (not recommended):
![\[
  i^X   = R^X F + R F^X
\]](form_311.png)
The orbital response i^X constructed in this way is automatically orthogonal to the occupied space because of the parallel term F^\parallel
References madness::apply(), calc, compute_cphf_parallel_term(), compute_fock_matrix(), madness::SCFProtocol::current_prec, madness::dot(), madness::SCFProtocol::end_prec, madness::SCFProtocol::finished(), get_calc(), get_calc_param(), madness::SCFProtocol::initialize(), load_function(), make_cphf_constant_term(), make_density(), make_incomplete_hessian(), make_incomplete_hessian_response_part(), molecule(), madness::mul(), madness::Molecule::natom(), p(), param, madness::print(), madness::NemoBase::R_square, madness::World::rank(), madness::save(), save_function(), set_protocol(), solve_cphf(), madness::stringify(), madness::timer::tag(), madness::truncate(), madness::wall_time(), madness::NemoBase::world, xc, and xi.
Referenced by hessian().
| vecfuncT madness::Nemo::compute_cphf_parallel_term | ( | const size_t | iatom, | 
| const int | iaxis | ||
| ) | const | 
this function computes that part of the orbital response that is parallel to the occupied space.
![\[
    F^X = F^\perp + F^\parallel
\]](form_293.png)
If no NCF's are used F^\parallel vanishes. If NCF's are used this term does not vanish because the derivatives of the NCF does not vanish, and it is given by
![\[
 F_i^\parallel = -\frac{1}{2}\sum_k|F_k ><F_k | (R^2)^X | F_i>
\]](form_312.png)
References calc, madness::matrix_inner(), madness::mul(), madness::transform(), madness::truncate(), and madness::NemoBase::world.
Referenced by compute_all_cphf().
| 
 | protected | 
given nemos, compute the HF energy using the regularized expressions for T and V
References madness::apply(), axis, calc, madness::NemoBase::compute_kinetic_energy(), madness::PCM::compute_pcm_energy(), madness::XCOperator< T, NDIM >::compute_xc_energy(), do_pcm(), madness::timer::end(), energy, get_calc_param(), madness::inner(), is_dft(), K(), pcm, madness::NemoBase::R_square, madness::World::rank(), madness::Tensor< T >::sum(), sum, madness::truncate(), and madness::NemoBase::world.
Referenced by solve().
compute the Fock matrix from scratch
References calc, compute_nemo_potentials(), do_pcm(), madness::matrix_inner(), madness::mul(), madness::NemoBase::R_square, T(), madness::transpose(), madness::truncate(), and madness::NemoBase::world.
Referenced by compute_all_cphf(), test_ethylene(), test_ne_boys(), and test_nemo().
| Tensor< double > madness::Nemo::compute_IR_intensities | ( | const Tensor< double > & | normalmodes, | 
| const vecfuncT & | dens_pt | ||
| ) | const | 
compute the IR intensities in the double harmonic approximation
use the projected normal modes; units are km/mol
| [in] | normalmodes | the normal modes | 
| [in] | dens_pt | the perturbed densities for each nuclear displacement | 
References madness::constants::atomic_mass_constant, madness::constants::atomic_mass_in_au, madness::constants::atomic_unit_of_electric_dipole_moment, madness::constants::atomic_unit_of_length, madness::constants::Avogadro_constant, c, madness::constants::dielectric_constant, e(), madness::Tensor< T >::emul(), madness::inner(), madness::Molecule::massweights(), molecule(), N, madness::Molecule::natom(), madness::Molecule::nuclear_dipole_derivative(), pi, madness::constants::pi, madness::MolecularOptimizer::projector_external_dof(), madness::Tensor< T >::scale(), and madness::constants::speed_of_light_in_vacuum.
Referenced by hessian().
| 
 | protected | 
compute the reconstructed orbitals, and all potentials applied on nemo
to use these potentials in the fock matrix computation they must be multiplied by the nuclear correlation factor
| [in] | nemo | the nemo orbitals | 
| [out] | Jnemo | Coulomb operator applied on the nemos | 
| [out] | Knemo | exchange operator applied on the nemos | 
| [out] | pcmnemo | PCM (solvent) potential applied on the nemos | 
| [out] | Unemo | regularized nuclear potential applied on the nemos | 
to use these potentials in the fock matrix computation they must be multiplied by the nuclear correlation factor
| [in] | nemo | the nemo orbitals | 
| [out] | psi | the reconstructed, full orbitals | 
| [out] | Jnemo | Coulomb operator applied on the nemos | 
| [out] | Knemo | exchange operator applied on the nemos | 
| [out] | Vnemo | nuclear potential applied on the nemos | 
| [out] | Unemo | regularized nuclear potential applied on the nemos | 
References ac, madness::AC< NDIM >::apply(), calc, madness::CalculationParameters::charge(), madness::PCM::compute_pcm_potential(), do_ac(), do_pcm(), madness::WorldGopInterface::fence(), get_calc_param(), madness::get_size(), madness::World::gop, K(), madness::XCOperator< T, NDIM >::make_xc_potential(), molecule(), pcm, madness::Coulomb< T, NDIM >::potential(), madness::print(), madness::World::rank(), madness::scale(), madness::Exchange< T, NDIM >::set_algorithm(), madness::Exchange< T, NDIM >::set_symmetric(), madness::Exchange< T, NDIM >::set_taskq(), madness::Coulomb< T, NDIM >::set_taskq(), madness::World::size(), madness::stringify(), madness::timer::tag(), madness::truncate(), and madness::NemoBase::world.
Referenced by compute_fock_matrix(), and solve().
| 
 | inline | 
References madness::CalculationParameters::ac_data(), and get_calc_param().
Referenced by compute_nemo_potentials(), and make_fock_operator().
| 
 | inline | 
References get_calc_param(), and madness::CalculationParameters::pcm_data().
Referenced by Nemo(), Nemo(), compute_energy_regularized(), compute_fock_matrix(), compute_nemo_potentials(), make_fock_operator(), and solve().
| 
 | inline | 
References madness::projector_irrep::get_pointgroup(), and symmetry_projector.
Referenced by madness::OEP::iterate(), and solve().
| 
 | inline | 
References calc.
Referenced by madness::Coulomb< T, NDIM >::Coulomb(), madness::Exchange< T, NDIM >::ExchangeImpl< T, NDIM >::ExchangeImpl(), madness::F12Potentials::F12Potentials(), madness::PNO::PNO(), madness::XCOperator< T, NDIM >::XCOperator(), madness::XCOperator< T, NDIM >::XCOperator(), compute_all_cphf(), madness::OEP::compute_exchange_potential(), madness::Coulomb< T, NDIM >::compute_potential(), madness::PNO::grow_rank(), madness::PNO::guess_virtuals(), hessian(), madness::F12Potentials::initialize_active_mos(), madness::PNO::initialize_pairs(), madness::OEP::iterate(), madness::PNO::iterate_pairs_internal(), kinetic_energy_potential(), madness::OEP::load_restartdata(), localize(), main(), madness::PNO::make_bsh_operators(), make_cphf_constant_term(), make_incomplete_hessian(), make_laplacian_density(), madness::PNO::nact(), madness::PNO::nocc(), madness::PNO::oit(), madness::ParametrizedExchange::operator()(), madness::F12Potentials::read_cabs_from_file(), madness::OEP::save_restartdata(), set_protocol(), solve(), solve_cphf(), test_ethylene(), test_fock(), test_ne_boys(), test_nemo(), and trantol().
| 
 | inline | 
References calc.
Referenced by Nemo(), Nemo(), madness::OEP::OEP(), madness::OEP::OEP(), compute_all_cphf(), madness::OEP::compute_and_print_final_energies(), compute_energy_regularized(), compute_nemo_potentials(), madness::OEP::compute_Pauli_kinetic_density(), madness::OEP::compute_slater_potential(), madness::OEP::compute_total_kinetic_density(), do_ac(), do_pcm(), get_coulomb_potential(), hessian(), madness::OEP::iterate(), make_cphf_constant_term(), make_fock_operator(), madness::OEP::output_calc_info_schema(), madness::OEP::print_parameters(), madness::OEP::selftest(), set_protocol(), solve(), solve_cphf(), value(), and madness::OEP::value().
| 
 | protected | 
return the Coulomb potential
References calc, get_calc_param(), MADNESS_ASSERT, make_density(), psi(), and madness::Function< T, NDIM >::scale().
| 
 | inline | 
References nemo_param.
Referenced by hessian(), and set_protocol().
| 
 | inline | 
return the symmetry_projector
References symmetry_projector.
compute the nuclear gradients
Reimplemented from madness::OptimizationTargetInterface.
References madness::Atom::atomic_number, madness::atomic_number_to_symbol(), calc, madness::NemoBase::compute_gradient(), madness::grad(), make_density(), madness::print(), madness::World::rank(), madness::Function< T, NDIM >::refine(), madness::Function< T, NDIM >::scale(), madness::NemoBase::world, madness::Atom::x, madness::Atom::y, and madness::Atom::z.
| 
 | inlinestatic | 
References madness::print(), and madness::print_header2().
Referenced by main().
returns the molecular hessian matrix at structure x
compute the nuclear hessian
References madness::Tensor< T >::absmax(), madness::constants::au2invcm, calc, compute_all_cphf(), compute_frequencies(), compute_IR_intensities(), compute_reduced_mass(), madness::timer::end(), get_calc(), get_calc_param(), get_nemo_param(), hessian(), madness::inner(), make_ddensity(), make_density(), molecule(), madness::mul(), madness::Molecule::natom(), madness::Molecule::nuclear_repulsion_hessian(), madness::print(), purify_hessian(), madness::NemoBase::R_square, madness::World::rank(), save_function(), madness::Tensor< T >::scale(), madness::BaseTensor::size(), sum, madness::transpose(), madness::truncate(), madness::wall_time(), madness::NemoBase::world, and xi.
Referenced by hessian(), and purify_hessian().
| 
 | inline | 
References calc.
Referenced by compute_energy_regularized(), make_cphf_constant_term(), and solve_cphf().
| real_function_3d madness::Nemo::kinetic_energy_potential | ( | const vecfuncT & | nemo | ) | const | 
compute the kinetic energy potential using Eq. (16) of R. A. King and N. C. Handy, “Kinetic energy functionals from the Kohn–Sham potential,” Phys. Chem. Chem. Phys., vol. 2, no. 22, pp. 5049–5056, 2000.
References madness::apply(), madness::binary_op(), calc, diff(), madness::dot(), madness::gaxpy(), get_calc(), madness::mul(), madness::norm2(), madness::constants::pi, poisson, madness::print(), madness::NemoBase::R, madness::NemoBase::R_square, madness::save(), madness::scale(), madness::SmoothingOperator3D(), madness::sub(), sum, and madness::NemoBase::world.
| 
 | protected | 
load a function
References madness::f, name(), madness::print(), madness::World::rank(), and madness::NemoBase::world.
Referenced by compute_all_cphf().
| 
 | inline | 
Referenced by madness::OEP::value().
| 
 | protected | 
localize the nemo orbitals
localize the nemo orbitals according to Pipek-Mezey or Foster-Boys
References aobasis, calc, madness::Localizer::compute_localization_matrix(), get_calc(), molecule(), madness::NemoBase::normalize(), madness::NemoBase::R, madness::Localizer::set_method(), madness::Localizer::set_metric(), madness::transform(), madness::transpose(), madness::truncate(), and madness::NemoBase::world.
Referenced by madness::OEP::iterate(), and solve().
| 
 | protected | 
compute the constant term for the CPHF equations
mainly all terms with the nuclear correlation factor's derivatives
References calc, madness::Coulomb< T, NDIM >::compute_potential(), get_calc(), get_calc_param(), is_dft(), lo, madness::Coulomb< T, NDIM >::potential(), Q(), madness::Exchange< T, NDIM >::set_bra_and_ket(), madness::BaseTensor::size(), madness::truncate(), and madness::NemoBase::world.
Referenced by compute_all_cphf().
| real_function_3d madness::Nemo::make_ddensity | ( | const real_function_3d & | rhonemo, | 
| const int | axis | ||
| ) | const | 
make the derivative of the density

| [in] | rhonemo | the regularized density | 
| [in] | axis | the component of the nabla operator | 
References axis, madness::copy(), madness::NemoBase::R_square, madness::Function< T, NDIM >::refine(), and madness::NemoBase::world.
Referenced by hessian().
| real_function_3d madness::Nemo::make_density | ( | const Tensor< double > & | occ, | 
| const vecfuncT & | nemo | ||
| ) | const | 
make the density (alpha or beta)
References calc, and madness::NemoBase::world.
Referenced by madness::XCOperator< T, NDIM >::XCOperator(), analyze(), compute_all_cphf(), madness::Coulomb< T, NDIM >::compute_potential(), get_coulomb_potential(), gradient(), hessian(), make_incomplete_hessian(), madness::ParametrizedExchange::operator()(), solve_cphf(), and value().
| real_function_3d madness::Nemo::make_density | ( | const tensorT & | occ, | 
| const vecfuncT & | bra, | ||
| const vecfuncT & | ket, | ||
| const bool | refine = false | ||
| ) | const | 
make the density using different bra and ket vectors
e.g. for computing the perturbed density \sum_i \phi_i \phi_i^X or when using nemos: \sum_i R2nemo_i nemo_i
References calc, madness::compress(), madness::Function< T, NDIM >::compressed, madness::WorldGopInterface::fence(), madness::Function< T, NDIM >::gaxpy(), madness::World::gop, madness::mul(), madness::refine(), and madness::NemoBase::world.
| 
 | virtual | 
construct the fock operator based on the calculation parameters (K or XC?)
Reimplemented from madness::NemoBase.
Reimplemented in madness::OEP.
References ac, madness::AC< NDIM >::apply(), calc, madness::CalculationParameters::charge(), madness::PCM::compute_pcm_potential(), do_ac(), do_pcm(), get_calc_param(), K(), MADNESS_CHECK, madness::XCOperator< T, NDIM >::make_xc_potential(), molecule(), pcm, madness::Coulomb< T, NDIM >::potential(), madness::LocalPotentialOperator< T, NDIM >::set_info(), madness::LocalPotentialOperator< T, NDIM >::set_potential(), madness::Exchange< T, NDIM >::set_symmetric(), and madness::NemoBase::world.
Referenced by madness::Fock< T, NDIM >::Fock(), and test_nemo().
| 
 | protected | 
compute the incomplete hessian
incomplete hessian is the nuclear-nuclear contribution, and the contribution from the second derivative of the nuclear potential, and also the derivative of the nuclear correlation factor. i.e. all contributions that do not contain the regularized perturbed density, but it will contain parts of the perturbed density
References get_calc(), madness::inner(), make_density(), molecule(), madness::Molecule::natom(), madness::Molecule::nuclear_repulsion_hessian(), madness::NemoBase::R_square, madness::refine(), and madness::NemoBase::world.
Referenced by compute_all_cphf().
| 
 | protected | 
compute the complementary incomplete hessian
| [in] | xi | the response functions including the parallel part | 
References madness::_(), calc, madness::dot(), h(), madness::inner(), molecule(), madness::Molecule::natom(), madness::NemoBase::R_square, madness::NemoBase::world, and xi.
Referenced by compute_all_cphf().
| real_function_3d madness::Nemo::make_laplacian_density | ( | const real_function_3d & | rhonemo | ) | const | 
the Laplacian of the density
The Laplacian should currently only be used for subsequent convolution with a Green's function (which is reasonably stable), but not on its own!
The Laplacian of the cuspy density is numerically fairly unstable:
It may be computed as
![\[
  \Delta \rho = \Delta (R^2 \rho_R)
         = \Delta (R^2) \rho_R + 2\nabla R \nabla \rho_R + R^2 \Delta \rho_R
         = 2 R^2 U1^2 \rho_R -4 R^2 ( U-V ) \rho_R + R^2 \Delta\rho_R
\]](form_314.png)
where we can use the identity
![\[
  U=V + R^{-1}[T,R]
  -2 R (U-V) = \Delta R + 2\nabla R\dot \nabla
\]](form_315.png)
first term comes from the definition of the U potential as the commutator over the kinetic energy (aka the Laplacian)
| [in] | rhonemo | the regularized density \rho_R | 
References axis, madness::Function< T, NDIM >::compressed, madness::copy(), get_calc(), madness::constants::pi, poisson, madness::NemoBase::R_square, madness::Function< T, NDIM >::refine(), madness::save(), smoothen(), madness::truncate(), and madness::NemoBase::world.
| 
 | inlineprotected | 
| real_function_3d madness::Nemo::make_sigma | ( | const real_function_3d & | rho1, | 
| const real_function_3d & | rho2 | ||
| ) | const | 
compute the reduced densities sigma (gamma) for GGA functionals
the reduced density is given by
![\[
  \sigma = \nabla\rho_1 \nabla\rho_2
         = (\nabla R^2 \rho_{R,1} + R^2 \nabla \rho_{R,1}) \cdot
             (\nabla R^2 \rho_{R,2} + R^2 \nabla \rho_{R,2})
         = 4R^4 U_1^2 \rho_{R,1}\rho_{R,2}
            + 2R^4 \vec U_1 \left(\rho_{R,1}\nabla \rho_{R,2} + \nabla\rho_{R,1} \rho_{R,2}\right)
            + R^4 \nabla \rho_{R,1}\cdot \nabla\rho_{R,2}
\]](form_310.png)
References madness::dot(), madness::FunctionDefaults< NDIM >::get_thresh(), madness::grad(), madness::NemoBase::R_square, madness::FunctionDefaults< NDIM >::set_thresh(), thresh, and madness::NemoBase::world.
| 
 | inlinevirtual | 
return a reference to the molecule
Reimplemented from madness::MolecularOptimizationTargetInterface.
References calc.
Referenced by Nemo(), Nemo(), analyze(), compute_all_cphf(), compute_IR_intensities(), compute_nemo_potentials(), hessian(), localize(), make_fock_operator(), make_incomplete_hessian(), make_incomplete_hessian_response_part(), solve_cphf(), test_ethylene(), and test_ne_boys().
| 
 | inline | 
return a reference to the molecule
References calc.
| 
 | inlinevirtual | 
Implements madness::QCPropertyInterface.
Reimplemented in madness::OEP.
Referenced by load_function(), make_plots(), and save_function().
| 
 | inlinestatic | 
References param, madness::print(), and madness::Molecule::print_parameters().
Referenced by main().
| 
 | inlinevirtual | 
Override this to return true if the derivative is implemented.
Reimplemented from madness::OptimizationTargetInterface.
purify and symmetrize the hessian
The hessian should be symmetric, but it is not, because
![\[
 \langle i^{Y_B}|H^{X_A}|i\rangle \neq \langle i|H^{X_A}|i^{Y_B}\rangle
\]](form_292.png)
does holds analytically, but not numerically. If the two numbers differ, pick the more trustworthy, which is the one with a heavy atom causing the perturbed density and the light atom being the nuclear singularity.
| [in] | hessian | the raw hessian | 
References calc, madness::copy(), diff(), hessian(), and madness::print().
Referenced by hessian().
| 
 | protected | 
save a function
References madness::f, name(), madness::print(), madness::World::rank(), and madness::NemoBase::world.
Referenced by compute_all_cphf(), and hessian().
| 
 | inlinevirtual | 
Implements madness::QCPropertyInterface.
Reimplemented in madness::OEP.
| 
 | inlineprotected | 
adapt the thresholds consistently to a common value
References calc, madness::NemoBase::construct_nuclear_correlation_factor(), madness::CoulombOperatorPtr(), madness::timer::end(), get_calc(), get_calc_param(), get_nemo_param(), lo, madness::NemoBase::ncf, madness::NemoBase::need_recompute_factors_and_potentials(), poisson, madness::set_thresh(), thresh, and madness::NemoBase::world.
Referenced by compute_all_cphf(), madness::OEP::selftest(), value(), and madness::OEP::value().
| 
 | inlinestatic | 
smooth a function by projecting it onto k-1 and then average with k
kept it here for further testing
References madness::f, k, and madness::project().
Referenced by make_laplacian_density().
| 
 | protected | 
solve the HF equations
References madness::Tensor< T >::absmax(), calc, madness::NemoBase::check_convergence(), madness::NemoBase::compute_density(), compute_energy_regularized(), compute_nemo_potentials(), madness::copy(), madness::CalculationParameters::dconv(), madness::SCFProtocol::dconv, madness::QCPropertyInterface::density(), madness::BaseTensor::dim(), madness::CalculationParameters::do_localize(), do_pcm(), do_symmetry(), madness::CalculationParameters::econv(), madness::SCFProtocol::econv, madness::timer::end(), energy, get_calc(), get_calc_param(), madness::BSHApply< T, NDIM >::levelshift, madness::BSHApply< T, NDIM >::lo, localize(), madness::matrix_inner(), madness::CalculationParameters::maxiter(), madness::BSHApply< T, NDIM >::metric, madness::mul(), madness::norm2(), madness::NemoBase::normalize(), madness::CalculationParameters::orbitalshift(), madness::NemoBase::orthonormalize(), madness::print(), madness::NemoBase::R, madness::NemoBase::R_square, madness::World::rank(), residual(), madness::BSHApply< T, NDIM >::ret_value, madness::save(), madness::World::size(), symmetry_projector, T(), madness::timer::tag(), madness::transform(), trantol(), madness::truncate(), update(), madness::wall_time(), and madness::NemoBase::world.
Referenced by value().
| vecfuncT madness::Nemo::solve_cphf | ( | const size_t | iatom, | 
| const int | iaxis, | ||
| const Tensor< double > | fock, | ||
| const vecfuncT & | guess, | ||
| const vecfuncT & | rhsconst, | ||
| const Tensor< double > | incomplete_hessian, | ||
| const vecfuncT & | parallel, | ||
| const SCFProtocol & | proto, | ||
| const std::string & | xc_data | ||
| ) | const | 
solve the CPHF equations for the nuclear displacements
solve the CPHF equation for the nuclear displacements
this function computes that part of the orbital response that is orthogonal to the occupied space. If no NCF's are used this corresponds to the normal response. If NCF's are used the part parallel to the occupied space must be added!
![\[
    F^X = F^\perp + F^\parallel
\]](form_293.png)
cf parallel_CPHF()
| [in] | iatom | the atom A to be moved | 
| [in] | iaxis | the coordinate X of iatom to be moved | 
| [in] | iatom | the atom A to be moved | 
| [in] | iaxis | the coordinate X of iatom to be moved | 
References madness::_(), madness::Tensor< T >::absmax(), madness::apply(), b, calc, madness::Coulomb< T, NDIM >::compute_potential(), madness::copy(), madness::SCFProtocol::dconv, madness::dot(), madness::gaxpy(), get_calc(), get_calc_param(), guess(), h(), madness::inner(), is_dft(), K(), lo, make_density(), molecule(), madness::mul(), madness::Molecule::natom(), norm(), madness::norm2(), madness::norm2s(), param, madness::Coulomb< T, NDIM >::potential(), madness::print(), Q(), madness::NemoBase::R_square, madness::World::rank(), residual(), madness::scale(), madness::Exchange< T, NDIM >::set_bra_and_ket(), madness::Exchange< T, NDIM >::set_symmetric(), madness::World::size(), madness::stringify(), madness::transform(), trantol(), madness::truncate(), V(), madness::NemoBase::world, xc, and xi.
Referenced by compute_all_cphf().
| 
 | inlineprotected | 
return the threshold for vanishing elements in orbital rotations
References calc, and get_calc().
Referenced by solve(), and solve_cphf().
| 
 | inlinevirtual | 
Reimplemented in madness::OEP.
Referenced by main(), test_fock(), test_nemo(), and value().
| 
 | virtual | 
Should return the value of the objective function.
Implements madness::OptimizationTargetInterface.
Reimplemented in madness::OEP.
References calc, check_converged(), coords_sum, get_calc_param(), madness::NemoBase::invalidate_factors_and_potentials(), make_density(), madness::CalculationParameters::no_compute(), p(), madness::print(), madness::print_header2(), madness::NemoBase::R_square, madness::World::rank(), madness::Tensor< T >::reshape(), madness::save(), set_protocol(), madness::FunctionDefaults< NDIM >::set_thresh(), solve(), madness::stringify(), madness::Tensor< T >::sumsq(), madness::truncate(), and madness::NemoBase::world.
| 
 | friend | 
| 
 | friend | 
| 
 | protected | 
asymptotic correction for DFT
Referenced by compute_nemo_potentials(), get_ac(), and make_fock_operator().
| 
 | protected | 
Referenced by madness::OEP::OEP(), madness::OEP::OEP(), analyze(), compute_all_cphf(), madness::OEP::compute_and_print_final_energies(), compute_cphf_parallel_term(), madness::OEP::compute_E_first(), madness::OEP::compute_energy(), compute_energy_regularized(), compute_fock_matrix(), compute_nemo_potentials(), get_calc(), get_calc_param(), get_coulomb_potential(), gradient(), hessian(), is_dft(), madness::OEP::iterate(), kinetic_energy_potential(), load_mos(), madness::OEP::load_restartdata(), localize(), make_cphf_constant_term(), make_density(), make_density(), make_fock_operator(), make_incomplete_hessian_response_part(), molecule(), molecule(), madness::OEP::output_calc_info_schema(), purify_hessian(), madness::OEP::save_restartdata(), madness::OEP::selftest(), set_protocol(), solve(), madness::OEP::solve(), solve_cphf(), trantol(), value(), madness::OEP::value(), value(), and madness::OEP::value().
| 
 | mutableprivate | 
sum of square of coords at last solved geometry
Referenced by check_converged(), and value().
| const NemoCalculationParameters madness::Nemo::nemo_param | 
Referenced by get_nemo_param().
| 
 | private | 
polarizable continuum model
Referenced by Nemo(), Nemo(), compute_energy_regularized(), compute_nemo_potentials(), get_pcm(), and make_fock_operator().
| 
 | protected | 
a poisson solver
Referenced by kinetic_energy_potential(), make_laplacian_density(), and set_protocol().
| 
 | protected | 
Referenced by Nemo(), Nemo(), do_symmetry(), get_symmetry_projector(), madness::OEP::iterate(), and solve().