Molecule

class psi4.core.Molecule

Bases: pybind11_object

Class to store the elements, coordinates, fragmentation pattern, basis sets, charge, multiplicity, etc. of a molecule.

Methods Summary

B787(ref_mol[, do_plot, verbose, atoms_map, ...])

Finds shift, rotation, and atom reordering of concern_mol that best aligns with ref_mol.

BFS([seed_atoms, bond_threshold, ...])

Detect fragments among real atoms through a breadth-first search (BFS) algorithm.

Z(self, arg0)

Nuclear charge of atom arg0 (0-indexed without dummies)

activate_all_fragments(self)

Sets all fragments in the molecule to be active

add_atom(self, Z, x, y, z, symbol, mass, ...)

Adds to self Molecule an atom with atomic number Z, Cartesian coordinates in Bohr (x, y, z), atomic symbol, mass, and charge, extended atomic label, and mass number A

atom_at_position(*args, **kwargs)

Overloaded function.

basis_on_atom(self, arg0)

Gets the label of the orbital basis set on a given atom arg0

center_of_mass(self)

Computes center of mass of molecule (does not translate molecule)

charge(self, atom)

Gets charge of atom (0-indexed without dummies)

clone(self)

Returns a new Molecule identical to arg1

com_fixed(self)

Gets whether or not center of mass is fixed

comment(self)

Gets molecule comment

connectivity(self)

Gets molecule connectivity

create_psi4_string_from_molecule(self)

Gets a string re-expressing in input format the current state of the molecule.Contains Cartesian geometry info, fragmentation, charges and multiplicities, and any frame restriction.

deactivate_all_fragments(self)

Sets all fragments in the molecule to be inactive

distance_matrix(self)

Returns Matrix of interatom distances

extract_subsets(*args, **kwargs)

Overloaded function.

fZ(self, arg0)

Nuclear charge of atom arg1 (0-indexed including dummies)

fcharge(self, atom)

Gets charge of atom (0-indexed including dummies)

find_highest_point_group(self[, tolerance])

Finds highest possible computational molecular point group

find_point_group(self[, tolerance])

Finds computational molecular point group, user can override this with the symmetry keyword

fix_com(self, arg0)

Sets whether to fix the Cartesian position, or to translate to the C.O.M.

fix_orientation(self, arg0)

Fix the orientation at its current frame.

flabel(self, atom)

Gets the original label of the atom arg0 as given in the input file (C2, H4)(0-indexed including dummies)

fmass(self, atom)

Gets mass of atom (0-indexed including dummies)

form_symmetry_information(self, arg0)

Uses the point group object obtain by calling point_group()

format_molecule_for_mol()

Returns a string of Molecule formatted for mol2.

from_arrays([geom, elea, elez, elem, mass, ...])

Construct Molecule from unvalidated arrays and variables.

from_dict(arg0)

Returns a new Molecule constructed from python dictionary.

from_schema(molschema[, return_dict, ...])

Construct Molecule from non-Psi4 schema.

from_string(molstr[, dtype, name, fix_com, ...])

fsymbol(self, atom)

Gets the cleaned up label of atom (C2 => C, H4 = H) (0-indexed including dummies)

ftrue_atomic_number(self, atom)

Gets atomic number of atom from element (0-indexed including dummies)

full_geometry(self)

Gets the geometry [Bohr] as a (Natom X 3) matrix of coordinates (including dummies)

full_pg_n(self)

Gets n in Cnv, etc.; If there is no n (e.g. Td) it's the highest-order rotation axis.

fx(self, arg0)

x position of atom arg0 (0-indexed including dummies in Bohr)

fy(self, arg0)

y position of atom arg0 (0-indexed including dummies in Bohr)

fz(self, arg0)

z position of atom arg0 (0-indexed including dummies in Bohr)

geometry(self)

Gets the geometry [Bohr] as a (Natom X 3) matrix of coordinates (excluding dummies)

get_fragment_charges(self)

Gets the charge of each fragment

get_fragment_multiplicities(self)

Gets the multiplicity of each fragment

get_fragment_types(self)

Returns a list describing how to handle each fragment {Real, Ghost, Absent}

get_fragments(self)

Returns list of pairs of atom ranges defining each fragment from parent molecule(fragments[frag_ind] = <Afirst,Alast+1>)

get_full_point_group(self)

Gets point group name such as C3v or S8

get_full_point_group_with_n(self)

Gets point group name such as Cnv or Sn

get_variable(self, arg0)

Returns the value of variable arg0 in the structural variables list

has_zmatrix(self)

Get whether or not this molecule has at least one zmatrix entry

inertia_tensor(self)

Returns intertial tensor

input_units_to_au(self)

Returns unit conversion to [a0] for geometry

irrep_labels(self)

Returns Irreducible Representation symmetry labels

is_variable(self, arg0)

Checks if variable arg0 is in the structural variables list

label(self, atom)

Gets the original label of the atom as given in the input file (C2, H4)(0-indexed without dummies)

mass(self, atom)

Returns mass of atom (0-indexed)

mass_number(self, arg0)

Mass number (A) of atom if known, else -1

molecular_charge(self)

Gets the molecular charge

move_to_com(self)

Moves molecule to center of mass

