MADNESS
0.10.1
|
#include <molecule.h>
Classes | |
struct | apply_c2 |
Apply to (x,y,z) a C2 rotation about an axis thru the origin and (xaxis,yaxis,zaxis) More... | |
struct | apply_inverse |
struct | apply_sigma |
Apply to (x,y,z) a reflection through a plane containing the origin with normal (xaxis,yaxis,zaxis) More... | |
struct | GeometryParameters |
Public Member Functions | |
Molecule () | |
Makes a molecule with zero atoms. More... | |
Molecule (std::vector< Atom > atoms, double eprec, CorePotentialManager core_pot={}, madness::Tensor< double > field=madness::Tensor< double >(3L)) | |
makes a molecule from a list of atoms More... | |
Molecule (World &world, const commandlineparser &parser) | |
makes a molecule using contents of parser More... | |
void | add_atom (double x, double y, double z, double q, int atn) |
void | add_atom (double x, double y, double z, double q, int atn, bool psat) |
double | atomic_attraction_potential (int iatom, double x, double y, double z) const |
nuclear attraction potential for a specific atom in the molecule More... | |
double | bounding_cube () const |
Returns the half width of the bounding cube. More... | |
void | center () |
Moves the center of nuclear charge to the origin. More... | |
Tensor< double > | center_of_mass () const |
compute the center of mass More... | |
double | core_derivative (int atom, int axis, unsigned int core, int m, double x, double y, double z) const |
double | core_eval (int atom, unsigned int core, int m, double x, double y, double z) const |
double | core_potential_derivative (int atom, int axis, double x, double y, double z) const |
std::vector< std::string > | cubefile_header () const |
print out a Gaussian cubefile header More... | |
madness::Tensor< double > | get_all_coords () const |
std::vector< madness::Vector< double, 3 > > | get_all_coords_vec () const |
const Atom & | get_atom (unsigned int i) const |
unsigned int | get_atom_charge (unsigned int i) const |
unsigned int | get_atomic_number (unsigned int i) const |
const std::vector< Atom > & | get_atoms () const |
double | get_core_bc (unsigned int atn, unsigned int c) const |
unsigned int | get_core_l (unsigned int atn, unsigned int c) const |
double | get_eprec () const |
std::string | get_pointgroup () const |
bool | get_pseudo_atom (unsigned int i) const |
std::vector< double > | get_rcut () const |
void | get_structure () |
std::string | guess_file () const |
hashT | hash () const |
double | inter_atomic_distance (unsigned int i, unsigned int j) const |
bool | is_potential_defined (unsigned int atn) const |
bool | is_potential_defined_atom (int i) const |
Tensor< double > | massweights () const |
compute the mass-weighting matrix for the hessian More... | |
double | mol_nuclear_charge_density (double x, double y, double z) const |
double | molecular_core_potential (double x, double y, double z) const |
Tensor< double > | moment_of_inertia () const |
unsigned int | n_core_orb (unsigned int atn) const |
unsigned int | n_core_orb_all () const |
size_t | natom () const |
double | nuclear_attraction_potential (double x, double y, double z) const |
nuclear attraction potential for the whole molecule More... | |
double | nuclear_attraction_potential_derivative (int atom, int axis, double x, double y, double z) const |
double | nuclear_attraction_potential_second_derivative (int atom, int iaxis, int jaxis, double x, double y, double z) const |
the second derivative of the (smoothed) nuclear potential Z/r More... | |
double | nuclear_charge_density (double x, double y, double z) const |
double | nuclear_dipole (int axis) const |
compute the dipole moment of the nuclei More... | |
Tensor< double > | nuclear_dipole_derivative (const int atom, const int axis) const |
compute the derivative of the nuclear dipole wrt a nuclear displacement More... | |
double | nuclear_repulsion_derivative (size_t iatom, int axis) const |
double | nuclear_repulsion_energy () const |
Tensor< double > | nuclear_repulsion_hessian () const |
return the hessian matrix of the second derivatives d^2/dxdy V More... | |
double | nuclear_repulsion_second_derivative (int iatom, int jatom, int iaxis, int jaxis) const |
compute the nuclear-nuclear contribution to the second derivatives More... | |
void | orient (bool verbose=false) |
Centers and orients the molecule in a standard manner. More... | |
void | print () const |
void | read (std::istream &f) |
void | read_core_file (const std::string &filename) |
void | read_file (const std::string &filename) |
void | read_structure_from_library (const std::string &name) |
void | read_xyz (const std::string filename) |
void | rotate (const Tensor< double > &D) |
rotates the molecule and the external field More... | |
template<typename Archive > | |
void | serialize (Archive &ar) |
void | set_all_coords (const madness::Tensor< double > &newcoords) |
void | set_atom_charge (unsigned int i, double zeff) |
void | set_atom_coords (unsigned int i, double x, double y, double z) |
void | set_core_eprec (double value) |
void | set_core_rcut (double value) |
void | set_pseudo_atom (unsigned int i, bool psat) |
void | set_rcut (double value) |
double | smallest_length_scale () const |
std::string | symmetrize_and_identify_point_group (const double symtol) |
json | to_json () const |
double | total_nuclear_charge () const |
void | translate (const Tensor< double > &translation) |
translate the molecule More... | |
void | update_rcut_with_eprec (double value) |
updates rcuts with given eprec More... | |
Static Public Member Functions | |
static std::string | get_structure_library_path () |
static std::istream & | position_stream_in_library (std::ifstream &f, const std::string &name) |
static void | print_parameters () |
Public Attributes | |
std::vector< double > | atomic_radii |
GeometryParameters | parameters |
Private Member Functions | |
template<typename opT > | |
int | find_symmetry_equivalent_atom (int iatom, opT op, const double symtol) const |
void | swapaxes (int ix, int iy) |
template<typename opT > | |
void | symmetrize_for_op (opT op, const double symtol) |
bool | test_for_c2 (double xaxis, double yaxis, double zaxis, const double symtol) const |
bool | test_for_inverse (const double symtol) const |
template<typename opT > | |
bool | test_for_op (opT op, const double symtol) const |
bool | test_for_sigma (double xaxis, double yaxis, double zaxis, const double symtol) const |
Private Attributes | |
std::vector< Atom > | atoms |
CorePotentialManager | core_pot |
madness::Tensor< double > | field |
std::string | pointgroup_ ="c1" |
The molecular point group is automatically assigned in the identify_pointgroup function. More... | |
std::vector< double > | rcut |
|
inline |
Makes a molecule with zero atoms.
madness::Molecule::Molecule | ( | std::vector< Atom > | atoms, |
double | eprec, | ||
CorePotentialManager | core_pot = {} , |
||
madness::Tensor< double > | field = madness::Tensor<double>(3L) |
||
) |
makes a molecule from a list of atoms
References atomic_radii, madness::constants::atomic_unit_of_length, madness::AtomicData::covalent_radius, e(), madness::get_atomic_data(), and update_rcut_with_eprec().
madness::Molecule::Molecule | ( | World & | world, |
const commandlineparser & | parser | ||
) |
makes a molecule using contents of parser
References madness::WorldGopInterface::broadcast_serializable(), madness::Molecule::GeometryParameters::field(), field, get_structure(), madness::World::gop, MADNESS_CHECK, parameters, madness::QCCalculationParametersBase::print(), and madness::World::rank().
void madness::Molecule::add_atom | ( | double | x, |
double | y, | ||
double | z, | ||
double | q, | ||
int | atn | ||
) |
References atomic_radii, madness::constants::atomic_unit_of_length, atoms, c, madness::AtomicData::covalent_radius, e(), madness::get_atomic_data(), get_eprec(), q(), rcut, and madness::smoothing_parameter().
Referenced by main(), read(), read_xyz(), and madness::InitParameters::readnw().
void madness::Molecule::add_atom | ( | double | x, |
double | y, | ||
double | z, | ||
double | q, | ||
int | atn, | ||
bool | psat | ||
) |
double madness::Molecule::atomic_attraction_potential | ( | int | iatom, |
double | x, | ||
double | y, | ||
double | z | ||
) | const |
nuclear attraction potential for a specific atom in the molecule
References atoms, madness::distance(), rcut, madness::smoothed_potential(), and madness::sum().
Referenced by madness::NuclearCorrelationFactor::square_times_V_functor::operator()().
double madness::Molecule::bounding_cube | ( | ) | const |
void madness::Molecule::center | ( | ) |
Moves the center of nuclear charge to the origin.
References atoms, and translate().
Referenced by orient().
Tensor< double > madness::Molecule::center_of_mass | ( | ) | const |
compute the center of mass
References get_atom(), madness::Atom::mass, natom(), madness::Atom::x, madness::Atom::y, and madness::Atom::z.
Referenced by projector_external_dof(), and madness::MolecularOptimizer::projector_external_dof().
double madness::Molecule::core_derivative | ( | int | atom, |
int | axis, | ||
unsigned int | core, | ||
int | m, | ||
double | x, | ||
double | y, | ||
double | z | ||
) | const |
References atoms, axis, madness::CorePotentialManager::core_derivative(), core_pot, m, and xi.
Referenced by madness::CoreOrbitalDerivativeFunctor::operator()().
double madness::Molecule::core_eval | ( | int | atom, |
unsigned int | core, | ||
int | m, | ||
double | x, | ||
double | y, | ||
double | z | ||
) | const |
References atoms, madness::CorePotentialManager::core_eval(), core_pot, and m.
Referenced by madness::CoreOrbitalFunctor::operator()().
double madness::Molecule::core_potential_derivative | ( | int | atom, |
int | axis, | ||
double | x, | ||
double | y, | ||
double | z | ||
) | const |
References atoms, axis, core_pot, madness::distance(), natom(), madness::CorePotentialManager::potential_derivative(), and xi.
Referenced by madchem::CorePotentialDerivativeFunctor::operator()().
std::vector< std::string > madness::Molecule::cubefile_header | ( | ) | const |
print out a Gaussian cubefile header
References get_atom(), madness::Atom::get_atomic_number(), madness::Atom::get_coords(), and natom().
Referenced by madness::cubefile_header(), and main().
|
private |
References atoms, madness::distance(), and op().
Referenced by symmetrize_for_op(), and test_for_op().
madness::Tensor< double > madness::Molecule::get_all_coords | ( | ) | const |
References c, get_atom(), natom(), madness::Atom::x, madness::Atom::y, and madness::Atom::z.
Referenced by madness::Znemo::gradient(), main(), madness::MolecularEnergy::output_calc_info_schema(), run_all_calculations(), test_read_restartaodata(), test_read_restartdata(), and madness::Znemo::value().
std::vector< madness::Vector< double, 3 > > madness::Molecule::get_all_coords_vec | ( | ) | const |
References c, get_atom(), natom(), madness::Atom::x, madness::Atom::y, and madness::Atom::z.
Referenced by GaussianNucleusFunctor::GaussianNucleusFunctor(), madness::Znemo::recompute_factors_and_potentials(), MolecularPotentialFunctor::special_points(), madchem::MolecularGuessDensityFunctor::special_points(), madness::MolecularPotentialFunctor::special_points(), madness::MolecularCorePotentialFunctor::special_points(), and madness::write_molecules_to_file().
const Atom & madness::Molecule::get_atom | ( | unsigned int | i | ) | const |
References atoms.
Referenced by madness::GTHPseudopotential< Q >::apply_potential(), madness::GTHPseudopotential< Q >::apply_potential_simple(), center_of_mass(), madness::PotentialManager::core_projection(), madness::PotentialManager::core_projector_derivative(), cubefile_header(), madness::SCF::derivatives(), get_all_coords(), get_all_coords_vec(), madness::Znemo::gradient(), madness::BasisFunctions::guess_virtuals_internal(), madness::SCF::initial_guess(), madness::GTHPseudopotential< Q >::load_pseudo_from_file(), madness::GTHPseudopotential< Q >::make_pseudo_potential(), massweights(), madchem::AtomicAttractionFunctor::operator()(), projector_external_dof(), madness::MolecularOptimizer::projector_external_dof(), madchem::AtomicAttractionFunctor::special_points(), madchem::MolecularDerivativeFunctor::special_points(), madchem::MolecularSecondDerivativeFunctor::special_points(), and test_ethylene().
unsigned int madness::Molecule::get_atom_charge | ( | unsigned int | i | ) | const |
References atoms.
Referenced by madness::SCF::initial_guess().
unsigned int madness::Molecule::get_atomic_number | ( | unsigned int | i | ) | const |
|
inline |
|
inline |
References c, core_pot, and madness::CorePotentialManager::get_core_bc().
Referenced by madness::PotentialManager::core_projection(), and madness::PotentialManager::core_projector_derivative().
|
inline |
References c, core_pot, and madness::CorePotentialManager::get_core_l().
Referenced by madness::PotentialManager::core_projection(), and madness::PotentialManager::core_projector_derivative().
|
inline |
References madness::Molecule::GeometryParameters::eprec(), and parameters.
Referenced by madness::GaussSlater::GaussSlater(), madness::GradientalGaussSlater::GradientalGaussSlater(), madness::LinearSlater::LinearSlater(), madness::poly4erfc::poly4erfc(), madness::Polynomial< N >::Polynomial(), madness::PseudoNuclearCorrelationFactor::PseudoNuclearCorrelationFactor(), madness::Slater::Slater(), add_atom(), read_core_file(), madness::NuclearCorrelationFactor::smoothed_unitvec(), and update_rcut_with_eprec().
|
inline |
References pointgroup_.
bool madness::Molecule::get_pseudo_atom | ( | unsigned int | i | ) | const |
References atoms.
Referenced by get_structure(), and madness::SCF::initial_guess().
|
inline |
References rcut.
Referenced by madchem::AtomicAttractionFunctor::operator()().
void madness::Molecule::get_structure | ( | ) |
References madness::Molecule::GeometryParameters::core_type(), madness::QCCalculationParametersBase::get(), get_atomic_number(), madness::get_charge_from_file(), get_pseudo_atom(), MADNESS_EXCEPTION, natom(), madness::Molecule::GeometryParameters::no_orient(), orient(), parameters, read(), read_core_file(), read_structure_from_library(), read_xyz(), set_atom_charge(), set_pseudo_atom(), madness::Molecule::GeometryParameters::source_name(), and madness::Molecule::GeometryParameters::source_type().
Referenced by Molecule().
|
static |
|
inline |
References core_pot, and madness::CorePotentialManager::guess_file().
|
inline |
References atoms, h(), madness::hash_combine(), madness::hash_range(), pointgroup_, and rcut.
double madness::Molecule::inter_atomic_distance | ( | unsigned int | i, |
unsigned int | j | ||
) | const |
References atoms, and madness::distance().
Referenced by nuclear_repulsion_derivative(), nuclear_repulsion_energy(), and nuclear_repulsion_second_derivative().
|
inline |
References core_pot, and madness::CorePotentialManager::is_defined().
|
inline |
References atoms, core_pot, and madness::CorePotentialManager::is_defined().
Referenced by madness::SCF::derivatives().
|
inline |
compute the mass-weighting matrix for the hessian
use as mass_weighted_hessian=inner(massweights,inner(hessian,massweights));
References get_atom(), and natom().
Referenced by compute_frequencies(), madness::Nemo::compute_IR_intensities(), and compute_reduced_mass().
double madness::Molecule::mol_nuclear_charge_density | ( | double | x, |
double | y, | ||
double | z | ||
) | const |
References atoms, madness::distance(), rcut, and madness::smoothed_density().
Referenced by NuclearDensityFunctor::operator()().
double madness::Molecule::molecular_core_potential | ( | double | x, |
double | y, | ||
double | z | ||
) | const |
Tensor< double > madness::Molecule::moment_of_inertia | ( | ) | const |
References atoms, I, k, and q().
Referenced by projector_external_dof(), and madness::MolecularOptimizer::projector_external_dof().
|
inline |
unsigned int madness::Molecule::n_core_orb_all | ( | ) | const |
References atoms, core_pot, madness::CorePotentialManager::is_defined(), madness::CorePotentialManager::n_core_orb(), natom(), and madness::sum().
Referenced by madness::SCF::initial_guess(), and madness::SCF::load_mos().
|
inline |
References atoms.
Referenced by NuclearVector::NuclearVector(), center_of_mass(), madness::Nemo::compute_all_cphf(), compute_frequencies(), madness::Nemo::compute_IR_intensities(), compute_reduced_mass(), core_potential_derivative(), madness::PotentialManager::core_projection(), cubefile_header(), madness::SCF::derivatives(), get_all_coords(), get_all_coords_vec(), get_structure(), madness::Znemo::gradient(), madness::BasisFunctions::guess_virtuals_internal(), madness::Nemo::hessian(), madness::SCF::initial_guess(), madness::GTHPseudopotential< Q >::load_pseudo_from_file(), madness::Localizer::localize_new(), madness::Nemo::make_incomplete_hessian(), madness::Nemo::make_incomplete_hessian_response_part(), madness::GTHPseudopotential< Q >::make_pseudo_potential(), massweights(), molecular_core_potential(), n_core_orb_all(), nuclear_generator(), nuclear_repulsion_hessian(), madness::MolecularEnergy::output_calc_info_schema(), print(), projector_external_dof(), madness::MolecularOptimizer::projector_external_dof(), read_core_file(), read_xyz(), set_all_coords(), madness::ResponseParameters::set_derived_values(), madness::SCF::solve(), madness::Nemo::solve_cphf(), to_json(), madness::MolecularEnergy::value(), madness::Znemo::value(), and madness::write_molecules_to_file().
double madness::Molecule::nuclear_attraction_potential | ( | double | x, |
double | y, | ||
double | z | ||
) | const |
nuclear attraction potential for the whole molecule
References atoms, madness::distance(), field, rcut, madness::smoothed_potential(), and madness::sum().
Referenced by MolecularPotentialFunctor::operator()(), and madness::MolecularPotentialFunctor::operator()().
double madness::Molecule::nuclear_attraction_potential_derivative | ( | int | atom, |
int | axis, | ||
double | x, | ||
double | y, | ||
double | z | ||
) | const |
double madness::Molecule::nuclear_attraction_potential_second_derivative | ( | int | atom, |
int | iaxis, | ||
int | jaxis, | ||
double | x, | ||
double | y, | ||
double | z | ||
) | const |
the second derivative of the (smoothed) nuclear potential Z/r
with
References atoms, madness::d2smoothed_potential(), madness::Vector< T, N >::normf(), rcut, madness::smoothed_potential(), and u().
Referenced by madchem::MolecularSecondDerivativeFunctor::operator()().
double madness::Molecule::nuclear_charge_density | ( | double | x, |
double | y, | ||
double | z | ||
) | const |
References atoms, madness::distance_sq(), rcut, and madness::smoothed_density().
Referenced by NuclearDensityFunctor::operator()().
double madness::Molecule::nuclear_dipole | ( | int | axis | ) | const |
compute the dipole moment of the nuclei
[in] | axis | the axis (x, y, z) |
References atoms, axis, core_pot, madness::CorePotentialManager::is_defined(), MADNESS_EXCEPTION, madness::CorePotentialManager::n_core_orb(), and madness::sum().
Referenced by madness::SCF::dipole().
Tensor< double > madness::Molecule::nuclear_dipole_derivative | ( | const int | atom, |
const int | axis | ||
) | const |
compute the derivative of the nuclear dipole wrt a nuclear displacement
[in] | atom | the atom which will be displaced |
[in] | axis | the axis where the atom will be displaced |
Referenced by madness::Nemo::compute_IR_intensities().
double madness::Molecule::nuclear_repulsion_derivative | ( | size_t | iatom, |
int | axis | ||
) | const |
References atoms, axis, core_pot, inter_atomic_distance(), madness::CorePotentialManager::is_defined(), madness::CorePotentialManager::n_core_orb(), and madness::sum().
Referenced by madness::SCF::derivatives().
double madness::Molecule::nuclear_repulsion_energy | ( | ) | const |
References atoms, core_pot, inter_atomic_distance(), madness::CorePotentialManager::is_defined(), madness::CorePotentialManager::n_core_orb(), and madness::sum().
Referenced by madness::Znemo::compute_energy(), main(), and madness::SCF::solve().
Tensor< double > madness::Molecule::nuclear_repulsion_hessian | ( | ) | const |
return the hessian matrix of the second derivatives d^2/dxdy V
compute the nuclear-nuclear contribution to the molecular hessian
no factor 0.5 included
References natom(), and nuclear_repulsion_second_derivative().
Referenced by madness::Nemo::hessian(), and madness::Nemo::make_incomplete_hessian().
double madness::Molecule::nuclear_repulsion_second_derivative | ( | int | iatom, |
int | jatom, | ||
int | iaxis, | ||
int | jaxis | ||
) | const |
compute the nuclear-nuclear contribution to the second derivatives
[in] | iatom | the i-th atom (row of the hessian) |
[in] | jatom | the j-th atom (column of the hessian) |
[in] | iaxis | the xyz axis of the i-th atom |
[in] | jaxis | the xyz axis of the j-th atom return the (3*iatom + iaxis, 3*jatom + jaxis) matix element of the hessian |
[in] | iatom | the i-th atom (row of the hessian) |
[in] | jatom | the j-th atom (column of the hessian) |
[in] | iaxis | the xyz axis of the i-th atom |
[in] | jaxis | the xyz axis of the j-th atom return the (3*iatom + iaxis, 3*jatom + jaxis) matix element |
References atoms, core_pot, inter_atomic_distance(), madness::CorePotentialManager::is_defined(), MADNESS_EXCEPTION, pow(), and madness::sum().
Referenced by nuclear_repulsion_hessian().
void madness::Molecule::orient | ( | bool | verbose = false | ) |
Centers and orients the molecule in a standard manner.
References atoms, center(), e(), I, k, MADNESS_CHECK, parameters, pointgroup_, madness::print(), q(), rotate(), madness::syev(), symmetrize_and_identify_point_group(), and madness::Molecule::GeometryParameters::symtol().
Referenced by madness::cubefile_header(), get_structure(), and main().
|
static |
void madness::Molecule::print | ( | ) | const |
|
static |
References param, and madness::print().
Referenced by madness::CC2::print_parameters(), madness::MP2::print_parameters(), madness::Nemo::print_parameters(), madness::OEP::print_parameters(), madness::PNO::print_parameters(), madness::SCF::print_parameters(), madness::TDHF::print_parameters(), madness::Zcis::print_parameters(), and madness::Znemo::print_parameters().
void madness::Molecule::read | ( | std::istream & | f | ) |
References add_atom(), madness::constants::atomic_unit_of_length, atoms, madness::check_if_pseudo_atom(), e(), madness::Molecule::GeometryParameters::eprec(), madness::f, madness::QCCalculationParametersBase::get_all_parameters(), p(), parameters, madness::position_stream(), rcut, madness::scale(), madness::symbol_to_atomic_number(), madness::Molecule::GeometryParameters::units(), and update_rcut_with_eprec().
Referenced by get_structure(), read_file(), and read_structure_from_library().
void madness::Molecule::read_core_file | ( | const std::string & | filename | ) |
References atoms, core_pot, madness::filename, get_eprec(), madness::CorePotentialManager::is_defined(), madness::CorePotentialManager::n_core_orb(), natom(), madness::print(), q(), rcut, madness::CorePotentialManager::read_file(), and madness::smoothing_parameter().
Referenced by madness::SCF::SCF(), and get_structure().
void madness::Molecule::read_file | ( | const std::string & | filename | ) |
References errmsg(), madness::f, madness::filename, MADNESS_EXCEPTION, and read().
Referenced by madness::cubefile_header(), and main().
void madness::Molecule::read_structure_from_library | ( | const std::string & | name | ) |
References errmsg(), madness::f, MADNESS_EXCEPTION, madness::name(), position_stream_in_library(), and read().
Referenced by get_structure().
void madness::Molecule::read_xyz | ( | const std::string | filename | ) |
References add_atom(), madness::constants::atomic_unit_of_length, atoms, madness::check_if_pseudo_atom(), e(), madness::Molecule::GeometryParameters::eprec(), errmsg(), madness::f, madness::filename, MADNESS_CHECK, MADNESS_EXCEPTION, natom(), parameters, rcut, madness::scale(), madness::symbol_to_atomic_number(), madness::Molecule::GeometryParameters::units(), and update_rcut_with_eprec().
Referenced by get_structure().
void madness::Molecule::rotate | ( | const Tensor< double > & | D | ) |
rotates the molecule and the external field
[in] | D | the rotation matrix |
References atoms, field, and madness::inner().
Referenced by madness::SCF::initial_guess(), and orient().
|
inline |
References atoms, core_pot, field, parameters, pointgroup_, and rcut.
void madness::Molecule::set_all_coords | ( | const madness::Tensor< double > & | newcoords | ) |
References atoms, c, MADNESS_ASSERT, and natom().
Referenced by madness::MolecularEnergy::value(), and madness::Znemo::value().
void madness::Molecule::set_atom_charge | ( | unsigned int | i, |
double | zeff | ||
) |
References atoms.
Referenced by get_structure().
void madness::Molecule::set_atom_coords | ( | unsigned int | i, |
double | x, | ||
double | y, | ||
double | z | ||
) |
References atoms.
|
inline |
References core_pot, and madness::CorePotentialManager::set_eprec().
|
inline |
References core_pot, and madness::CorePotentialManager::set_rcut().
void madness::Molecule::set_pseudo_atom | ( | unsigned int | i, |
bool | psat | ||
) |
References atoms.
Referenced by get_structure().
double madness::Molecule::smallest_length_scale | ( | ) | const |
References atoms, max, and rcut.
Referenced by NuclearDensityFunctor::operator()(), and madness::ResponseParameters::set_ground_state_calculation_data().
|
private |
References atoms, field, and madness::swap().
Referenced by symmetrize_and_identify_point_group().
std::string madness::Molecule::symmetrize_and_identify_point_group | ( | const double | symtol | ) |
References madness::inverse(), madness::print(), swapaxes(), symmetrize_for_op(), test_for_c2(), test_for_inverse(), and test_for_sigma().
Referenced by orient().
|
private |
References atoms, madness::distance(), find_symmetry_equivalent_atom(), MADNESS_CHECK, op(), and test_for_op().
Referenced by symmetrize_and_identify_point_group().
|
private |
References test_for_op().
Referenced by symmetrize_and_identify_point_group().
|
private |
References test_for_op().
Referenced by symmetrize_and_identify_point_group().
|
private |
References atoms, find_symmetry_equivalent_atom(), and op().
Referenced by symmetrize_for_op(), test_for_c2(), test_for_inverse(), and test_for_sigma().
|
private |
References test_for_op().
Referenced by symmetrize_and_identify_point_group().
nlohmann::json madness::Molecule::to_json | ( | ) | const |
References atoms, madness::get_atomic_data(), natom(), and madness::AtomicData::symbol.
Referenced by madness::MolecularEnergy::output_calc_info_schema().
double madness::Molecule::total_nuclear_charge | ( | ) | const |
References atoms, and madness::sum().
Referenced by MiniDFT::doit(), and main().
void madness::Molecule::translate | ( | const Tensor< double > & | translation | ) |
translate the molecule
References atoms.
Referenced by center(), madness::SCF::initial_guess(), projector_external_dof(), and madness::MolecularOptimizer::projector_external_dof().
void madness::Molecule::update_rcut_with_eprec | ( | double | value | ) |
updates rcuts with given eprec
References atoms, core_pot, get_eprec(), parameters, rcut, madness::CorePotentialManager::set_eprec(), madness::QCCalculationParametersBase::set_user_defined_value(), and madness::smoothing_parameter().
Referenced by Molecule(), dnuclear_anchor_test(), main(), read(), and read_xyz().
std::vector<double> madness::Molecule::atomic_radii |
Referenced by Molecule(), and add_atom().
|
private |
Referenced by add_atom(), atomic_attraction_potential(), bounding_cube(), center(), core_derivative(), core_eval(), core_potential_derivative(), find_symmetry_equivalent_atom(), get_atom(), get_atom_charge(), get_atomic_number(), get_atoms(), get_pseudo_atom(), hash(), inter_atomic_distance(), is_potential_defined_atom(), mol_nuclear_charge_density(), molecular_core_potential(), moment_of_inertia(), n_core_orb_all(), natom(), nuclear_attraction_potential(), nuclear_attraction_potential_derivative(), nuclear_attraction_potential_second_derivative(), nuclear_charge_density(), nuclear_dipole(), nuclear_dipole_derivative(), nuclear_repulsion_derivative(), nuclear_repulsion_energy(), nuclear_repulsion_second_derivative(), orient(), print(), read(), read_core_file(), read_xyz(), rotate(), serialize(), set_all_coords(), set_atom_charge(), set_atom_coords(), set_pseudo_atom(), set_rcut(), smallest_length_scale(), swapaxes(), symmetrize_for_op(), test_for_op(), to_json(), total_nuclear_charge(), translate(), and update_rcut_with_eprec().
|
private |
Referenced by core_derivative(), core_eval(), core_potential_derivative(), get_core_bc(), get_core_l(), guess_file(), is_potential_defined(), is_potential_defined_atom(), molecular_core_potential(), n_core_orb(), n_core_orb_all(), nuclear_dipole(), nuclear_repulsion_derivative(), nuclear_repulsion_energy(), nuclear_repulsion_second_derivative(), read_core_file(), serialize(), set_core_eprec(), set_core_rcut(), and update_rcut_with_eprec().
|
private |
Referenced by Molecule(), nuclear_attraction_potential(), nuclear_attraction_potential_derivative(), rotate(), serialize(), and swapaxes().
GeometryParameters madness::Molecule::parameters |
Referenced by Molecule(), madness::Nuclear< T, NDIM >::Nuclear(), madness::SCF::SCF(), madness::SCF::apply_potential(), madness::SCF::derivatives(), get_eprec(), get_structure(), madness::SCF::initial_guess(), madness::SCF::initial_load_bal(), madness::SCF::loadbal(), madness::SCF::make_nuclear_potential(), orient(), madness::MolecularEnergy::output_calc_info_schema(), print(), read(), read_xyz(), madness::Znemo::recompute_factors_and_potentials(), serialize(), madness::SCF::solve(), and update_rcut_with_eprec().
|
private |
The molecular point group is automatically assigned in the identify_pointgroup function.
Referenced by get_pointgroup(), hash(), orient(), and serialize().
|
private |
Referenced by add_atom(), atomic_attraction_potential(), get_rcut(), hash(), mol_nuclear_charge_density(), nuclear_attraction_potential(), nuclear_attraction_potential_derivative(), nuclear_attraction_potential_second_derivative(), nuclear_charge_density(), read(), read_core_file(), read_xyz(), serialize(), set_rcut(), smallest_length_scale(), and update_rcut_with_eprec().