.. : PyBEST: Pythonic Black-box Electronic Structure Tool : Copyright (C) 2016-- The PyBEST Development Team : : This file is part of PyBEST. : : PyBEST is free software; you can redistribute it and/or : modify it under the terms of the GNU General Public License : as published by the Free Software Foundation; either version 3 : of the License, or (at your option) any later version. : : PyBEST is distributed in the hope that it will be useful, : but WITHOUT ANY WARRANTY; without even the implied warranty of : MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the : GNU General Public License for more details. : : You should have received a copy of the GNU General Public License : along with this program; if not, see : : -- .. _user_estruct_rcc: The coupled cluster module ########################## Before reading this section, please read the general introduction mentioned in :ref:`user_posthf_intro`. This part of the Documentation builds upon it. With PyBEST, you can perform a coupled cluster correction on top of an RHF (see :ref:`user_hf`) reference function. Both modules are explained in greater detail below. .. _getstartedcc: Quick Guide: RCC ================ This is a short introduction to the RCC module. More information on the input and output structure can be found in the following sections. Similar to the previous modules, we will assume the following names for all PyBEST objects :lf: A :py:class:`~pybest.linalg.base.LinalgFactory` instance (see :ref:`user_linalg_intro`). :occ_model: An Aufbau occupation model of the :py:class:`~pybest.scf.occ.AufbauOccModel` class :kin: the kinetic energy integrals :ne: the nucleus-electron attraction integrals :eri: the two-electron repulsion integrals .. _qg_rcc: RHF reference function ---------------------- We assume that you have performed a restricted Hartree-Fock calculation, whose results are stored in the :py:class:`~pybest.io.iodata.IOData` container ``hf_`` (see :ref:`user_hf`). The code snippet below shows how to perform an RCCSD calculation with an RHF reference function in PyBEST, .. code-block:: python ccsd = RCCSD(lf, occ_model) ccsd_ = ccsd(kin, ne, eri, hf_) The results are returned as an :py:class:`~pybest.io.iodata.IOData` container, while all results are saved to the ``pybest-results/checkpoint_RCCSD.h5`` file. Specifically, the :py:class:`~pybest.io.iodata.IOData` container contains (amongst others) the following attributes :e_tot: The total RHFLCCD energy (including all external terms) :e_corr: The RHFLCCD correlation energy :e_ref: The energy of the reference determinant (here, the Hartree-Fock determinant) :t_1: The single excitation cluster amplitudes :t_2: The double excitation cluster amplitudes .. note:: The order of arguments does not matter because PyBEST assigns them in an automatic fashion. To exclude single excitations, you have to use the :py:class:`~pybest.cc.rccsd.RCCD` class. The overall input structure is equivalent to :py:class:`~pybest.cc.rccsd.RCCSD`, .. code-block:: python ccd = RCCD(lf, occ_model) ccd_ = ccd(kin, ne, eri, hf_) All results are saved to the ``pybest-results/checkpoint_RCCD.h5`` checkpoint file. Defining a frozen core ---------------------- By default, all (occupied and virtual) orbitals are active. To freeze some (occupied) orbitals, the number of frozen cores has to be specified when an instance of some RCC class is created. For instance, if you wish to freeze the first 4 (occupied) orbitals in a RHF-CCSD calculation, specify the ``ncore`` argument during the initialization of :py:class:`~pybest.cc.rccsd.RCCSD` .. code-block:: python # Select 4 frozen core orbital #----------------------------- ccsd = RCCSD(lf, occ_model, ncore=4) ccsd_ = ccsd(kin, ne, eri, hf_) This syntax is working for all CC modules mentioned above. .. _cc_restart: Restart options --------------- The :py:class:`~pybest.cc.rccsd.RCC` module supports restarts from PyBEST's internal checkpoint files. If your :py:class:`~pybest.cc.rccsd.RCC` calculation did not converge and you wish to restart it and take the amplitudes from a previous calculation as guess, specify the path for checkpoint file using ``initguess`` keyword. .. code-block:: python ccd = RCCD(lf, occ_model) # Restart only amplitudes from checkpoint_RCCD_d.h5 ccd_ = ccd( kin, ne, eri, hf_, initguess='pybest-results/checkpoint_RCCD.h5' ) PyBEST will read the amplitudes stored in ``pybest-results/checkpoint_RCCD.h5`` and provide those as guess amplitudes to the RCC solver. .. note:: Currently, PyBEST does not dump the Hamiltonian to disk that is used in the CC calculation. Thus, if the user modifies the Hamiltonian, the calculation might be corrupted and special care is advised. .. _rccsd_keywords: Summary of keyword arguments ---------------------------- The :py:class:`~pybest.cc.rccsd.RCC` module supports various keyword arguments that allow us to steer the optimization of the cluster amplitudes. In the following, all supported keyword arguments are listed together with their default values. Please note, that for most cases the default values should be sufficient to reach convergence. :checkpoint: (int) frequency of checkpointing :checkpoint_path: (str) path to file where current solution is dumped (default pybest-results/checkpoint_method.h5). :e_ref: (float) reference energy (default value 0.0). :initguess: (str or numpy.ndarray) initial guess for amplitudes, one of: "random" - random numbers, "mp2" - performs MP2 calculations before, "const" - uses constant number, used for debugging, path to file - path to .h5 file containing amplitudes (saved as checkpoint_method.h5 by default), numpy.ndarray - 1d vector with amplitudes. :maxiter: (int) maximum number of iterations (default: 100). :solver: (str) one of scipy.optimize.root solvers (default: "krylov"). :threshold: (float) tolerance for amplitudes (default: 1e-8). RCC calculation on the water molecule using a frozen core --------------------------------------------------------- This is a basic example on how to perform a CCD and CCSD calculation in PyBEST. This script performs a CCD calculation followed by a CCSD calculation on the water molecule using the cc-pVDZ basis set and freezing the lowest orbital. .. literalinclude:: ../src/pybest/data/examples/rccsd/rccsd_frozen_core_water_cc-pvdz.py :caption: data/examples/rccsd/rccsd_frozen_core_water_cc-pvdz.py :lines: 3- RCC calculation on the water molecule ------------------------------------- This is a basic example on how to perform a CCSD and CCSD calculation in PyBEST. This script performs a CCD calculation followed by a CCSD calculation on the water molecule using the cc-pVDZ basis set. .. literalinclude:: ../src/pybest/data/examples/rccsd/rccsd_water_cc-pvdz.py :caption: data/examples/rccsd/rccsd_water_cc-pvdz.py :lines: 3- RCCSD calculation using data from file as an initial guess for amplitudes ------------------------------------------------------------------------- This is a basic example on how to perform a CCSD calculation in PyBEST. This script performs a CCSD calculation where the maximum number of iterations is too small to converge. This calculation is followed up by the CCSD calculations that use the previous results that were dumped to file. .. literalinclude:: ../src/pybest/data/examples/rccsd/rccsd_restart_water_cc-pvdz.py :caption: data/examples/rccsd/rccsd_frozen_core_water_cc-pvdz.py :lines: 3-