multiplicity(self)

Gets the multiplicity (defined as 2Ms + 1)

nallatom(self)

Number of real and dummy atoms

name(self)

Gets molecule name

natom(self)

Number of real atoms

nfragments(self)

Gets the number of fragments in the molecule

nuclear_dipole(*args, **kwargs)

Overloaded function.

nuclear_repulsion_energy(self[, dipole_field])

Computes nuclear repulsion energy

nuclear_repulsion_energy_deriv1(self[, ...])

Returns first derivative of nuclear repulsion energy as a matrix (natom, 3)

nuclear_repulsion_energy_deriv2(self)

Returns second derivative of nuclear repulsion energy as a matrix (natom X 3, natom X 3)

orientation_fixed(self)

Get whether or not orientation is fixed

point_group(self)

Returns the current point group object

print_bond_angles(self)

Print the bond angle geometrical parameters

print_cluster(self)

Prints the molecule in Cartesians in input units adding fragment separators

print_distances(self)

Print the interatomic distance geometrical parameters

print_in_input_format(self)

Prints the molecule as Cartesian or ZMatrix entries, just as inputted.

print_out(self)

Prints the molecule in Cartesians in input units to output file

print_out_in_angstrom(self)

Prints the molecule in Cartesians in Angstroms to output file

print_out_in_bohr(self)

Prints the molecule in Cartesians in Bohr to output file

print_out_of_planes(self)

Print the out-of-plane angle geometrical parameters to output file

print_rotational_constants(self)

Print the rotational constants to output file

provenance(self)

Gets molecule provenance

reinterpret_coordentry(self, arg0)

Do reinterpret coordinate entries during update_geometry().

reset_point_group(self, arg0)

Overrides symmetry from outside the molecule string

rotational_constants(self)

Returns the rotational constants [cm^-1] of the molecule

rotational_symmetry_number(self)

Returns number of unique orientations of the rigid molecule that only interchange identical atoms

rotor_type(self)

Returns rotor type, e.g. 'RT_ATOM' or 'RT_SYMMETRIC_TOP'.

run_dftd3([func, dashlvl, dashparam, ...])

Compute dispersion correction via Grimme's DFTD3 program.

run_dftd4([func, dashlvl, dashparam, ...])

Compute dispersion correction via Grimme's DFTD4 program.

run_gcp([func, dertype, verbose])

Compute geometrical BSSE correction via Grimme's GCP program.

save_string_xyz(self)

Saves the string of an XYZ file to arg2

save_string_xyz_file(self)

Saves an XYZ file to arg2

save_xyz_file(self, arg0, arg1)

Saves an XYZ file to arg0

schoenflies_symbol(self)

Returns the Schoenflies symbol

scramble([do_shift, do_rotate, do_resort, ...])

Tester for B787 by shifting, rotating, and atom shuffling ref_mol and checking that the aligner returns the opposite transformation.

set_active_fragment(self, arg0)

Sets the specified fragment arg0 to be Real

set_active_fragments(self, arg0)

Sets the specified list arg0 of fragments to be Real

set_basis_all_atoms(self, arg0, arg1)

Sets basis set arg0 to all atoms

set_basis_by_label(self, arg0, arg1, arg2)

Sets basis set arg1 to all atoms with label (e.g., H4) arg0

set_basis_by_number(self, arg0, arg1, arg2)

Sets basis set arg1 to all atoms with number arg0

set_basis_by_symbol(self, arg0, arg1, arg2)

Sets basis set arg1 to all atoms with symbol (e.g., H) arg0

set_comment(self, arg0)

Sets molecule comment

set_connectivity(self, arg0)

Sets molecule connectivity

set_full_geometry(self, arg0)

