21#ifndef VOTCA_XTP_PMLOCALIZATION_H
22#define VOTCA_XTP_PMLOCALIZATION_H
35 method_ = options.
get(
".method").as<std::string>();
47 double cost(
const Eigen::MatrixXd &W,
48 const std::vector<Eigen::MatrixXd> &Sat_all,
49 const Index nat)
const;
51 const Eigen::MatrixXd &W,
const std::vector<Eigen::MatrixXd> &Sat_all,
52 const Index nat)
const;
55 const Eigen::VectorXd &y)
const;
59 Eigen::MatrixXd
rotate_W(
const double step,
const Eigen::MatrixXd &W,
60 const Eigen::VectorXcd &eval,
61 const Eigen::MatrixXcd &evec)
const;
64 const Eigen::MatrixXd &occ_orbitals);
66 double inner_prod(
const Eigen::MatrixXd &A,
const Eigen::MatrixXd &B)
const {
67 return (0.5 * A.transpose() * B).trace();
78 std::pair<Eigen::MatrixXd, Eigen::VectorXd>
sort_lmos(
79 const Eigen::VectorXd &energies);
81 Eigen::VectorXd
pop_per_atom(
const Eigen::VectorXd &orbital);
Container to hold Basisfunctions for all atoms.
Logger is used for thread-safe output of messages.
container for molecular orbitals
Eigen::VectorXcd find_complex_roots(const Eigen::VectorXcd &coeff) const
double inner_prod(const Eigen::MatrixXd &A, const Eigen::MatrixXd &B) const
Eigen::MatrixXd MullikenPop_orb_per_atom_
void check_orthonormality()
void computePML(Orbitals &orbitals)
Eigen::VectorXd calculate_lmo_energies(const Orbitals &orbitals)
void computePML_UT(Orbitals &orbitals)
Eigen::Vector2d offdiag_penalty_elements(const Eigen::MatrixXd &s_overlap, Index s, Index t)
std::pair< Eigen::MatrixXd, Eigen::VectorXd > sort_lmos(const Eigen::VectorXd &energies)
Eigen::MatrixXd PM_penalty_
void update_penalty(Index s, Index t)
double find_smallest_step(const Eigen::VectorXd &coeff) const
void computePML_JS(Orbitals &orbitals)
Eigen::MatrixXcd companion_matrix(const Eigen::VectorXcd &coeff) const
Eigen::VectorXd fit_polynomial(const Eigen::VectorXd &x, const Eigen::VectorXd &y) const
PMLocalization(Logger &log, const tools::Property &options)
Eigen::MatrixX2d rotateorbitals(const Eigen::MatrixX2d &maxorbs, Index s, Index t)
Eigen::MatrixXd localized_orbitals_
std::vector< Eigen::MatrixXd > setup_pop_matrices(const Eigen::MatrixXd &occ_orbitals)
Eigen::VectorXd pop_per_atom(const Eigen::VectorXd &orbital)
double convergence_limit_
std::vector< Index > numfuncpatom_
std::pair< double, Eigen::MatrixXd > cost_derivative(const Eigen::MatrixXd &W, const std::vector< Eigen::MatrixXd > &Sat_all, const Index nat) const
Eigen::MatrixXd rotate_W(const double step, const Eigen::MatrixXd &W, const Eigen::VectorXcd &eval, const Eigen::MatrixXcd &evec) const
double cost(const Eigen::MatrixXd &W, const std::vector< Eigen::MatrixXd > &Sat_all, const Index nat) const
base class for all analysis tools