.. : 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_pccd: The pCCD module ############### Two-electron functions, called geminals, can be used to incorporate electron correlation effects into the many-particle wave function. PyBEST supports a special type of geminal-based wave function models, the antisymmetric product of 1-reference orbital geminals (pCCD), which is equivalent to pair-coupled cluster doubles. The pCCD wave function effectively parameterizes the doubly occupied configuration interaction wave function (DOCI), but requires only mean-field computational cost in contrast to the factorial scaling of traditional DOCI implementations [limacher2013]_. Currently, the pCCD module is limited to closed-shell systems only. .. _user_pccd_model: The pCCD model ================ The pCCD wave function ansatz can be rewritten in terms of one-particle functions as a fully general pair-coupled-cluster wave function, .. math:: :label: pccd \vert \textrm{pCCD}\rangle = \exp(\sum_{ia} c_i^a a^\dagger_a a^\dagger_{\bar{a}} a_{\bar{i}} a_i) \vert \Phi_0 \rangle where :math:`a_p^{\dagger}`, :math:`a_{\bar{p}}^{\dagger}`, and :math:`a_p`, :math:`a_{\bar{p}}` are the electron creation and annihilation operators and :math:`p` and :math:`\bar{p}` denote :math:`\alpha` and :math:`\beta` spins, respectively. :math:`\vert \Phi_0 \rangle` is some independent-particle wave function (for instance, the Hartree-Fock determinant). Indices :math:`i` and :math:`a` correspond to virtual and occupied orbitals with respect to :math:`\vert \Phi_0 \rangle`, :math:`P` and :math:`K` denote the number of electron pairs (:math:`P = N/2` with :math:`N` being the total number of electrons) and orbitals, respectively. The geminal coefficient matrix (:math:`\mathbf{C}`) of pCCD links the geminals with the underlying one-particle basis functions and has the following form, .. math:: :label: cia \mathbf{C} = \begin{pmatrix} 1 & 0 & \cdots & 0 & c_{1;P+1} & c_{1;P+2}&\cdots &c_{1;K}\\ 0 & 1 & \cdots & 0 & c_{2;P+1} & c_{2;P+2}&\cdots &c_{2;K}\\ \vdots & \vdots & \ddots & \vdots & \vdots & \vdots &\ddots &\vdots\\ 0 & 0 & \cdots & 1 & c_{P;P+1} & c_{P;P+2}&\cdots & c_{P;K} \end{pmatrix} The exponential form of eq. :eq:`pccd` assures the size extensivity of the geminal wave function, however, to ensure the size consistency, one has to optimize the orbitals (see [boguslawski2014a]_ and [boguslawski2014b]_). The simplest and most robust way is to use the variational orbital optimization (vOO-pCCD) method implemented in PyBEST (see :ref:`oopccd`). Currently supported features ============================ If not mentioned otherwise, the pCCD module supports spin-restricted orbitals and the ``DenseLinalgFactory`` and ``CholeskyLinalgFactory``. Specifically, the following features are provided: 1. Optimization of pCCD (eq. :eq:`pccd`) with (see :ref:`oopccd` and [boguslawski2014a]_) and without orbital optimization (see :ref:`hfpccd`) for a given Hamiltonian (see :ref:`preamble`). 2. Variational orbital optimization and PS2c orbital optimization (see :ref:`keywords-oopccd` to choose the orbital optimizer and [boguslawski2014b]_) 3. Calculation of response one- and two-particle reduced density matrices (see :ref:`responsedms`) 4. Determination of pCCD natural orbitals and occupation numbers (see :ref:`responsedms` and :ref:`natorb`) 5. Calculation of the exact orbital Hessian (see :ref:`exacthessian`). Note that the orbital optimizer uses only a diagonal approximation to the orbital Hessian. Since the exact orbital Hessian is internally evaluated for the ``DenseLinalgFactory``, PyBEST might require large amounts of memory to store the matrix representation of the orbital Hessian. Thus, it is recommended to evaluate the orbital Hessian only for small molecules or basis sets. The pCCD wave function and its response density matrices can then be used for post-processing. This version of PyBEST offers the following post-pCCD features: 1. A posteriori addition of dynamic electron correlation using perturbation theory (see :ref:`pt2ap1rog` and [boguslawski2017a]_) or coupled cluster corrections (see :ref:`user_estruct_lcc` and [boguslawski2015b]_) 2. Analysis of orbital entanglement and correlations in the pCCD wave function using the orbital entanglement module (see :ref:`orbital_entanglementseniorityzero`) .. _getstartedpccd: Quick Guide: pCCD and orbital-optimized pCCD ============================================ If you use this part of the code, please cite [boguslawski2014a]_ and [boguslawski2014b]_. 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_output`` (see :ref:`user_hf`). Furthermore, 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 How to: pCCD ------------ The code snippet below shows how to perform a restricted pCCD calculation in PyBEST, .. code-block:: python pccd = RpCCD(lf, occ_model) pccd_output = pccd(kin, ne, eri, hf_output) The results are returned as an :py:class:`~pybest.io.iodata.IOData` container, while all results are saved to the ``pybest-results/checkpoint_pccd.h5`` file. Specifically, the :py:class:`~pybest.io.iodata.IOData` container contains the following attributes :e_tot: The total pCCD energy (including all external terms) :e_corr: The pCCD correlation energy :e_ref: The energy of the reference determinant (here, the Hartree-Fock determinant) :t_p: The pair amplitudes The pCCD module features additional options, which are discussed in greater detail below. How to: orbital-optimized pCCD ------------------------------ The code snippet below shows how to perform an orbital-optimized pCCD calculation in PyBEST, .. code-block:: python pccd = ROOpCCD(lf, occ_model) pccd_output = pccd(kin, ne, eri, hf_output) The results are returned as an :py:class:`~pybest.io.iodata.IOData` container, while all results are saved to the ``pybest-resutls/checkpoint_pccd.h5`` file. Similar to above, Similar to pCCD (without orbital optimization), the :py:class:`~pybest.io.iodata.IOData` container contains the following attributes :e_tot: The total pCCD energy (including all external terms) :e_corr: The pCCD correlation energy :e_ref: The energy of the reference determinant (here, the Hartree-Fock determinant) :t_p: The pair amplitudes If the orbitals are optimized, the :py:class:`~pybest.io.iodata.IOData` container also includes the :math:`\lambda_p` (electron-pair de-excitation) amplitudes, the optimized pCCD natural orbitals, and the response 1- and 2-particle reduced density matrices (RDMs). :orb_a: The pCCD natural orbitals :l_p: The pair de-excitation amplitudes :dm_1: The response 1-RDM :dm_2: The response 2-RDM The OOpCCD module features additional options, which are discussed in greater detail below. .. _pccd: Detailed input structure of pCCD ================================ .. _preamble: Getting started --------------- To optimize a pCCD wave function, the module requires some one- and two-electron integrals (that is, some Hamiltonian) and an initial guess for the orbitals as input arguments. PyBEST provides different options for specifying the Hamiltonian and an orbital guess. .. note:: See :ref:`user_hamiltonian` for how to define the molecular geometry, basis set, or the Hamiltonian. - The Hamiltonian is divided into three contributions: the one- and two-electron integrals as well as an external term (also referred to as core energy). Possible choices are (see below for examples): 1. In-house calculation of the quantum chemical Hamiltonian expressed in the AO basis (all terms are calculated separately in PyBEST; see :ref:`user_molecularham_matrix_elements` for documentation). 2. In-house calculation of model Hamiltonians. Supported model Hamiltonians are summarized in :ref:`modphysham`. 3. External (one- and two-electron) integrals (in an orthonormal basis) and core energy can be read from a file. The integral file must use the Molpro file format (see :ref:`hamiltonian_io` for more details). - A set of initial guess orbitals can be either generated in PyBEST (including the AO overlap matrix) or read from disk. Examples for initial guess orbitals are: 1. Restricted canonical Hartree-Fock orbitals (see :ref:`user_hf`) 2. Localized orbitals (see :ref:`localization` for documentation) .. _oopccd: pCCD with orbital optimization -------------------------------- If you use this part of the code, please cite [boguslawski2014a]_ and [boguslawski2014b]_. .. _setup-oo-pccd: How to set up a calculation ^^^^^^^^^^^^^^^^^^^^^^^^^^^ After having specified a Hamiltonian and initial guess orbitals, a pCCD calculation with variational orbital optimization can be initiated as follows .. literalinclude:: ../src/pybest/data/examples/pccd/water_oopccd_cc-pvdz.py :lines: 45-49 All arguments above have been introduced in :ref:`getstartedpccd`. .. note:: The core energy (also referred to as external energy) is automatically determined from the :py:class:`~pybest.io.iodata.IOData` container ``hf_output`` and hence does not need to be passed explicitly. The final results of the calculation are returned as an :py:class:`~pybest.io.iodata.IOData` container (here ``oopccd_output``) .. note:: The number of electron pairs is automatically determined using the :py:class:`~pybest.scf.occ.OccModel` instance ``occ_model``. A set of keyword arguments is optional and contain optimization-specific options (like the number of orbital optimization steps, etc.). Their default values are chosen to give a reasonable performance and can be adjusted if convergence difficulties are encountered (see also :ref:`keywords-oopccd`). After the pCCD calculation is finished (because pCCD converged or the maximum number of iterations was reached), all final results are, by default, stored in a checkpoint file ``checkpoint_pccd.h5`` of the ``pybest-results`` directory and can be used for a subsequent restart. .. _oopccd_core: Defining a frozen core ^^^^^^^^^^^^^^^^^^^^^^ The pCCD module in PyBEST supports frozen core orbitals. These orbitals are not optimized during the orbital optimization, neither are they included in the amplitude equations. The frozen (core) orbitals are thus by construction doubly occupied. The number of frozen orbitals can be defined when creating an instance of the :py:class:`~pybest.occ_model.OccupationModel` class, .. literalinclude:: ../src/pybest/data/examples/pccd/water_oopccd_frozen_core_cc-pvdz.py :lines: 26,44-49 .. note:: If ``ncore`` is set to some nonzero number, the orbitals are frozen with respect to their order in the ``orb_a`` attribute passed to :py:class:`~pybest.geminal.rpccd.ROOpCCD`. In general, these are the energetically lowest-lying orbitals or those with the largest occupation numbers. The number of virtual orbitals cannot be changed and hence PyBEST assumes all virtual orbitals to be active. .. _restart-pccd: How to restart ^^^^^^^^^^^^^^ PyBEST checkpoint files +++++++++++++++++++++++ To restart a pCCD calculation (for instance, using the orbitals from a different molecular geometry as initial guess or from a previous calculation using the same molecular geometry), you can use the restart keyword in the function call providing the filename (or its full path) .. literalinclude:: ../src/pybest/data/examples/pccd/water_oopccd_cc-pvdz.py :lines: 51-57 ``pybest-results/checkpoint_pccd.h5`` is the default checkpoint file generated py PyBEST if not specified otherwise by the user. Perturbed restart orbitals ++++++++++++++++++++++++++ Sometimes, it might be advantageous to slightly perturb some of the orbitals used for restarts (for instance, to overcome local minima or to facilitate orbital localization). This can be achieved similar to the SCF module as explained in :ref:`user_scf_restart_perturbed`. .. literalinclude:: ../src/pybest/data/examples/pccd/water_oopccd_cc-pvdz.py :lines: 16,59- After some restart file ``pybest-results/checkpoint_pccd.h5`` has been read in using the :py:meth:`~pybest.io.iodata.IOData.from_file` method of the :py:class:`~pybest.io.iodata.IOData` container class, you can perturb the orbitals using the swapping feature (see :ref:`orbitals_swapping`) or 2x2 rotations (see :ref:`orbitals_rotations`). The modified/updated/perturbed :py:class:`~pybest.io.iodata.IOData` container can then be simply passed as an argument to the :py:class:`~pybest.geminal.rpccd.ROOpCCD` function call (as done above with the :py:class:`~pybest.wrapper.hf.RHF` container). .. _responsedms: Response density matrices ^^^^^^^^^^^^^^^^^^^^^^^^^ PyBEST supports the calculation of the response 1- and 2-particle reduced density matrices (1-RDM and 2-RDM), :math:`\gamma_{pq}` and :math:`\Gamma_{pqrs}`, respectively. Since pCCD is a product of natural geminals, the 1-RDM is diagonal and is calculated from .. math:: \gamma_p = \langle \Psi_0| (1+\hat{\Lambda}) a^\dagger_p a_p | \textrm{pCCD} \rangle, where :math:`\hat{\Lambda}` contains the de-excitation operator, .. math:: \hat{\Lambda} = \sum_{ia} \lambda_i^a (a^\dagger_i a^\dagger_{\bar{i}} a_{\bar{a}} a_a - c_i^a). The most recent response 1-RDM (a ``OneIndex`` instance) of OO-pCCD is saved in the return value of pCCD (that is the :py:class:`~pybest.io.iodata.IOData` container). The 1-RDM is stored as the ``dm_1`` attribute of the :py:class:`~pybest.io.iodata.IOData` container and can be accessed as follows .. code-block:: python one_dm = oopccd_output.dm_1 The response 2-RDM is defined as .. math:: \Gamma_{pqrs} = \langle \Psi_0| (1+\hat{\Lambda})a^\dagger_p a^\dagger_{q} a_{s} a_r| \textrm{pCCD} \rangle. In PyBEST, only the non-zero elements of the response 2-RDM are calculated, which are :math:`\Gamma_{pqpq}=\Gamma_{p\bar{q}p\bar{q}}` and :math:`\Gamma_{p\bar{p}q\bar{q}}`. Specifically, the non-zero elements :math:`\Gamma_{pqpq}` and :math:`\Gamma_{ppqq}` (where we have omitted the information about electron spin) are calculated separately and stored as ``TwoIndex`` objects. Note that :math:`\gamma_p=\Gamma_{p\bar{p}p\bar{p}}`. Similar to the 1-RDM, the most recent 2-RDM is stored in the :py:class:`~pybest.io.iodata.IOData` container as the attribute ``dm_2``. Since most of elements of the 2-RDM are zero, only its non-zero blocks are calculated and stored in PyBEST. To distinguish between the :math:`\Gamma_{pqpq}=\Gamma_{p\bar{q}p\bar{q}}` and :math:`\Gamma_{p\bar{p}q\bar{q}}` block, the 2-RDM is a dictionary with attributes ``ppqq`` and ``pqpq``. They can be accessed as follows .. code-block:: python two_dm = oopccd_output.dm_2 # Access the ppqq block two_dm_ppqq = two_dm['ppqq'] # Access the pqpq block two_dm_pqpq = two_dm['pqpq'] .. note:: In PyBEST, we have :math:`\Gamma_{p\bar{p}q\bar{q}} = 0 \, \forall \, p=q \in \textrm{occupied}` and :math:`\Gamma_{p\bar{q}p\bar{q}} = 0 \, \forall \, p=q \in \textrm{virtual}`. .. _natorb: Natural orbitals and occupation numbers ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ If pCCD converges, the final orbitals are the pCCD natural orbitals and are stored in the ``orb_a`` attribute of the :py:class:`~pybest.io.iodata.IOData` container. The natural orbitals can be exported to the molden file format (see :ref:`user_orbitals_molden`) and visualized using, for instance, `Jmol `_ or `VESTA `_. The natural occupation numbers, the eigenvalues of the response 1-RDM (see :ref:`responsedms`) are stored in the ``occupations`` attribute (a 1-dim np.array) of ``orb_a`` and can be directly accessed after a pCCD calculation using .. code-block:: python oopccd_output.orb_a.occupations .. _exacthessian: The exact orbital Hessian ^^^^^^^^^^^^^^^^^^^^^^^^^ Although the orbital optimizer uses a diagonal approximation to the exact orbital Hessian, the exact orbital Hessian can be evaluated after an pCCD calculation. Note that this feature uses :py:class:`~pybest.linalg.dense.DenseLinalgFactory`. If :py:class:`~pybest.linalg.cholesky.CholeskyLinalgFactory` is chosen, the corresponding 2-electron integrals are transformed back to a dense object prior the calculation of the exact orbital Hessian. Thus, this feature is limited by the memory bottleneck of the 4-index transformation of :py:class:`~pybest.linalg.dense.DenseLinalgFactory` (see also ``indextrans`` in :ref:`keywords-oopccd`). To calculate the exact orbital Hessian, all one-electron and two-electron integrals (``two``) need to be transformed into the pCCD MO basis first (see also :ref:`hamiltonian_transformation`), .. code-block:: python # transform integrals for restricted orbitals oopccd_output.orb_a t_ints = transform_integrals(kin, ne, er, oopccd_output) # transformed one-electron integrals: attribute 'one' (list) (one,) = t_ints.one # or: one = t_ints.one[0] # transformed two-electron integrals: attribute 'two' (list) (two,) = t_ints.two # or: two = t_ints.two[0] where ``kin`` and ``na`` (``er``) are the one- (two-)electron integrals expressed in the AO basis. In the above example, the molecular orbitals (here, the optimized pCCD orbitals) are passed as an :py:class:`~pybest.io.iodata.IOData` container ``oopccd_output``. PyBEST automatically assigns all arguments. Thus their order does not matter. This step can be skipped if the one- and two-electron integrals are already expressed in the (optimized) MO basis. The transformed one- and two-electron integrals (first element in each list) are passed as function arguments to the :py:meth:`~pybest.geminal.rpccd.get_exact_hessian` method of :py:class:`~pybest.geminal.rpccd.ROOpCCD`, which returns a 2-dim np.array with elements :math:`H_{pq,rs} = H_{p,q,r,s}`, .. code-block:: python # calculate the exact orbital hessian of OOpCCD hessian = oopccd.get_exact_hessian(one, two) where ``oopccd`` is an instance of :py:class:`~pybest.geminal.rpccd.ROOpCCD` (see above :ref:`oopccd`). The exact orbital Hessian can be diagonalized using, for instance, the `np.linalg.eigvalsh routine `_, .. code-block:: python # diagonalize the exact orbital Hessian eigv = np.linalg.eigvalsh(hessian) If frozen core orbitals have been defined for the pCCD optimization, the one- and two-electron integrals need to be transformed using the :py:func:`~pybest.orbital_utils.split_core_active` function (see also :ref:`hamiltonian_activespace` for more details) .. code-block:: python # Define an active space with ncore=1 frozen core orbitals and all virtual # active orbitals from a pCCD wave function stored in oopccd_ cas = split_core_active(kin, ne, er, oopccd_, ncore=1) # all one-electron integrals one = cas.one # all two-electron integrals two = cas.two # the core energy (not required for the orbital Hessian) e_core = cas.e_core The exact orbital Hessian can then be calculated and diagonalized using the transformed one- and two-electron integrals. .. _keywords-oopccd: Summary of keyword arguments ^^^^^^^^^^^^^^^^^^^^^^^^^^^^ The :py:class:`~pybest.geminal.rpccd.ROOpCCD` module supports various keyword arguments that allow us to steer the optimization of the pCCD wave function, including orbital optimization. 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. For convergence difficulties see :ref:`pccd_troubleshooting`. :restart: (str) if provided, a restart from a given checkpoint file is performed (see :ref:`restart-pccd`) :indextrans: (str) 4-index transformation. The choice between ``cupy``, ``tensordot`` (default) and ``einsum``. ``tensordot`` is faster than ``einsum``. If ``DenseLinalgFactory`` is used, the memory requirement scales roughly as :math:`3N^4`. Due to the storage of the two-electron integrals, the total amount of memory increases to :math:`4N^4` .. note:: If Cupy is not available or unsuccessful, ``td`` is selected instead. :warning: (boolean) if ``True``, (scipy) solver-specific warnings are printed (default ``False``) :guess: (dictionary) initial guess specifications: :type: (str) guess type. One of ``mp2`` (default), ``random`` (random numbers), ``const`` (``1.0`` scaled by **factor**) :factor: (float) a scaling factor for the initial guess of type ``type`` (default ``-0.1``, applies only to ``random`` and ``const``) :geminal: (1-dim np.array) external guess for geminal coefficients (default ``None``). If provided, **type** and **factor** are ignored. The elements of the geminal matrix of eq. :eq:`cia` have to be indexed in C-like order. Note that the identity block is not required. The size of the 1-dim np.array is thus equal to the number of unknowns, that is, :math:`n_{\rm pairs}*n_{\rm virtuals}`. :lagrange: (1-dim np.array) external guess for Lagrange multipliers (default ``None``). If provided, **type** and **factor** are ignored. The elements have to be indexed in C-like order. The size of the 1-dim np.array is equal to the number of unknowns, that is, :math:`n_{\rm pairs}*n_{\rm virtuals}`. :solver: (dictionary) scipy wave function/Lagrange solver: :wfn: (str) wave function solver (default ``krylov``) :lagrange: (str) Lagrange multiplier solver (default ``krylov``) .. note:: The exact Jacobian of **wfn** and **lagrange** is not supported. Thus, scipy solvers that need the exact Jacobian cannot be used. See `scipy root-solvers `_ for more details. :maxiter: (dictionary) a maximum number of iterations: :wfniter: (int) maximum number of iterations for the **wfn/lagrange** solver (default ``200``) :orbiter: (int) maximum number of orbital optimization steps (default ``100``) :thresh: (dictionary) optimization thresholds: :wfn: (float) optimization threshold for geminal coefficients and Lagrange multipliers (default ``1e-12``) :energy: (float) convergence threshold for energy (default ``1e-8``) :gradientnorm: (float) convergence threshold for norm of orbital gradient (default ``1e-4``) :gradientmax: (float) threshold for maximum absolute value of the orbital gradient (default ``5e-5``) :printoptions: (dictionary) print level: :geminal: (boolean) if True, geminal matrix is printed (default ``False``). Note that the identity block is omitted. :ci: (float) threshold for CI coefficients (requires evaluation of a permanent). All coefficients (for a given excitation order) larger than **ci** are printed (default ``0.01``) :excitationlevel: (int) a number of excited pairs w.r.t. the reference determinant for which the wave function amplitudes are reconstructed (default ``1``). At most, the coefficients corresponding to hextuply excited Slater determinants w.r.t the reference determinant can be calculated. .. note:: The reconstruction of the wave function amplitudes requires evaluating a permanent which is in general slow (in addition to the factorial number of determinants in the active space). :dumpci: (dictionary) dump Slater determinants and corresponding CI coefficients to file: :amplitudestofile: (boolean) write wave function amplitudes to file (default ``False``) :amplitudesfilename: (str) filename (default ``pccd_amplitudes.dat``) :stepsearch: (dictionary) optimizes an orbital rotation step: :method: (str) step search method used. One of ``trust-region`` (default), ``None``, ``backtracking`` :optimizer: (str) optimizes step to a boundary of trust radius in ``trust-region``. One of ``pcg`` (preconditioned conjugate gradient), ``dogleg`` (Powell's single dogleg step), ``ddl`` (Powell's double-dogleg step) (default ``ddl``) :alpha: (float) scaling factor for the Newton step. Used in ``backtracking`` and ``None`` method (default ``1.00``) :c1: (float) a parameter used in the Armijo condition of ``backtracking`` (default ``1e-4``) :minalpha: (float) minimum step length used in ``backracking`` (default ``1e-6``). If step length falls below **minalpha**, the ``backtracking`` line search is terminated and the most recent step is accepted :maxiterouter: (int) maximum number of iterations to optimize orbital rotation step (default ``10``) :maxiterinner: (int) maximum number of optimization steps in each step search (used only in ``pcg``, default ``500``) :maxeta: (float) upper bound for estimated vs. actual change in ``trust-region`` (default ``0.75``) :mineta: (float) lower bound for estimated vs. actual change in ``trust-region`` (default ``0.25``) :upscale: (float) scaling factor to increase trust radius in ``trust-region`` (default ``2.0``) :downscale: (float) scaling factor to decrease trust radius in ``trust-region`` (default ``0.25``) :trustradius: (float) initial trust radius (default ``0.75``) :maxtrustradius: (float) maximum trust radius (default ``0.75``) :reset: (int) number of iterations after which the trust radius is reset to maxtrustradius (default ``10``) :threshold: (float) trust-region optimization threshold, only used in ``pcg`` (default ``1e-8``) :checkpoint: (int) frequency of checkpointing. If **checkpoint** > 0, a checkpoint file is written every int iteration (default ``1``) :checkpoint_fn: (str) name of checkpoint file (default ``checkpoint_pccd.h5``) :levelshift: (float) a level shift of Hessian (default ``1e-8``). The absolute value of elements of the orbital Hessian smaller than **levelshift** are shifted by **levelshift** :absolute: (boolean), if ``True``, the absolute value of the orbital Hessian is taken (default ``False``) :sort: (boolean), if ``True``, orbitals are sorted according to their natural occupation numbers. This requires re-solving for the wave function after each orbital optimization step. Works only if **orbitaloptimizer** is set to ``variational`` (default ``True``) :orbitaloptimizer: (str) switch between variational orbital optimization (``variational``) and PS2c orbital optimization (``ps2c``) (default ``variational``) .. _hfpccd: pCCD without orbital optimization ----------------------------------- If you use this part of the code, please cite [boguslawski2014a]_ and [boguslawski2014b]_. How to set-up a calculation ^^^^^^^^^^^^^^^^^^^^^^^^^^^ If you only want to optimize the pCCD amplitudes (that is, to skip the orbital optimization step), use the :py:class:`~pybest.geminal.rpccd.RpCCD` class instead of the :py:class:`~pybest.geminal.rpccd.ROOpCCD` class. The remaining syntax is equivalent to the orbital optimized pCCD module, .. code-block:: python pccd = RpCCD(lf, occ_model) pccd_output = pccd(kin, ne, eri, hf_output) The function call again returns an :py:class:`~pybest.io.iodata.IOData` container. .. note:: In the :py:class:`~pybest.geminal.rpccd.RpCCD` module, the :math:`\Lambda` equations are not solved. Thus, no response density matrices are available. If response density matrices are required, the :py:class:`~pybest.geminal.rpccd.ROOpCCD` class has to be used. To skip the orbital optimization procedure, the ``maxiter`` keyword has to be adjusted as follows (see also :ref:`keywords-oopccd`) .. code-block:: python pccd = ROOpCCD(lf, occ_model) pccd_output = pccd(kin, ne, eri, hf_output, maxiter={'orbiter': 0}) This will force PyBEST to solve the :math:`\Lambda` equations and calculate by default the 1- and 2-RDMs. .. _pccd_core: Defining a frozen core ^^^^^^^^^^^^^^^^^^^^^^ Similar to the orbital-optimized version, frozen core orbitals are also supported. These orbitals are not included in the amplitude equations and are thus doubly occupied. The number of frozen orbitals can be defined when creating an instance of the :py:class:`~pybest.occ_model.OccupationModel` class, .. literalinclude:: ../src/pybest/data/examples/pccd/water_pccd_frozen_core_cc-pvdz.py :lines: 27,47- .. note:: Similar to the ``ncore`` feature of :py:class:`~pybest.geminal.rpccd.ROOpCCD`, the orbitals are frozen with respect to their order in the ``orb_a`` attribute passed to :py:class:`~pybest.geminal.rpccd.RpCCD`. In general, these are the energetically lowest-lying orbitals or those with the largest occupation numbers. The number of virtual orbitals cannot be changed and hence PyBEST assumes all virtual orbitals to be active. .. _keywordspccd: Summary of keyword arguments ^^^^^^^^^^^^^^^^^^^^^^^^^^^^ Similar to the :py:class:`~pybest.geminal.rpccd.ROOpCCD` module, :py:class:`~pybest.geminal.rpccd.RpCCD` also supports various keyword arguments that allow us to steer the optimization of the pCCD wave function. In the following, all supported keyword arguments are listed together with their default values. :indextrans: (str) 4-index transformation. The choice between ``cupy``, ``tensordot`` (default) and ``einsum``. ``tensordot`` is faster than ``einsum``. If ``DenseLinalgFactory`` is used, the memory requirement scales roughly as :math:`3N^4`. Due to the storage of the two-electron integrals, the total amount of memory increases to :math:`4N^4` .. note:: If Cupy is not available or unsuccessful, ``td`` is selected instead. :warning: (boolean) if ``True``, (scipy) solver-specific warnings are printed (default ``False``) :guess: (dictionary) initial guess specifications: :type: (str) guess type. One of ``mp2`` (default), ``random`` (random numbers), ``const`` (``1.0`` scaled by **factor**) :factor: (float) a scaling factor for the initial guess of type ``type`` (default ``-0.1``, applies only to ``random`` and ``const``) :geminal: (1-dim np.array) external guess for geminal coefficients (default ``None``). If provided, **type** and **factor** are ignored. The elements of the geminal matrix of eq. :eq:`cia` have to be indexed in C-like order. Note that the identity block is not required. The size of the 1-dim np.array is thus equal to the number of unknowns, that is, :math:`n_{\rm pairs}*n_{\rm virtuals}`. :solver: (dictionary) scipy wave function solver: :wfn: (str) wave function solver (default ``krylov``) .. note:: The exact Jacobian of **wfn** is not supported. Thus, scipy solvers that need the exact Jacobian cannot be used. See `scipy root-solvers `_ for more details. :maxiter: (dictionary) maximum number of iterations: :wfniter: (int) maximum number of iterations for the **wfn** solver (default ``200``) :thresh: (dictionary) optimization thresholds: :wfn: (float) optimization threshold for geminal coefficients and Lagrange multipliers (default ``1e-12``) :printoptions: (dictionary) print level: :geminal: (boolean) if True, geminal matrix is printed (default ``False``). Note that the identity block is omitted. :ci: (float) threshold for CI coefficients (requires evaluation of a permanent). All coefficients (for a given excitation order) larger than **ci** are printed (default ``0.01``) :excitationlevel: (int) number of excited pairs w.r.t. the reference determinant for which the wave function amplitudes are reconstructed (default ``1``). At most, the coefficients corresponding to hextuply excited Slater determinants w.r.t the reference determinant can be calculated. .. note:: The reconstruction of the wave function amplitudes requires evaluating a permanent which is in general slow (in addition to the factorial number of determinants in the active space). :dumpci: (dictionary) dump Slater determinants and corresponding CI coefficients to file: :amplitudestofile: (boolean) write wave function amplitudes to file (default ``False``) :amplitudesfilename: (str) filename (default ``pccd_amplitudes.dat``) .. _pccd_troubleshooting: Troubleshooting in OOpCCD calculations ====================================== - **How to change the number of orbital optimization steps:** To increase the number of iterations in the orbital optimization, adjust the keyword ``maxiter`` (see :ref:`keywords-oopccd`): .. code-block:: python maxiter={'orbiter': int} where ``int`` is the desired number of iterations - **The energy oscillates during orbital optimization:** The occupied-virtual separation breaks down and the reference determinant cannot be optimized. In some cases, fixing the reference determinant might accelerate convergence. However, the final solution might not be reasonable if the optimized geminal coefficient matrix contains elements that are significantly larger than 1.0 in absolute value. To fix the reference determinant, the ``sort`` keyword has to be set to ``False`` (see :ref:`keywords-oopccd`): .. code-block:: python sort=False - **The orbital optimization converges very, very slowly:** Usually, the orbital optimization converges fast around the equilibrium. For stretched distances (in the vicinity of dissociation, etc.) convergence can be very slow, especially if the final solution results in symmetry-broken orbitals. In such cases, the diagonal approximation to the Hessian is not optimal. However, the current version of PyBEST does not support orbital optimization with the exact Hessian nor Hessian updates. You can either try to start with localized guess orbitals, or some perturbed orbital guess (2x2 rotations of, for instance, s and p orbitals), or simply wait until the optimizer converges. - **How to scan a potential energy surface:** To accelerate convergence, restart from adjacent points on the potential energy surface (see :ref:`restart-pccd`). Using Hartree-Fock orbitals as initial guess might result in convergence difficulties and optimization problems of the reference determinant. .. note:: Different guess orbitals might result in different natural orbitals optimized within the OOpCCD module. This is especially true for stretched bond distances, where (partial) localization might occur rather quickly. Please check the final pCCD natural orbitals by, for instance, dumping them to Molden files (see :ref:`user_orbitals_molden`). - **How to perturb the orbitals:** The initial guess orbitals can be perturbed using a sequence of Givens rotations (see also :ref:`orbitals_rotations`). Givens rotations between orbital pairs can be used if, for instance, the orbital optimizer converges to a saddle point. Example Python scripts ====================== Several complete examples can be found in the directory ``data/examples/pccd``. The water molecule ------------------ This is a basic example on how to perform an orbital-optimized pCCD calculation in PyBEST. This script performs an orbital-optimized pCCD calculation on the water molecule using the cc-pVDZ basis set and canonical RHF orbitals as initial orbitals. .. literalinclude:: ../src/pybest/data/examples/pccd/water_oopccd_cc-pvdz.py :caption: data/examples/pccd/water_oopccd_cc-pvdz.py :lines: 3-50 The water molecule including restart for the same geometry ---------------------------------------------------------- This is a basic example on how to perform an orbital-optimized pCCD calculation in PyBEST. This script performs an orbital-optimized pCCD calculation on the water molecule using the cc-pVDZ basis set and canonical RHF orbitals as initial orbitals. Afterwards, the calculation is restarted from the default checkpoint file generated by the previous iteration. .. literalinclude:: ../src/pybest/data/examples/pccd/water_oopccd_cc-pvdz.py :caption: data/examples/pccd/water_oopccd_cc-pvdz.py :lines: 3-58 The next example shows how to perform a simple restart from some checkpoint file without any RHF calculation. .. literalinclude:: ../src/pybest/data/examples/pccd/water_oopccd_restart_cc-pvdz.py :caption: data/examples/pccd/water_oopccd_restart_cc-pvdz.py :lines: 3- .. note:: This will only work, if the restart is performed for the same molecular geometry. The water molecule including restart for a different geometry ------------------------------------------------------------- This is a basic example on how to perform an orbital-optimized pCCD calculation in PyBEST, where the orbitals from a **different** molecular geometry are taken as guess orbitals. .. literalinclude:: ../src/pybest/data/examples/pccd/water_oopccd_restart_perturbed_cc-pvdz.py :caption: data/examples/pccd/water_oopccd_restart_perturbed_cc-pvdz.py :lines: 3- .. note:: If the molecular structure is perturbed, we need to explicitly pass the external energy and the overlap matrix (as well as an empty set of orbitals). If not provided, the restart will result in wrong energies. pCCD with external integrals (FCIDUMP) -------------------------------------- This is a basic example on how to perform an orbital-optimized pCCD calculation using one- and two-electron integrals from an external file. The number of doubly-occupied orbitals is ``5`` and has to be specified, while the total number of basis functions is ``28`` (it is set automatically). .. literalinclude:: ../src/pybest/data/examples/pccd/water_oopccd_fcidump.py :caption: data/examples/pccd/water_oopccd_fcidump.py :lines: 3- pCCD using model Hamiltonians: Hubbard -------------------------------------- This is a basic example on how to perform an orbital-optimized pCCD calculation using the half-filled 1-D Hubbard model Hamiltonian with 6 sites. The number of doubly-occupied sites is thus ``3``. The ``t`` parameter is set to -1, the ``U`` parameter is set to 1. Periodic boundary conditions are employed. .. literalinclude:: ../src/pybest/data/examples/pccd/hubbard_oopccd.py :caption: data/examples/pccd/hubbard_oopccd.py :lines: 3-