.. : 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_tcc: The restricted tailored 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 single-reference tailored coupled-cluster correction on top of a CASCI-type wave function. The module is explained in greater detail below. If you use this module, please cite [leszczyk2022]_. .. _getstartedtcc: Quick Guide: RtCC ================= This is a short introduction to the RtCC 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 electron repulsion integrals .. _qg_tcc: :math:`e^{T_{CAS}}` reference function in the canonical RHF orbital basis ------------------------------------------------------------------------- Currently, CASCI-type calculations are not available in PyBEST and must be performed with external software, e.g., the QC-DMRG-Budapest program. To ensure the correctness of calculations, the one- and two-body Hamiltonian integrals must be stored in the same order at all stages of the calculations. Thus, we suggest either generating integrals in the FCIDUMP format using some external software or generating integrals using PyBEST and read them in the external software. The following snippet allows us to save one- and two-body integrals to file. The integrals from the cas.FCIDUMP file can be used for CASCI-type calculations with some external software. .. code-block:: python from pybest.orbital_utils import transform_integrals, split_core_active # Transform Hamiltonian to MO basis mo_ints = transform_integrals(kin, ne, eri, hf_output.orb_a) one = mo_ints.one[0] two = mo_ints.two[0] # Write all integrals to a FCIDUMP file data = IOData(one=one, two=two, e_core=e_core, ms2=0, nelec=14) data.to_file("all.FCIDUMP") # Get CAS Hamiltonian mo_ints_cas = split_core_active(one, two, ncore=4, nactive=6, e_core=e_core) one_cas = mo_ints_cas.one two_cas = mo_ints_cas.two # Write CAS integrals to a FCIDUMP file data = IOData(one=one_cas, two=two_cas, e_core=e_core, ms2=0, nelec=6) data.to_file("cas.FCIDUMP") See :ref:`hamiltonian_io_fcidump` for more details about dumping integrals. We assume that you have performed some CASCI calculations using ``cas.FCIDUMP`` with some external software, and you obtained a file with the CC amplitudes reconstructed from CASCI result called ``T_DUMP``. One example of such a file (obtained with the QC-DMRG-Budapest program) is available in ``data/examples/dmrg-tccsd/data/T_DUMP``. The following code snippet below shows how to perform an RtCCSD calculation with a DMRG reference function in PyBEST. .. code-block:: python # Read integrals data = IOData.from_file("all.FCIDUMP") # Do tailored coupled cluster calculations tcc = RtCCSD(lf, occ_model) tcc_out = tcc( data.one, data.two, data.orb_a, external_file="T_DUMP", e_ref=energy_of_reference_determinant, ) .. note:: The order of arguments does not matter because PyBEST recognizes them in an automatic fashion. The results are returned as an :py:class:`~pybest.io.iodata.IOData` container, while all results are saved to the ``pybest-results/checkpoint_RtCCSD.h5`` file. Specifically, the :py:class:`~pybest.io.iodata.IOData` container contains (amongst others) the following attributes :e_tot: the total RtCCSD energy (sum of reference energy and RCC energy) :e_corr: the RtCCSD correlation energy :e_corr_s: the single excitation contribution in RtCCSD correlation energy :e_corr_d: the double excitation contribution in RtCCSD correlation energy :e_ref: the energy of the reference determinant :t_1: the single excitation cluster amplitudes :t_2: the double excitation cluster amplitudes :method: the name of the model (here: RCCSD) :nocc: the number of all occupied orbitals (including core orbitals) :nvirt: the number of all virtual (unoccupied in RHF model) orbitals :ncore: the number of (frozen) core orbitals (occupied but not correlated) :nact: the number of correlated occupied and unoccupied orbitals :converged: boolean value indicating if CC converged :math:`e^{T_{CAS}}` reference function in ROOpCCD orbital basis --------------------------------------------------------------- To perform analogous calculations in the ROOpCCD basis, follow the steps from :ref:`qg_tcc` but use the ROOpCCD output instead of the RHF output while transforming the integrals. .. _tcc_keywords: Summary of keyword arguments ---------------------------- The :py:class:`~pybest.cc.rtcc.RtCCSD` and :py:class:`~pybest.cc.rtcc.RtCCSD` classes support 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. :filename: (str) path to file where current solution is dumped (default pybest-results/checkpoint_method.h5). :e_ref: (float) reference energy (default value 0.0). :external_file: (str) path to file containing external amplitudes. :initguess: (str or numpy.ndarray) initial guess for the amplitudes, one of: "random", "mp2", "const", path to .h5 file, or IOData instance. :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). RtCCSD calculation on the carbon nitride cation ----------------------------------------------- This is a basic example on how to perform a RtCCSD calculation in PyBEST. These scripts perform a RtCCSD calculation on the carbon nitride action using the CC amplitudes from DMRG calculations. In the 0th step, we generate one- and two-body integrals in ROOpCCD orbital basis. .. literalinclude:: ../src/pybest/data/examples/dmrg-tccsd/step0.py :caption: data/examples/dmrg-tccsd/step0.py :lines: 3- The generated complete active space integrals (saved as ``cas.FCIDUMP``) can be utilized to compute a matrix product state (MPS) wave function. We extract CC amplitudes from the MPS and save them as ``data/T_DUMP``. In the second step, we perform DMRG-tCCSD calculations based on the data obtained from the previous steps. .. literalinclude:: ../src/pybest/data/examples/dmrg-tccsd/step2.py :caption: data/examples/dmrg-tccsd/step2.py :lines: 3-