Sets the geometry, given a (Natom X 3) matrix arg0 of coordinates (in Bohr) (including dummies

set_geometry(self, arg0)

Sets the geometry, given a (Natom X 3) matrix arg0 of coordinates [a0] (excluding dummies)

set_ghost_fragment(self, arg0)

Sets the specified fragment arg0 to be Ghost

set_ghost_fragments(self, arg0)

Sets the specified list arg0 of fragments to be Ghost

set_input_units_to_au(self, arg0)

Sets unit conversion to [a0] for geometry

set_mass(self, atom, mass)

Sets mass of atom (0-indexed) to mass (good for isotopic substitutions)

set_molecular_charge(self, arg0)

Change the overall molecular charge.

set_multiplicity(self, arg0)

Change the multiplicity (defined as 2S + 1).

set_name(self, arg0)

Sets molecule name

set_nuclear_charge(self, arg0, arg1)

Set the nuclear charge of the given atom arg0 to the value arg1 (primarily for ECP).

set_point_group(self, arg0)

Sets the molecular point group to the point group object arg0

set_provenance(self, arg0)

Sets molecule provenance

set_units(self, arg0)

Sets units (Angstrom or Bohr) used to define the geometry.

set_variable(self, arg0, arg1)

Sets the value arg1 to the variable arg0 in the list of structure variables, then calls update_geometry()

symbol(self, atom)

Gets the cleaned up label of atom (C2 => C, H4 = H) (0-indexed without dummies)

symmetrize(self, arg0)

Finds the highest point Abelian point group within the specified tolerance, and forces the geometry to have that symmetry.

symmetry_from_input(self)

Returns the symmetry specified in the input

to_arrays([dummy, ghost_as_dummy])

Exports coordinate info into NumPy arrays.

to_dict([force_c1, force_units, np_out])

Serializes instance into Molecule dictionary.

to_schema(dtype[, units])

Serializes instance into dictionary according to schema dtype.

to_string(dtype[, units, atom_format, ...])

Format a string representation of QM molecule.

translate(self, arg0)

Translates molecule by arg0

true_atomic_number(self, atom)

Gets atomic number of atom from element (0-indexed without dummies)

units(self)

Returns units used to define the geometry, i.e. 'Angstrom' or 'Bohr'.

update_geometry(self)

Reevaluates the geometry with current variable values, orientation directives, etc.

x(self, arg0)

x position [Bohr] of atom arg0 (0-indexed without dummies)

xyz(self, i)

Return the Vector3 for atom i (0-indexed without dummies)

y(self, arg0)

y position [Bohr] of atom arg0 (0-indexed without dummies)

z(self, arg0)

z position [Bohr] of atom arg0 (0-indexed without dummies)

Methods Documentation

B787(ref_mol, do_plot=False, verbose=1, atoms_map=False, run_resorting=False, mols_align=False, run_to_completion=False, uno_cutoff=0.001, run_mirror=False)[source]

Finds shift, rotation, and atom reordering of concern_mol that best aligns with ref_mol.

Wraps qcelemental.molutil.B787() for psi4.driver.qcdb.Molecule or psi4.core.Molecule. Employs the Kabsch, Hungarian, and Uno algorithms to exhaustively locate the best alignment for non-oriented, non-ordered structures.

Parameters:
  • concern_mol (Union[Molecule, Molecule]) – Molecule of concern, to be shifted, rotated, and reordered into best coincidence with ref_mol.

  • ref_mol (Union[Molecule, Molecule]) – Molecule to match.

  • atoms_map (bool) – Whether atom1 of ref_mol corresponds to atom1 of concern_mol, etc. If true, specifying True can save much time.

  • mols_align (bool) – Whether ref_mol and concern_mol have identical geometries by eye (barring orientation or atom mapping) and expected final RMSD = 0. If True, procedure is truncated when RMSD condition met, saving time.

  • do_plot (bool) – Pops up a mpl plot showing before, after, and ref geometries.

  • run_to_completion (bool) – Run reorderings to completion (past RMSD = 0) even if unnecessary because mols_align=True. Used to test worst-case timings.

  • run_resorting (bool) – Run the resorting machinery even if unnecessary because atoms_map=True.

  • uno_cutoff (float) – TODO

  • run_mirror (bool) – Run alternate geometries potentially allowing best match to ref_mol from mirror image of concern_mol. Only run if system confirmed to be nonsuperimposable upon mirror reflection.

  • verbose (int) –

Returns:

First item is RMSD [A] between ref_mol and the optimally aligned geometry computed. Second item is a AlignmentMill namedtuple with fields (shift, rotation, atommap, mirror) that prescribe the transformation from concern_mol and the optimally aligned geometry. Third item is a crude charge-, multiplicity-, fragment-less Molecule at optimally aligned (and atom-ordered) geometry. Return type determined by concern_mol type.

Return type:

float, tuple, qcdb.Molecule or psi4.core.Molecule

BFS(seed_atoms=None, bond_threshold=1.2, return_arrays=False, return_molecules=False, return_molecule=False)[source]

Detect fragments among real atoms through a breadth-first search (BFS) algorithm.

Parameters:
  • self (qcdb.Molecule or psi4.core.Molecule) –

  • seed_atoms (Optional[List]) – List of lists of atoms (0-indexed) belonging to independent fragments. Useful to prompt algorithm or to define intramolecular fragments through border atoms. Example: [[1, 0], [2]]

  • bond_threshold (float) – Factor beyond average of covalent radii to determine bond cutoff.

  • return_arrays (bool) – If True, also return fragments as list of arrays.

  • return_molecules (bool) – If True, also return fragments as list of Molecules.

  • return_molecule (bool) – If True, also return one big Molecule with fragmentation encoded.

Returns:

  • bfs_map (list of lists) – Array of atom indices (0-indexed) of detected fragments.

  • bfs_arrays (tuple of lists of ndarray, optional) – geom, mass, elem info per-fragment. Only provided if return_arrays is True.

  • bfs_molecules (list of qcdb.Molecule or psi4.core.Molecule, optional) – List of molecules, each built from one fragment. Center and orientation of fragments is fixed so orientation info from self is not lost. Loses chgmult and ghost/dummy info from self and contains default chgmult. Only provided if return_molecules is True. Returned are of same type as self.

  • bfs_molecule (qcdb.Molecule or psi4.core.Molecule, optional) – Single molecule with same number of real atoms as self with atoms reordered into adjacent fragments and fragment markers inserted. Loses ghost/dummy info from self; keeps total charge but not total mult. Only provided if return_molecule is True. Returned is of same type as self.

  • Authors

  • ——-

  • Original code from Michael S. Marshall, linear-scaling algorithm from

  • Trent M. Parker, revamped by Lori A. Burns

Notes

Relies upon van der Waals radii and so faulty for close (especially hydrogen-bonded) fragments. See` seed_atoms`. Any existing fragmentation info/chgmult encoded in self is lost.

Z(self: psi4.core.Molecule, arg0: int) float

Nuclear charge of atom arg0 (0-indexed without dummies)

activate_all_fragments(self: psi4.core.Molecule) None

Sets all fragments in the molecule to be active

add_atom(self: psi4.core.Molecule, Z: float, x: float, y: float, z: float, symbol: str, mass: float, charge: float, label: str, A: int) None

Adds to self Molecule an atom with atomic number Z, Cartesian coordinates in Bohr (x, y, z), atomic symbol, mass, and charge, extended atomic label, and mass number A

atom_at_position(*args, **kwargs)

Overloaded function.

  1. atom_at_position(self: psi4.core.Molecule, coord: float, tol: float) -> int

Returns the index of the atom inside tol radius around coord. Returns -1 for no atoms, throws an exception if more than one is found.

  1. atom_at_position(self: psi4.core.Molecule, coord: Annotated[List[float], FixedSize(3)], tol: float) -> int

Returns the index of the atom inside tol radius around coord. Returns -1 for no atoms, throws an exception if more than one is found.

basis_on_atom(self: psi4.core.Molecule, arg0: int) str

Gets the label of the orbital basis set on a given atom arg0

center_of_mass(self: psi4.core.Molecule) psi4.core.Vector3

Computes center of mass of molecule (does not translate molecule)

charge(self: psi4.core.Molecule, atom: int) float

Gets charge of atom (0-indexed without dummies)

clone(self: psi4.core.Molecule) psi4.core.Molecule

Returns a new Molecule identical to arg1

com_fixed(self: psi4.core.Molecule) bool

Gets whether or not center of mass is fixed

comment(self: psi4.core.Molecule) str

Gets molecule comment

connectivity(self: psi4.core.Molecule) List[Tuple[int, int, float]]

Gets molecule connectivity

create_psi4_string_from_molecule(self: psi4.core.Molecule) str

Gets a string re-expressing in input format the current state of the molecule.Contains Cartesian geometry info, fragmentation, charges and multiplicities, and any frame restriction.

deactivate_all_fragments(self: psi4.core.Molecule) None

Sets all fragments in the molecule to be inactive

distance_matrix(self: psi4.core.Molecule) psi4.core.Matrix

Returns Matrix of interatom distances

extract_subsets(*args, **kwargs)

Overloaded function.

  1. extract_subsets(self: psi4.core.Molecule, arg0: List[int], arg1: List[int]) -> psi4.core.Molecule

Returns copy of self with arg0 fragments Real and arg1 fragments Ghost

  1. extract_subsets(self: psi4.core.Molecule, arg0: List[int], arg1: int) -> psi4.core.Molecule

Returns copy of self with arg0 fragments Real and arg1 fragment Ghost

  1. extract_subsets(self: psi4.core.Molecule, arg0: int, arg1: List[int]) -> psi4.core.Molecule

Returns copy of self with arg0 fragment Real and arg1 fragments Ghost

  1. extract_subsets(self: psi4.core.Molecule, arg0: int, arg1: int) -> psi4.core.Molecule

Returns copy of self with arg0 fragment Real and arg1 fragment Ghost

  1. extract_subsets(self: psi4.core.Molecule, arg0: List[int]) -> psi4.core.Molecule

Returns copy of self with arg0 fragments Real

  1. extract_subsets(self: psi4.core.Molecule, arg0: int) -> psi4.core.Molecule

Returns copy of self with arg0 fragment Real

fZ(self: psi4.core.Molecule, arg0: int) float

Nuclear charge of atom arg1 (0-indexed including dummies)

fcharge(self: psi4.core.Molecule, atom: int) float

Gets charge of atom (0-indexed including dummies)

find_highest_point_group(self: psi4.core.Molecule, tolerance: float = 1e-08) psi4.core.PointGroup

Finds highest possible computational molecular point group

find_point_group(self: psi4.core.Molecule, tolerance: float = 1e-08) psi4.core.PointGroup

Finds computational molecular point group, user can override this with the symmetry keyword

fix_com(self: psi4.core.Molecule, arg0: bool) None

Sets whether to fix the Cartesian position, or to translate to the C.O.M. Expert use only; use before molecule finalized by update_geometry.

fix_orientation(self: psi4.core.Molecule, arg0: bool) None

Fix the orientation at its current frame. Expert use only; use before molecule finalized by update_geometry.

flabel(self: psi4.core.Molecule, atom: int) str

Gets the original label of the atom arg0 as given in the input file (C2, H4)(0-indexed including dummies)

fmass(self: psi4.core.Molecule, atom: int) float

Gets mass of atom (0-indexed including dummies)

form_symmetry_information(self: psi4.core.Molecule, arg0: float) None

Uses the point group object obtain by calling point_group()

format_molecule_for_mol()

Returns a string of Molecule formatted for mol2.

Written by Trent M. Parker 9 Jun 2014

classmethod from_arrays(geom=None, elea=None, elez=None, elem=None, mass=None, real=None, elbl=None, name=None, units='Angstrom', input_units_to_au=None, fix_com=None, fix_orientation=None, fix_symmetry=None, fragment_separators=None, fragment_charges=None, fragment_multiplicities=None, molecular_charge=None, molecular_multiplicity=None, comment=None, provenance=None, connectivity=None, missing_enabled_return='error', tooclose=0.1, zero_ghost_fragments=False, nonphysical=False, mtol=0.001, verbose=1, return_dict=False)

Construct Molecule from unvalidated arrays and variables.

Light wrapper around from_arrays() that is a full-featured constructor to dictionary representa- tion of Molecule. This follows one step further to return Molecule instance.

:param See from_arrays().:

Return type:

psi4.core.Molecule

static from_dict(arg0: dict) psi4.core.Molecule

Returns a new Molecule constructed from python dictionary. In progress: name and capabilities should not be relied upon

classmethod from_schema(molschema, return_dict=False, nonphysical=False, verbose=1)

Construct Molecule from non-Psi4 schema.

Light wrapper around from_arrays().

Parameters:
  • molschema (Dict) – Dictionary form of Molecule following known schema.

  • return_dict (bool) – Additionally return Molecule dictionary intermediate.

  • nonphysical (bool) – Do allow masses outside an element’s natural range to pass validation?

  • verbose (int) – Amount of printing.

Return type:

Union[Molecule, Tuple[Molecule, Dict]]

Returns:

  • mol (psi4.core.Molecule)

  • molrec (dict) – Dictionary representation of instance. Only provided if return_dict is True.

classmethod from_string(molstr, dtype=None, name=None, fix_com=None, fix_orientation=None, fix_symmetry=None, return_dict=False, enable_qm=True, enable_efp=True, missing_enabled_return_qm='none', missing_enabled_return_efp='none', verbose=1)
fsymbol(self: psi4.core.Molecule, atom: int) str

Gets the cleaned up label of atom (C2 => C, H4 = H) (0-indexed including dummies)

ftrue_atomic_number(self: psi4.core.Molecule, atom: int) int

Gets atomic number of atom from element (0-indexed including dummies)

full_geometry(self: psi4.core.Molecule) psi4.core.Matrix

Gets the geometry [Bohr] as a (Natom X 3) matrix of coordinates (including dummies)

full_pg_n(self: psi4.core.Molecule) int

Gets n in Cnv, etc.; If there is no n (e.g. Td) it’s the highest-order rotation axis

fx(self: psi4.core.Molecule, arg0: int) float

x position of atom arg0 (0-indexed including dummies in Bohr)

fy(self: psi4.core.Molecule, arg0: int) float

y position of atom arg0 (0-indexed including dummies in Bohr)

fz(self: psi4.core.Molecule, arg0: int) float

z position of atom arg0 (0-indexed including dummies in Bohr)

geometry(self: psi4.core.Molecule) psi4.core.Matrix

Gets the geometry [Bohr] as a (Natom X 3) matrix of coordinates (excluding dummies)

get_fragment_charges(self: psi4.core.Molecule) List[int]

Gets the charge of each fragment

get_fragment_multiplicities(self: psi4.core.Molecule) List[int]

Gets the multiplicity of each fragment

get_fragment_types(self: psi4.core.Molecule) List[str]

Returns a list describing how to handle each fragment {Real, Ghost, Absent}

get_fragments(self: psi4.core.Molecule) List[Tuple[int, int]]

Returns list of pairs of atom ranges defining each fragment from parent molecule(fragments[frag_ind] = <Afirst,Alast+1>)

get_full_point_group(self: psi4.core.Molecule) str

Gets point group name such as C3v or S8

get_full_point_group_with_n(self: psi4.core.Molecule) str

Gets point group name such as Cnv or Sn

get_variable(self: psi4.core.Molecule, arg0: str) float

Returns the value of variable arg0 in the structural variables list

has_zmatrix(self: psi4.core.Molecule) bool

Get whether or not this molecule has at least one zmatrix entry

inertia_tensor(self: psi4.core.Molecule) psi4.core.Matrix

Returns intertial tensor

input_units_to_au(self: psi4.core.Molecule) float

Returns unit conversion to [a0] for geometry

irrep_labels(self: psi4.core.Molecule) List[str]

Returns Irreducible Representation symmetry labels

is_variable(self: psi4.core.Molecule, arg0: str) bool

Checks if variable arg0 is in the structural variables list

label(self: psi4.core.Molecule, atom: int) str

Gets the original label of the atom as given in the input file (C2, H4)(0-indexed without dummies)

mass(self: psi4.core.Molecule, atom: int) float

Returns mass of atom (0-indexed)

mass_number(self: psi4.core.Molecule, arg0: int) int

Mass number (A) of atom if known, else -1

molecular_charge(self: psi4.core.Molecule) int

Gets the molecular charge

move_to_com(self: psi4.core.Molecule) None

Moves molecule to center of mass

multiplicity(self: psi4.core.Molecule) int

Gets the multiplicity (defined as 2Ms + 1)

nallatom(self: psi4.core.Molecule) int

Number of real and dummy atoms

name(self: psi4.core.Molecule) str

Gets molecule name

natom(self: psi4.core.Molecule) int

Number of real atoms

nfragments(self: psi4.core.Molecule) int

Gets the number of fragments in the molecule

nuclear_dipole(*args, **kwargs)

Overloaded function.

  1. nuclear_dipole(self: psi4.core.Molecule, arg0: psi4.core.Vector3) -> psi4.core.Vector3

Gets the nuclear contribution to the dipole, with respect to a specified origin atg0

  1. nuclear_dipole(self: psi4.core.Molecule) -> psi4.core.Vector3

Gets the nuclear contribution to the dipole, with respect to the origin

nuclear_repulsion_energy(self: psi4.core.Molecule, dipole_field: Annotated[List[float], FixedSize(3)] = [0.0, 0.0, 0.0]) float

Computes nuclear repulsion energy

nuclear_repulsion_energy_deriv1(self: psi4.core.Molecule, dipole_field: Annotated[List[float], FixedSize(3)] = [0.0, 0.0, 0.0]) psi4.core.Matrix

Returns first derivative of nuclear repulsion energy as a matrix (natom, 3)

nuclear_repulsion_energy_deriv2(self: psi4.core.Molecule) psi4.core.Matrix

Returns second derivative of nuclear repulsion energy as a matrix (natom X 3, natom X 3)

orientation_fixed(self: psi4.core.Molecule) bool

Get whether or not orientation is fixed

point_group(self: psi4.core.Molecule) psi4.core.PointGroup

Returns the current point group object

print_bond_angles(self: psi4.core.Molecule) None

Print the bond angle geometrical parameters

print_cluster(self: psi4.core.Molecule) None

Prints the molecule in Cartesians in input units adding fragment separators

print_distances(self: psi4.core.Molecule) None

Print the interatomic distance geometrical parameters

print_in_input_format(self: psi4.core.Molecule) None

Prints the molecule as Cartesian or ZMatrix entries, just as inputted.

print_out(self: psi4.core.Molecule) None

Prints the molecule in Cartesians in input units to output file

print_out_in_angstrom(self: psi4.core.Molecule) None

Prints the molecule in Cartesians in Angstroms to output file

print_out_in_bohr(self: psi4.core.Molecule) None

Prints the molecule in Cartesians in Bohr to output file

print_out_of_planes(self: psi4.core.Molecule) None

Print the out-of-plane angle geometrical parameters to output file

print_rotational_constants(self: psi4.core.Molecule) None

Print the rotational constants to output file

provenance(self: psi4.core.Molecule) Dict[str, str]

Gets molecule provenance

reinterpret_coordentry(self: psi4.core.Molecule, arg0: bool) None

Do reinterpret coordinate entries during update_geometry().

reset_point_group(self: psi4.core.Molecule, arg0: str) None

Overrides symmetry from outside the molecule string

rotational_constants(self: psi4.core.Molecule) psi4.core.Vector

Returns the rotational constants [cm^-1] of the molecule

rotational_symmetry_number(self: psi4.core.Molecule) int

Returns number of unique orientations of the rigid molecule that only interchange identical atoms

rotor_type(self: psi4.core.Molecule) str

Returns rotor type, e.g. ‘RT_ATOM’ or ‘RT_SYMMETRIC_TOP’

run_dftd3(func=None, dashlvl=None, dashparam=None, dertype=None, verbose=1)[source]

Compute dispersion correction via Grimme’s DFTD3 program.

Parameters:
  • func (Optional[str]) – Name of functional (func only, func & disp, or disp only) for which to compute dispersion (e.g., blyp, BLYP-D2, blyp-d3bj, blyp-d3(bj), hf+d). Any or all parameters initialized from dashcoeff[dashlvl][func] can be overwritten via dashparam.

  • dashlvl (Optional[str]) – Name of dispersion correction to be applied (e.g., d, D2, d3(bj), das2010). Must be key in dashcoeff or “alias” or “formal” to run.

  • dashparam (Optional[Dict]) – Values for the same keys as dashcoeff[dashlvl][‘default’] used to override any or all values initialized by func. Extra parameters will error.

  • dertype (Union[int, str, None]) – Maximum derivative level at which to run DFTD3. For large molecules, energy-only calculations can be significantly more efficient. Influences return values, see below.

  • verbose (int) – Amount of printing.

Returns:

  • energy (float) – When dertype=0, energy [Eh].

  • gradient (~numpy.ndarray) – When dertype=1, (nat, 3) gradient [Eh/a0].

  • (energy, gradient) (tuple of float and ~numpy.ndarray) – When dertype=None, both energy [Eh] and (nat, 3) gradient [Eh/a0].

run_dftd4(func=None, dashlvl=None, dashparam=None, dertype=None, verbose=1)[source]

Compute dispersion correction via Grimme’s DFTD4 program.

Parameters:
  • func (Optional[str]) – Name of functional (func only, func & disp, or disp only) for which to compute dispersion (e.g., blyp, BLYP-D2, blyp-d3bj, blyp-d3(bj), hf+d). Unlike run_dftd3, func overwrites any parameter initialized via dashparam.

  • dashlvl (Optional[str]) – Name of dispersion correction to be applied (e.g., d, D2, d3(bj), das2010). Must be key in dashcoeff or “alias” or “formal” to run.

  • dashparam (Optional[Dict]) – Values for the same keys as dashcoeff[dashlvl][‘default’] used to provide custom values. Unlike run_dftd3, will not have effect if func given. Must provide all parameters. Extra parameters will error.

  • dertype (Union[int, str, None]) – Maximum derivative level at which to run DFTD3. For large molecules, energy-only calculations can be significantly more efficient. Influences return values, see below.

  • verbose (int) – Amount of printing.

Returns:

  • energy (float) – When dertype=0, energy [Eh].

  • gradient (ndarray) – When dertype=1, (nat, 3) gradient [Eh/a0].

  • (energy, gradient) (tuple of float and ndarray) – When dertype=None, both energy [Eh] and (nat, 3) gradient [Eh/a0].

Notes

This function wraps the QCEngine dftd4 harness which wraps the internal DFTD4 Python API. As such, the upstream convention of func trumping dashparam holds, rather than the run_dftd3() behavior of dashparam extending or overriding func.

run_gcp(func=None, dertype=None, verbose=1)[source]

Compute geometrical BSSE correction via Grimme’s GCP program.

Function to call Grimme’s GCP program https://www.chemie.uni-bonn.de/pctc/mulliken-center/software/gcp/gcp to compute an a posteriori geometrical BSSE correction to self for several HF, generic DFT, and specific HF-3c and PBEh-3c method/basis combinations, func. Returns energy if dertype is 0, gradient if dertype is 1, else tuple of energy and gradient if dertype unspecified. The gcp executable must be independently compiled and found in PATH or PSIPATH. self may be either a qcdb.Molecule (sensibly) or a psi4.Molecule (works b/c psi4.Molecule has been extended by this method py-side and only public interface fns used) or a string that can be instantiated into a qcdb.Molecule.

Parameters:
  • func (Optional[str]) – Name of method/basis combination or composite method for which to compute the correction (e.g., HF/cc-pVDZ, DFT/def2-SVP, HF3c, PBEh3c).

  • dertype (Union[int, str, None]) – Maximum derivative level at which to run GCP. For large molecules, energy-only calculations can be significantly more efficient. Influences return values, see below.

  • verbose (int) – Amount of printing. Unused at present.

Returns:

  • energy (float) – When dertype=0, energy [Eh].

  • gradient (ndarray) – When dertype=1, (nat, 3) gradient [Eh/a0].

  • (energy, gradient) (tuple of float and ndarray) – When dertype=None, both energy [Eh] and (nat, 3) gradient [Eh/a0].

save_string_xyz(self: psi4.core.Molecule) str

Saves the string of an XYZ file to arg2

save_string_xyz_file(self: psi4.core.Molecule) str

Saves an XYZ file to arg2

save_xyz_file(self: psi4.core.Molecule, arg0: str, arg1: bool) None

Saves an XYZ file to arg0

schoenflies_symbol(self: psi4.core.Molecule) str

Returns the Schoenflies symbol

scramble(do_shift=True, do_rotate=True, do_resort=True, deflection=1.0, do_mirror=False, do_plot=False, run_to_completion=False, run_resorting=False, verbose=1)[source]

Tester for B787 by shifting, rotating, and atom shuffling ref_mol and checking that the aligner returns the opposite transformation.

Parameters:
  • ref_mol (Molecule) – Molecule to perturb.

  • do_shift (Union[bool, ndarray, List]) – Whether to generate a random atom shift on interval [-3, 3) in each dimension (True) or leave at current origin. To shift by a specified vector, supply a 3-element list.

  • do_rotate (Union[bool, ndarray, List[List]]) – Whether to generate a random 3D rotation according to algorithm of Arvo. To rotate by a specified matrix, supply a 9-element list of lists.

  • do_resort (Union[bool, List]) – Whether to shuffle atoms (True) or leave 1st atom 1st, etc. (False). To specify shuffle, supply a nat-element list of indices.

  • deflection (float) – If do_rotate, how random a rotation: 0.0 is no change, 0.1 is small perturbation, 1.0 is completely random.

  • do_mirror (bool) – Whether to construct the mirror image structure by inverting y-axis.

  • do_plot (bool) – Pops up a mpl plot showing before, after, and ref geometries.

  • run_to_completion (bool) – By construction, scrambled systems are fully alignable (final RMSD=0). Even so, True turns off the mechanism to stop when RMSD reaches zero and instead proceed to worst possible time.

  • run_resorting (bool) – Even if atoms not shuffled, test the resorting machinery.

  • verbose (int) – Print level.

Return type:

None

set_active_fragment(self: psi4.core.Molecule, arg0: int) None

Sets the specified fragment arg0 to be Real

set_active_fragments(self: psi4.core.Molecule, arg0: List[int]) None

Sets the specified list arg0 of fragments to be Real

set_basis_all_atoms(self: psi4.core.Molecule, arg0: str, arg1: str) None

Sets basis set arg0 to all atoms

set_basis_by_label(self: psi4.core.Molecule, arg0: str, arg1: str, arg2: str) None

Sets basis set arg1 to all atoms with label (e.g., H4) arg0

set_basis_by_number(self: psi4.core.Molecule, arg0: int, arg1: str, arg2: str) None

Sets basis set arg1 to all atoms with number arg0

set_basis_by_symbol(self: psi4.core.Molecule, arg0: str, arg1: str, arg2: str) None

Sets basis set arg1 to all atoms with symbol (e.g., H) arg0

set_comment(self: psi4.core.Molecule, arg0: str) None

Sets molecule comment

set_connectivity(self: psi4.core.Molecule, arg0: List[Tuple[int, int, float]]) None

Sets molecule connectivity

set_full_geometry(self: psi4.core.Molecule, arg0: psi4.core.Matrix) None

Sets the geometry, given a (Natom X 3) matrix arg0 of coordinates (in Bohr) (including dummies

set_geometry(self: psi4.core.Molecule, arg0: psi4.core.Matrix) None

Sets the geometry, given a (Natom X 3) matrix arg0 of coordinates [a0] (excluding dummies)

set_ghost_fragment(self: psi4.core.Molecule, arg0: int) None

Sets the specified fragment arg0 to be Ghost

set_ghost_fragments(self: psi4.core.Molecule, arg0: List[int]) None

Sets the specified list arg0 of fragments to be Ghost

set_input_units_to_au(self: psi4.core.Molecule, arg0: float) None

Sets unit conversion to [a0] for geometry

set_mass(self: psi4.core.Molecule, atom: int, mass: float) None

Sets mass of atom (0-indexed) to mass (good for isotopic substitutions)

set_molecular_charge(self: psi4.core.Molecule, arg0: int) None

Change the overall molecular charge. Setting in initial molecule string or constructor preferred.

set_multiplicity(self: psi4.core.Molecule, arg0: int) None

Change the multiplicity (defined as 2S + 1). Setting in initial molecule string or constructor preferred.

set_name(self: psi4.core.Molecule, arg0: str) None

Sets molecule name

set_nuclear_charge(self: psi4.core.Molecule, arg0: int, arg1: float) None

Set the nuclear charge of the given atom arg0 to the value arg1 (primarily for ECP).

set_point_group(self: psi4.core.Molecule, arg0: psi4.core.PointGroup) None

Sets the molecular point group to the point group object arg0

set_provenance(self: psi4.core.Molecule, arg0: Dict[str, str]) None

Sets molecule provenance

set_units(self: psi4.core.Molecule, arg0: psi4.core.GeometryUnits) None

Sets units (Angstrom or Bohr) used to define the geometry. Imposes Psi4 physical constants conversion for input_units_to_au.

set_variable(self: psi4.core.Molecule, arg0: str, arg1: float) None

Sets the value arg1 to the variable arg0 in the list of structure variables, then calls update_geometry()

symbol(self: psi4.core.Molecule, atom: int) str

Gets the cleaned up label of atom (C2 => C, H4 = H) (0-indexed without dummies)

symmetrize(self: psi4.core.Molecule, arg0: float) None

Finds the highest point Abelian point group within the specified tolerance, and forces the geometry to have that symmetry.

symmetry_from_input(self: psi4.core.Molecule) str

Returns the symmetry specified in the input

to_arrays(dummy=False, ghost_as_dummy=False)[source]

Exports coordinate info into NumPy arrays.

Parameters:
  • dummy (bool) – Whether or not to include dummy atoms in returned arrays.

  • ghost_as_dummy (bool) – Whether or not to treat ghost atoms as dummies.

Return type:

Tuple[ndarray, ndarray, ndarray, ndarray, ndarray]

Returns:

  • geom, mass, elem, elez, uniq (numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray) – (nat, 3) geometry [a0]. (nat,) mass [u]. (nat,) element symbol. (nat,) atomic number. (nat,) hash of element symbol and mass. Note that coordinate, orientation, and element information is preserved but fragmentation, chgmult, and dummy/ghost is lost.

  • Usage

  • —–

  • geom, mass, elem, elez, uniq = molinstance.to_arrays()

to_dict(force_c1=False, force_units=False, np_out=True)[source]

Serializes instance into Molecule dictionary.

to_schema(dtype, units='Bohr')[source]

Serializes instance into dictionary according to schema dtype.

to_string(dtype, units=None, atom_format=None, ghost_format=None, width=17, prec=12)[source]

Format a string representation of QM molecule.

translate(self: psi4.core.Molecule, arg0: psi4.core.Vector3) None

Translates molecule by arg0

true_atomic_number(self: psi4.core.Molecule, atom: int) int

Gets atomic number of atom from element (0-indexed without dummies)

units(self: psi4.core.Molecule) str

Returns units used to define the geometry, i.e. ‘Angstrom’ or ‘Bohr’

update_geometry(self: psi4.core.Molecule) None

Reevaluates the geometry with current variable values, orientation directives, etc. by clearing the atoms list and rebuilding it. Idempotent. Use liberally.Must be called after initial Molecule definition by string.

x(self: psi4.core.Molecule, arg0: int) float

x position [Bohr] of atom arg0 (0-indexed without dummies)

xyz(self: psi4.core.Molecule, i: int) psi4.core.Vector3

Return the Vector3 for atom i (0-indexed without dummies)

y(self: psi4.core.Molecule, arg0: int) float

y position [Bohr] of atom arg0 (0-indexed without dummies)

z(self: psi4.core.Molecule, arg0: int) float

z position [Bohr] of atom arg0 (0-indexed without dummies)