{ "cells": [ { "cell_type": "code", "execution_count": 1, "id": "corrected-cheese", "metadata": { "tags": [ "remove_cell" ] }, "outputs": [], "source": [ "\"\"\"First steps with VeloxChem\"\"\"\n", "\n", "__author__ = \"Roberto Di Remigio\"\n", "__credit__ = [\"Roberto Di Remigio\", \"Xin Li\"]\n", "\n", "__copyright__ = \"(c) 2021, ENCSS and PDC\"\n", "__license__ = \"MIT\"\n", "__date__ = \"2021-04-13\"" ] }, { "cell_type": "markdown", "id": "organized-toolbox", "metadata": { "tags": [ "remove_cell" ] }, "source": [ "
\n", " \n", "
" ] }, { "cell_type": "markdown", "id": "peripheral-proportion", "metadata": { "tags": [] }, "source": [ "# First steps with VeloxChem" ] }, { "cell_type": "markdown", "id": "wicked-madison", "metadata": { "tags": [ "remove_cell" ] }, "source": [ "
\n", " \n", " \n", "
\n", " To run the selected code cell, hit
Shift + Enter
\n", "
\n", "
" ] }, { "cell_type": "markdown", "id": "heavy-mortgage", "metadata": {}, "source": [ "We want to run a restricted Hartree-Fock test calculation on a water molecule. We will use the following structure." ] }, { "cell_type": "code", "execution_count": 2, "id": "crude-norwegian", "metadata": {}, "outputs": [ { "data": { "application/3dmoljs_load.v0": "
\n

You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n jupyter labextension install jupyterlab_3dmol

\n
\n", "text/html": [ "
\n", "

You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n", " jupyter labextension install jupyterlab_3dmol

\n", "
\n", "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "import py3Dmol as p3d\n", "\n", "h2o_xyz = \"\"\"3\n", "water \n", "O 0.000000000000 0.000000000000 0.000000000000 \n", "H 0.000000000000 0.740848095288 0.582094932012 \n", "H 0.000000000000 -0.740848095288 0.582094932012\n", "\"\"\"\n", "\n", "v = p3d.view(width=200, height=200)\n", "\n", "v.addModel(h2o_xyz, \"xyz\")\n", "v.setStyle({'stick':{}})\n", "v.show()" ] }, { "cell_type": "markdown", "id": "fixed-variance", "metadata": {}, "source": [ "## Running with an input file\n", "\n", "Let us save the following as `water.inp`:\n", "\n", "```\n", "@jobs\n", "task: scf\n", "@end\n", "\n", "@method settings\n", "basis: aug-cc-pvdz\n", "@end\n", "\n", "@molecule\n", "charge: 0\n", "multiplicity: 1\n", "xyz: \n", "O 0.000000000000 0.000000000000 0.000000000000 \n", "H 0.000000000000 0.740848095288 0.582094932012 \n", "H 0.000000000000 -0.740848095288 0.582094932012\n", "@end \n", "```\n", "\n", "A VeloxChem input file has:\n", "- *sections*. Each section starts with`@` and is closed with `@end`,\n", "- *keywords*. You specify them with `: `. \n", "\n", "Documentation for all input options is available at: https://docs.veloxchem.org/inputs/keywords.html\n", "\n", "### The command-line interface\n", "\n", "Possibly the most common way of running a quantum chemistry code is by using a command-line interface (CLI). You can access the VeloxChem CLI with:\n", "```\n", "$ vlx --help\n", "```\n", "this will print the following help text:\n", "```\n", "usage: \n", "================= VeloxChem =================\n", "\n", "vlx input_file [output_file]\n", "\n", "or\n", "\n", "python -m veloxchem input_file [output_file]\n", " \n", "\n", "positional arguments:\n", " input_file Input file\n", " output_file Output file (default: STDOUT)\n", "\n", "optional arguments:\n", " -h, --help show this help message and exit\n", " -v, --version show program's version number and exit\n", "```\n", "\n", "Let's run our water calculation with:\n", "```\n", "$ vlx water.inp water.out\n", "```\n", "such that the output will be saved to the plain-text file `water.out`. Congratulations! You have completed your first quantum chemical calculation with VeloxChem!\n", "\n", "### The output file\n", " \n", "Let us look now at the information contained in the output file. VeloxChem logs in the output stream the operations performed during a calculation. You will notice:\n", "\n", "- The molecular geometry\n", "```\n", " Molecular Geometry (Angstroms)\n", " ================================\n", "\n", " Atom Coordinate X Coordinate Y Coordinate Z\n", "\n", " O 0.000000000000 0.000000000000 0.000000000000\n", " H 0.000000000000 0.740848095288 0.582094932012\n", " H 0.000000000000 -0.740848095288 0.582094932012\n", "\n", " Molecular charge : 0\n", " Spin multiplicity : 1\n", " Number of atoms : 3\n", " Number of alpha electrons : 5\n", " Number of beta electrons : 5\n", "```\n", "- The molecular basis and self-consistent field setup\n", "\n", "```\n", " Molecular Basis (Atomic Basis)\n", " ================================\n", "\n", " Basis: AUG-CC-PVDZ\n", "\n", " Atom Contracted GTOs Primitive GTOs\n", "\n", " H (3S,2P) (5S,2P)\n", " O (4S,3P,2D) (18S,5P,2D)\n", "\n", " Contracted Basis Functions : 41\n", " Primitive Basis Functions : 65\n", "\n", "\n", " Self Consistent Field Driver Setup\n", " ====================================\n", "\n", " Wave Function Model : Spin-Restricted Hartree-Fock\n", " Initial Guess Model : Superposition of Atomic Densities\n", " Convergence Accelerator : Two Level Direct Inversion of Iterative Subspace\n", " Max. Number of Iterations : 50\n", " Max. Number of Error Vectors : 10\n", " Convergence Threshold : 1.0e-06\n", " ERI Screening Scheme : Cauchy Schwarz + Density\n", " ERI Screening Mode : Dynamic\n", " ERI Screening Threshold : 1.0e-12\n", " Linear Dependence Threshold : 1.0e-06\n", "```\n", "- The report on the iterative SCF procedure\n", "\n", "```\n", " Iter. | Hartree-Fock Energy | Energy Change | Gradient Norm | Max. Gradient | Density Change\n", " --------------------------------------------------------------------------------------------\n", " 1 -76.040225772729 0.0000000000 0.10527360 0.01410687 0.00000000\n", " 2 -76.041581009463 -0.0013552367 0.02396700 0.00364108 0.03764551\n", " 3 -76.041683167061 -0.0001021576 0.00872358 0.00122472 0.00801861\n", " 4 -76.041694355980 -0.0000111889 0.00512350 0.00094645 0.00386334\n", " 5 -76.041697479316 -0.0000031233 0.00046249 0.00008300 0.00125597\n", " 6 -76.041697547982 -0.0000000687 0.00007707 0.00001089 0.00028877\n", " 7 -76.041697549783 -0.0000000018 0.00000995 0.00000158 0.00004453\n", " 8 -76.041697549817 -0.0000000000 0.00000212 0.00000032 0.00000530\n", " 9 -76.041697549819 -0.0000000000 0.00000039 0.00000005 0.00000148\n", "```\n", "- The final summary of SCF energies\n", "\n", "```\n", " Spin-Restricted Hartree-Fock:\n", " -----------------------------\n", " Total Energy : -76.0416975498 au\n", " Electronic Energy : -85.3853357078 au\n", " Nuclear Repulsion Energy : 9.3436381580 au\n", " ------------------------------------\n", " Gradient Norm : 0.0000003901 au\n", "\n", "\n", " Ground State Information\n", " ------------------------\n", " Charge of Molecule : 0.0\n", " Multiplicity (2S+1) : 1.0\n", " Magnetic Quantum Number (M_S) : 0.0\n", "```\n", "- The molecular orbitals and the electric dipole moment\n", "- Last, but not least, the preferred citation for VeloxChem\n", "\n", "```\n", "!========================================================================================================================!\n", "! Rinkevicius, Z.; Li, X.; Vahtras, O.; Ahmadzadeh, K.; Brand, M.; Ringholm, M.; !\n", "! List, N. H.; Scheurer, M.; Scott, M.; Dreuw, A.; Norman, P. !\n", "! VeloxChem: A Python-driven Density-functional Theory Program for Spectroscopy !\n", "! Simulations in High-performance Computing Environments. !\n", "! WIREs Comput Mol Sci 2020, 10 (5), e1457. !\n", "!========================================================================================================================!\n", "```" ] }, { "cell_type": "markdown", "id": "trying-finance", "metadata": {}, "source": [ "## Running with a Python script\n", "\n", "Extracting information from a plain-text output file is not always desirable, for a number of reasons:\n", "1. The format of the output is not guaranteed to stay the same between different versions of the code. Data-harvesting scripts might fail unexpectedly.\n", "2. Reading floating-point numbers from their text representation might introduce subtle round-off errors.\n", "\n", "VeloxChem offers a powerful Python application programming interface (API) to *programmatically* access computed data, for example, SCF ground state electric dipole moments, **and** computational primitives, such as AO basis integrals. Let's see how we can write a Python script to perform a calculation on the water molecule without using plain-text input and output files. First and foremost, we import the VeloxChem API:" ] }, { "cell_type": "code", "execution_count": 3, "id": "olive-region", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "* Warning * Environment variable OMP_NUM_THREADS not set.\n", "* Warning * Setting OMP_NUM_THREADS to 8.\n" ] } ], "source": [ "import veloxchem as vlx" ] }, { "cell_type": "markdown", "id": "bacterial-characteristic", "metadata": {}, "source": [ "### The `Molecule` object\n", "\n", "The next task is to create a `Molecule` object. We can use the `h2o_xyz` string defined above directly." ] }, { "cell_type": "code", "execution_count": 4, "id": "welcome-restaurant", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Molecular Geometry (Angstroms)\n", "================================\n", "\n", " Atom Coordinate X Coordinate Y Coordinate Z \n", "\n", " O 0.000000000000 0.000000000000 0.000000000000\n", " H 0.000000000000 0.740848095288 0.582094932012\n", " H 0.000000000000 -0.740848095288 0.582094932012\n", "\n", "\n" ] } ], "source": [ "h2o = vlx.Molecule.from_xyz_string(h2o_xyz)\n", "print(h2o.get_string())" ] }, { "cell_type": "markdown", "id": "circular-ribbon", "metadata": {}, "source": [ "The `help` native method of Python works as well and we can explore the VeloxChem API easily. For example, `h2o.get_coordinates()` will return the molecular structure as a [NumPy array](https://numpy.org/doc/stable/reference/arrays.html)." ] }, { "cell_type": "code", "execution_count": 5, "id": "subjective-garbage", "metadata": { "tags": [ "hide_output" ] }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Help on class Molecule in module veloxchem.veloxchemlib:\n", "\n", "class Molecule(pybind11_builtins.pybind11_object)\n", " | Method resolution order:\n", " | Molecule\n", " | pybind11_builtins.pybind11_object\n", " | builtins.object\n", " | \n", " | Methods defined here:\n", " | \n", " | __eq__(...)\n", " | __eq__(self: veloxchem.veloxchemlib.Molecule, arg0: veloxchem.veloxchemlib.Molecule) -> bool\n", " | \n", " | __init__(...)\n", " | __init__(*args, **kwargs)\n", " | Overloaded function.\n", " | \n", " | 1. __init__(self: veloxchem.veloxchemlib.Molecule) -> None\n", " | \n", " | 2. __init__(self: veloxchem.veloxchemlib.Molecule, arg0: veloxchem.veloxchemlib.Molecule) -> None\n", " | \n", " | 3. __init__(self: veloxchem.veloxchemlib.Molecule, arg0: veloxchem.veloxchemlib.Molecule, arg1: veloxchem.veloxchemlib.Molecule) -> None\n", " | \n", " | 4. __init__(self: veloxchem.veloxchemlib.Molecule, symbols: List[str], coordinates: numpy.ndarray[numpy.float64], units: str = 'angstrom') -> None\n", " | \n", " | 5. __init__(self: veloxchem.veloxchemlib.Molecule, Zs: List[int], coordinates: numpy.ndarray[numpy.float64], units: str = 'angstrom') -> None\n", " | \n", " | broadcast(...)\n", " | broadcast(self: veloxchem.veloxchemlib.Molecule, arg0: int, arg1: object) -> None\n", " | \n", " | center_of_mass = _Molecule_center_of_mass(self)\n", " | Computes center of mass of a molecule.\n", " | \n", " | :return:\n", " | The center of mass.\n", " | \n", " | check_multiplicity(...)\n", " | check_multiplicity(self: veloxchem.veloxchemlib.Molecule) -> None\n", " | \n", " | check_proximity(...)\n", " | check_proximity(self: veloxchem.veloxchemlib.Molecule, arg0: float) -> None\n", " | \n", " | coordination_numbers(...)\n", " | coordination_numbers(self: veloxchem.veloxchemlib.Molecule) -> numpy.ndarray[numpy.float64]\n", " | \n", " | elem_ids_to_numpy(...)\n", " | elem_ids_to_numpy(self: veloxchem.veloxchemlib.Molecule) -> numpy.ndarray[numpy.int32]\n", " | \n", " | get_charge(...)\n", " | get_charge(self: veloxchem.veloxchemlib.Molecule) -> float\n", " | \n", " | get_coordinates = _Molecule_get_coordinates(self)\n", " | Returns atom coordinates in atomic units.\n", " | \n", " | :return:\n", " | A numpy array of atom coordinates (nx3) in atomic units.\n", " | \n", " | get_elemental_composition(...)\n", " | get_elemental_composition(self: veloxchem.veloxchemlib.Molecule) -> list\n", " | \n", " | get_ic_rmsd = _Molecule_get_ic_rmsd(self, ref_mol)\n", " | Gets statistical deviation of bonds, angles and dihedral angles between\n", " | self and reference geometry.\n", " | \n", " | :param ref_mol:\n", " | The reference molecule (or xyz filename).\n", " | \n", " | get_labels = _Molecule_get_labels(self)\n", " | Returns atom labels.\n", " | \n", " | :return:\n", " | A list of atom labels.\n", " | \n", " | get_multiplicity(...)\n", " | get_multiplicity(self: veloxchem.veloxchemlib.Molecule) -> int\n", " | \n", " | get_string(...)\n", " | get_string(self: veloxchem.veloxchemlib.Molecule) -> str\n", " | \n", " | get_sub_molecule(...)\n", " | get_sub_molecule(self: veloxchem.veloxchemlib.Molecule, arg0: int, arg1: int) -> veloxchem.veloxchemlib.Molecule\n", " | \n", " | masses_to_numpy(...)\n", " | masses_to_numpy(self: veloxchem.veloxchemlib.Molecule) -> numpy.ndarray[numpy.float64]\n", " | \n", " | more_info = _Molecule_more_info(self)\n", " | Returns more information about the molecule.\n", " | \n", " | :return:\n", " | Molecular information in plain text.\n", " | \n", " | nuclear_repulsion_energy(...)\n", " | nuclear_repulsion_energy(self: veloxchem.veloxchemlib.Molecule) -> float\n", " | \n", " | number_of_alpha_electrons(...)\n", " | number_of_alpha_electrons(self: veloxchem.veloxchemlib.Molecule) -> int\n", " | \n", " | number_of_atoms(...)\n", " | number_of_atoms(*args, **kwargs)\n", " | Overloaded function.\n", " | \n", " | 1. number_of_atoms(self: veloxchem.veloxchemlib.Molecule) -> int\n", " | \n", " | 2. number_of_atoms(self: veloxchem.veloxchemlib.Molecule, idElemental: int) -> int\n", " | \n", " | 3. number_of_atoms(self: veloxchem.veloxchemlib.Molecule, iatom: int, natoms: int, idElemental: int) -> int\n", " | \n", " | number_of_beta_electrons(...)\n", " | number_of_beta_electrons(self: veloxchem.veloxchemlib.Molecule) -> int\n", " | \n", " | number_of_electrons(...)\n", " | number_of_electrons(self: veloxchem.veloxchemlib.Molecule) -> int\n", " | \n", " | partial_charges(...)\n", " | partial_charges(self: veloxchem.veloxchemlib.Molecule) -> numpy.ndarray[numpy.float64]\n", " | \n", " | set_charge(...)\n", " | set_charge(self: veloxchem.veloxchemlib.Molecule, arg0: float) -> None\n", " | \n", " | set_multiplicity(...)\n", " | set_multiplicity(self: veloxchem.veloxchemlib.Molecule, arg0: int) -> None\n", " | \n", " | vdw_radii_to_numpy(...)\n", " | vdw_radii_to_numpy(self: veloxchem.veloxchemlib.Molecule) -> numpy.ndarray[numpy.float64]\n", " | \n", " | write_xyz = _Molecule_write_xyz(self, xyz_filename)\n", " | Writes molecular geometry to xyz file.\n", " | \n", " | :param xyz_filename:\n", " | The name of the xyz file.\n", " | \n", " | x_to_numpy(...)\n", " | x_to_numpy(self: veloxchem.veloxchemlib.Molecule) -> numpy.ndarray[numpy.float64]\n", " | \n", " | y_to_numpy(...)\n", " | y_to_numpy(self: veloxchem.veloxchemlib.Molecule) -> numpy.ndarray[numpy.float64]\n", " | \n", " | z_to_numpy(...)\n", " | z_to_numpy(self: veloxchem.veloxchemlib.Molecule) -> numpy.ndarray[numpy.float64]\n", " | \n", " | ----------------------------------------------------------------------\n", " | Static methods defined here:\n", " | \n", " | from_dict = _Molecule_from_dict(mol_dict)\n", " | Reads molecule from a dictionary.\n", " | \n", " | :param mol_dict:\n", " | The molecule dictionary.\n", " | \n", " | :return:\n", " | The molecule.\n", " | \n", " | from_xyz_string = _Molecule_from_xyz_string(xyz)\n", " | Generate molecule from string in XYZ format.\n", " | \n", " | :param xyz:\n", " | String with XYZ structure.\n", " | \n", " | :return:\n", " | The molecule.\n", " | \n", " | read_str = _Molecule_read_str(xyzstr, units='angstrom')\n", " | Reads molecule from xyz string.\n", " | \n", " | :param xyzstr:\n", " | The xyz string.\n", " | :param units:\n", " | The unit of coordinates.\n", " | \n", " | :return:\n", " | The molecule.\n", " | \n", " | read_xyz = _Molecule_read_xyz(xyzfile)\n", " | Reads molecule from file in XYZ format.\n", " | \n", " | :param xyzfile:\n", " | File with molecular structure in XYZ format.\n", " | \n", " | :return:\n", " | The molecule.\n", " | \n", " | ----------------------------------------------------------------------\n", " | Data and other attributes defined here:\n", " | \n", " | __hash__ = None\n", " | \n", " | ----------------------------------------------------------------------\n", " | Static methods inherited from pybind11_builtins.pybind11_object:\n", " | \n", " | __new__(*args, **kwargs) from pybind11_builtins.pybind11_type\n", " | Create and return a new object. See help(type) for accurate signature.\n", "\n" ] }, { "data": { "text/plain": [ "array([[ 0. , 0. , 0. ],\n", " [ 0. , 1.4, 1.1],\n", " [ 0. , -1.4, 1.1]])" ] }, "execution_count": 5, "metadata": {}, "output_type": "execute_result" } ], "source": [ "help(vlx.Molecule)\n", "\n", "h2o.get_coordinates()" ] }, { "cell_type": "markdown", "id": "announced-manor", "metadata": {}, "source": [ "### Gaussian molecular basis sets\n", "\n", "Once we have a molecule, we can define a Gaussian molecular basis set. `vlx.MolecularBasis.get_avail_basis` shows the list of already available basis sets in VeloxChem for a given chemical element." ] }, { "cell_type": "code", "execution_count": 6, "id": "sixth-pipeline", "metadata": { "tags": [ "hide_output" ] }, "outputs": [ { "data": { "text/plain": [ "['6-31++G',\n", " '6-31+G',\n", " '6-31+G(D,P)',\n", " '6-311++G',\n", " '6-311++G(2D,2P)',\n", " '6-311++G(D,P)',\n", " '6-311+G',\n", " '6-311+G(2D,P)',\n", " '6-311G',\n", " '6-311G(2DF,2PD)',\n", " '6-31G',\n", " '6-31G(2DF,P)',\n", " 'AUG-CC-PCVDZ',\n", " 'AUG-CC-PCVQZ',\n", " 'AUG-CC-PCVTZ',\n", " 'AUG-CC-PVDZ',\n", " 'AUG-CC-PVTZ',\n", " 'CC-PCVDZ',\n", " 'CC-PVDZ',\n", " 'CC-PVDZ-RI',\n", " 'CC-PVTZ',\n", " 'CC-PVTZ-RI',\n", " 'D-AUG-CC-PVDZ',\n", " 'D-AUG-CC-PVQZ',\n", " 'D-AUG-CC-PVTZ',\n", " 'DEF2-QZVP',\n", " 'DEF2-QZVPD',\n", " 'DEF2-QZVPP',\n", " 'DEF2-QZVPPD',\n", " 'DEF2-RI-J',\n", " 'DEF2-SV(P)',\n", " 'DEF2-SVP',\n", " 'DEF2-SVPD',\n", " 'DEF2-TZVP',\n", " 'DEF2-TZVPP',\n", " 'DEF2-TZVPPD',\n", " 'MIN-CC-PVDZ',\n", " 'SADLEJ-PVTZ',\n", " 'STO-3G',\n", " 'STO-6G']" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "vlx.MolecularBasis.get_avail_basis(\"C\")" ] }, { "cell_type": "markdown", "id": "after-mistress", "metadata": {}, "source": [ "We generate a basis object for our molecule using one of the bases included in VeloxChem with:" ] }, { "cell_type": "code", "execution_count": 7, "id": "manufactured-canada", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Molecular Basis (Atomic Basis)\n", "================================\n", "\n", "Basis: STO-3G \n", "\n", " Atom Contracted GTOs Primitive GTOs \n", "\n", " H (1S) (3S) \n", " O (2S,1P) (6S,3P) \n", "\n", "Contracted Basis Functions : 7 \n", "Primitive Basis Functions : 21 \n", "\n", "\n" ] } ], "source": [ "basis = vlx.MolecularBasis.read(h2o, \"sto-3g\")\n", "\n", "print(basis.get_string(h2o))" ] }, { "cell_type": "markdown", "id": "turkish-promise", "metadata": {}, "source": [ "### Self-consistent field calculation drivers\n", "\n", "You launch restricted or unrestricted SCF calculations through their respective driver object: `vlx.ScfRestrictedDriver` and `vlx.ScfUnrestrictedDriver`." ] }, { "cell_type": "code", "execution_count": 8, "id": "wanted-southeast", "metadata": {}, "outputs": [], "source": [ "scfdrv = vlx.ScfRestrictedDriver()" ] }, { "cell_type": "markdown", "id": "confident-twelve", "metadata": {}, "source": [ "You can configure the drivers with options passed as a Python dictionary." ] }, { "cell_type": "code", "execution_count": 9, "id": "fitting-entrepreneur", "metadata": {}, "outputs": [], "source": [ "settings = {\"conv_thresh\": 1e-6}\n", "scfdrv.update_settings(settings)" ] }, { "cell_type": "markdown", "id": "constant-doctor", "metadata": {}, "source": [ "Finally, we run the calculation with the `compute` method. The output will be exactly identical to the one obtained from running the input file." ] }, { "cell_type": "code", "execution_count": 10, "id": "practical-bench", "metadata": { "tags": [ "hide_output" ] }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ " \n", " Self Consistent Field Driver Setup \n", " ==================================== \n", " \n", " Wave Function Model : Spin-Restricted Hartree-Fock \n", " Initial Guess Model : Superposition of Atomic Densities \n", " Convergence Accelerator : Two Level Direct Inversion of Iterative Subspace \n", " Max. Number of Iterations : 50 \n", " Max. Number of Error Vectors : 10 \n", " Convergence Threshold : 1.0e-06 \n", " ERI Screening Scheme : Cauchy Schwarz + Density \n", " ERI Screening Mode : Dynamic \n", " ERI Screening Threshold : 1.0e-12 \n", " Linear Dependence Threshold : 1.0e-06 \n", " \n", "* Info * Nuclear repulsion energy: 9.3436381580 au \n", " \n", "* Info * Overlap matrix computed in 0.18 sec. \n", " \n", "* Info * Kinetic energy matrix computed in 0.00 sec. \n", " \n", "* Info * Nuclear potential matrix computed in 0.00 sec. \n", " \n", "* Info * Orthogonalization matrix computed in 0.09 sec. \n", " \n", "* Info * SAD initial guess computed in 0.00 sec. \n", " \n", "* Info * Starting Reduced Basis SCF calculation... \n", "* Info * ...done. SCF energy in reduced basis set: -74.960337067134 au. Time: 0.48 sec. \n", " \n", "* Info * Overlap matrix computed in 0.09 sec. \n", " \n", "* Info * Kinetic energy matrix computed in 0.01 sec. \n", " \n", "* Info * Nuclear potential matrix computed in 0.00 sec. \n", " \n", "* Info * Orthogonalization matrix computed in 0.00 sec. \n", " \n", " \n", " Iter. | Hartree-Fock Energy | Energy Change | Gradient Norm | Max. Gradient | Density Change \n", " -------------------------------------------------------------------------------------------- \n", " 1 -74.960337068856 0.0000000000 0.00002391 0.00000577 0.00000000 \n", " 2 -74.960337069026 -0.0000000002 0.00000785 0.00000178 0.00001947 \n", " 3 -74.960337069049 -0.0000000000 0.00000050 0.00000012 0.00001040 \n", " \n", " *** SCF converged in 3 iterations. Time: 0.49 sec. \n", " \n", " Spin-Restricted Hartree-Fock: \n", " ----------------------------- \n", " Total Energy : -74.9603370690 au \n", " Electronic Energy : -84.3039752270 au \n", " Nuclear Repulsion Energy : 9.3436381580 au \n", " ------------------------------------ \n", " Gradient Norm : 0.0000004994 au \n", " \n", " \n", " Ground State Information \n", " ------------------------ \n", " Charge of Molecule : 0.0 \n", " Multiplicity (2S+1) : 1.0 \n", " Magnetic Quantum Number (M_S) : 0.0 \n", " \n", " \n", " Spin Restricted Orbitals \n", " ------------------------ \n", " \n", " Molecular Orbital No. 1: \n", " -------------------------- \n", " Occupation: 2.0 Energy: -20.24119 au \n", " ( 1 O 1s : -0.99) \n", " \n", " Molecular Orbital No. 2: \n", " -------------------------- \n", " Occupation: 2.0 Energy: -1.27737 au \n", " ( 1 O 1s : 0.23) ( 1 O 2s : -0.83) ( 2 H 1s : -0.16) \n", " ( 3 H 1s : -0.16) \n", " \n", " Molecular Orbital No. 3: \n", " -------------------------- \n", " Occupation: 2.0 Energy: -0.62456 au \n", " ( 1 O 1p-1: -0.61) ( 2 H 1s : -0.44) ( 3 H 1s : 0.44) \n", " \n", " Molecular Orbital No. 4: \n", " -------------------------- \n", " Occupation: 2.0 Energy: -0.45684 au \n", " ( 1 O 2s : 0.54) ( 1 O 1p0 : -0.78) ( 2 H 1s : -0.27) \n", " ( 3 H 1s : -0.27) \n", " \n", " Molecular Orbital No. 5: \n", " -------------------------- \n", " Occupation: 2.0 Energy: -0.39296 au \n", " ( 1 O 1p+1: 1.00) \n", " \n", " Molecular Orbital No. 6: \n", " -------------------------- \n", " Occupation: 0.0 Energy: 0.62228 au \n", " ( 1 O 2s : -0.91) ( 1 O 1p0 : -0.75) ( 2 H 1s : 0.81) \n", " ( 3 H 1s : 0.81) \n", " \n", " Molecular Orbital No. 7: \n", " -------------------------- \n", " Occupation: 0.0 Energy: 0.75803 au \n", " ( 1 O 1p-1: -1.00) ( 2 H 1s : 0.85) ( 3 H 1s : -0.85) \n", " \n" ] } ], "source": [ "scfdrv.compute(h2o, basis)" ] }, { "cell_type": "markdown", "id": "dangerous-gauge", "metadata": {}, "source": [ "The driver object contains a wealth of useful information:\n", "\n", "- The final SCF energy: `scfdrv.get_scf_energy()`\n", "- Various tensors, such as the molecular orbital coefficients, AO density and Fock matrices: `scfdrv.scf_tensors`.\n", "\n", "The `scf_tensors` property is a Python dictionary with keys:\n", "\n", "- `C`, for the molecular orbital coefficients.\n", "- `E`, for the orbital energies.\n", "- `S`, for the overlap integral.\n", "- `D`, for the $\\alpha$ and $\\beta$ density matrices as a tuple of NumPy arrays.\n", "- `F`, for the $\\alpha$ and $\\beta$ Fock matrices as a tuple of NumPy arrays." ] }, { "cell_type": "code", "execution_count": 11, "id": "proof-malpractice", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "dict_keys(['C', 'E', 'S', 'D', 'F'])" ] }, "execution_count": 11, "metadata": {}, "output_type": "execute_result" } ], "source": [ "scfdrv.scf_tensors.keys()" ] }, { "cell_type": "code", "execution_count": 12, "id": "mature-stephen", "metadata": { "tags": [ "hide_output" ] }, "outputs": [ { "data": { "text/plain": [ "(array([[ 1.05303213e+00, -2.22302631e-01, -1.46417671e-02,\n", " -1.46417671e-02, -2.04228488e-16, 5.41902851e-02,\n", " -7.42025150e-18],\n", " [-2.22302631e-01, 9.82076521e-01, -1.64989375e-02,\n", " -1.64989375e-02, -2.16731540e-16, -3.10480924e-01,\n", " 6.28860220e-17],\n", " [-1.46417671e-02, -1.64989375e-02, 2.97409605e-01,\n", " -9.77472101e-02, 2.70379139e-01, 2.33638805e-01,\n", " 1.91065463e-16],\n", " [-1.46417671e-02, -1.64989375e-02, -9.77472101e-02,\n", " 2.97409605e-01, -2.70379139e-01, 2.33638805e-01,\n", " -3.43153996e-16],\n", " [-2.04228488e-16, -2.16731540e-16, 2.70379139e-01,\n", " -2.70379139e-01, 3.70004394e-01, -1.39439910e-15,\n", " 1.06935025e-16],\n", " [ 5.41902851e-02, -3.10480924e-01, 2.33638805e-01,\n", " 2.33638805e-01, -1.39439910e-15, 6.22944826e-01,\n", " -7.67046565e-17],\n", " [-7.42025150e-18, 6.28860220e-17, 1.91065463e-16,\n", " -3.43153996e-16, 1.06935025e-16, -7.67046565e-17,\n", " 1.00000000e+00]]),\n", " array([[ 1.05303213e+00, -2.22302631e-01, -1.46417671e-02,\n", " -1.46417671e-02, -2.04228488e-16, 5.41902851e-02,\n", " -7.42025150e-18],\n", " [-2.22302631e-01, 9.82076521e-01, -1.64989375e-02,\n", " -1.64989375e-02, -2.16731540e-16, -3.10480924e-01,\n", " 6.28860220e-17],\n", " [-1.46417671e-02, -1.64989375e-02, 2.97409605e-01,\n", " -9.77472101e-02, 2.70379139e-01, 2.33638805e-01,\n", " 1.91065463e-16],\n", " [-1.46417671e-02, -1.64989375e-02, -9.77472101e-02,\n", " 2.97409605e-01, -2.70379139e-01, 2.33638805e-01,\n", " -3.43153996e-16],\n", " [-2.04228488e-16, -2.16731540e-16, 2.70379139e-01,\n", " -2.70379139e-01, 3.70004394e-01, -1.39439910e-15,\n", " 1.06935025e-16],\n", " [ 5.41902851e-02, -3.10480924e-01, 2.33638805e-01,\n", " 2.33638805e-01, -1.39439910e-15, 6.22944826e-01,\n", " -7.67046565e-17],\n", " [-7.42025150e-18, 6.28860220e-17, 1.91065463e-16,\n", " -3.43153996e-16, 1.06935025e-16, -7.67046565e-17,\n", " 1.00000000e+00]]))" ] }, "execution_count": 12, "metadata": {}, "output_type": "execute_result" } ], "source": [ "scfdrv.scf_tensors[\"D\"]" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.8.9" } }, "nbformat": 4, "nbformat_minor": 5 }