{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Overview of adcc\n", "\n", "ADC-connect – or adcc in short – is a Python-based framework to connect to arbitrary programs and perform calculations based on the algebraic-diagrammatic construction approach (ADC) on top of their existing self-consistent field (SCF) procedures. Four SCF codes can be used with adcc out of the box, namely [molsturm](https://molsturm.org),\n", "[psi4](https://psicode.org/), [PySCF](https://github.com/pyscf/pyscf), and [veloxchem](https://veloxchem.org/). Without expressing any particular preference these notebooks will focus on psi4 and pyscf to obtain an SCF reference, but the other codes work very similar.\n", "\n", "The purpose of these tutorial notebooks is to introduce adcc from the practitioners perspective and give an idea of supported features. For a more structured documentation and an API reference, see https://adc-connect.org. In particular for installation instructions see https://adc-connect.org/installation.html. If you cannot be bothered to install adcc just to try these notebooks, just go to https://try.adc-connect.org.\n", "\n", "These notebooks do not put a strong emphasis on the theory and numerical methods behind adcc, with the exception of a small introduction in [02_Theory.ipynb](02_Theory.ipynb). If have the desire to dig deeper, a review of ADC literature is provided in https://adc-connect.org/theory.html." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## A first calculation" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Running a simple ADC calculation on top of an HF reference obtained with one of the four packages mentioned above takes pretty much a single line of code as we will see now.\n", "\n", "Let's say we want to calculate an excitation spectrum of water in the UV/vis range. First we prepare our reference in psi4 in a cc-pVTZ basis:" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "import psi4\n", "import adcc\n", "\n", "# By default adcc uses all cores on the machnine. For the binder VMs we need\n", "# to reduce this a little and tell adcc to use two threads\n", "adcc.set_n_threads(2)\n", "\n", "# Run SCF in Psi4 using a cc-pVTZ basis\n", "mol = psi4.geometry(\"\"\"\n", " O 0 0 0\n", " H 0 0 1.795239827225189\n", " H 1.693194615993441 0 -0.599043184453037\n", " symmetry c1\n", " units au\n", "\"\"\")\n", "psi4.core.be_quiet()\n", "psi4.set_options({'basis': \"cc-pvtz\", })\n", "scf_e, wfn = psi4.energy('SCF', return_wfn=True)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Next we perform an ADC(2) calculation on top, asking for 10 singlet states:" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "state = adcc.adc2(wfn, n_singlets=10)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We are greeted with a nice convergence table and the result of the calculation is stored in the `state`. Now let's actually look at it:" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "state" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We get a nice table, summarising excitation energies, oscillator strengths as well as `|v1|^2` and `|v2|^2`, the norms of the singles and doubles blocks of the ADC vectors (to be explained in more detail in [02_Theory.ipynb](02_Theory.ipynb)).\n", "\n", "Note: If you are running adcc from a script, than just `print(state)` will not give you the same answer, instead you will have to explicitly call the `state.describe()` function, like so:" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "print(state.describe())" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "With 10 states nicely computed, we want to plot a simulated excitation spectrum, which again can be done in a single function call:" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "state.plot_spectrum()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The plotted spectrum shows both the exact excitation energies with dashed lines as well as a Lorenzian envelope in solid. Notice that plot spectrum is pretty flexible and takes for example typical arguments you can pass to a `plt.plot` function call from matplotlib, for example:" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "from matplotlib import pyplot as plt\n", "\n", "state.plot_spectrum(color=\"green\", label=\"ADC(2)\")\n", "plt.legend()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Another important way to get insight into a result is by looking at the dominating amplitudes (dominating orbital excitations) for each of the states:" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "print(state.describe_amplitudes())" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The default output is pretty detailed. You can get less printed by increasing the tolerance for printing an amplitude. For example:" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "print(state.describe_amplitudes(0.1))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Comparing ADC methods" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Let us get a feeling for adcc as well as the different ADC methods implemented by comparing the excitation spectra for our water molecule in a single plot. Since visual comparison can certainly not distinguish five digits in the energies we will also lower the convergence tolerance a little to speed up the calculation:" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "n_singlets = 10\n", "state_1 = adcc.adc1(wfn, n_singlets=n_singlets, conv_tol=1e-3) # ADC(1)\n", "print()\n", "state_2x = adcc.adc2x(state.ground_state, n_singlets=n_singlets, # ADC(2)-x, as guesses\n", " guesses=state.excitation_vectors, conv_tol=1e-3) # use ADC(2) result\n", "print()\n", "state_3 = adcc.adc3(state_2x.ground_state, n_singlets=n_singlets,\n", " guesses=state_2x.excitation_vectors, conv_tol=1e-3) # ADC(3)\n", "\n", "# Plot all results in one plot\n", "state_1.plot_spectrum(label=\"ADC(1)\")\n", "state.plot_spectrum(label=\"ADC(2)\")\n", "state_2x.plot_spectrum(label=\"ADC(2)-x\")\n", "state_3.plot_spectrum(label=\"ADC(3)\")\n", "plt.legend()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Getting help\n", "\n", "Documentation in adcc is still a little sparser than it should. Still most functions are documented and their docstrings can be easily obtained, for example:" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "adcc.adc2?" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Last but not least, if you find a bug or need help you are always welcome to file [an issue](https://github.com/adc-connect/adcc/issues) or contact us by mail: developers@adc-connect.org." ] } ], "metadata": { "@webio": { "lastCommId": null, "lastKernelId": null }, "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.2" } }, "nbformat": 4, "nbformat_minor": 4 }