Defining Fragments

Fragments are added within a fragmentation context manager, which determines in which way orbitals are constructed for particular fragment choices.

Fragmentation

In order to define fragments in an ab initio system, one needs to specify a set of orthogonal projectors onto subspaces. There is no unique way to define such projectors—the question where one atom or orbital ‘ends’ and the next one begins within a molecule or solid is fundamentally ill-posed.

Nevertheless, we do not require the projectors to be uniquely defined in order to perform quantum embedding calculations. Within the embedded wave function (EWF) framework, the systematic improvability to the exact result is guaranteed via an expansion of the bath space and appropriate choice of property functionals [1]. The only requirement on the fragments is that they are orthogonal and (when taken together) they span the complete occupied space. For DMET calculations, it is additionally desirable for the union of the fragment spaces to span the entire virtual space.

Note

To check if a user-defined fragmentation is orthogonal and complete in the occupied/virtual space, lock for the following line in the output:

Fragmentation: orthogonal= True, occupied-complete= True, virtual-complete= False

For ab initio systems atomic projectors can be defined in terms of three different types of local atomic orbitals:

For lattice models, on the other hand, the mapping between the site basis and the sites is clearer and more natural to define.

Table 1 shows a comparison of the the different fragmentation types:

Table 1 Comparison of fragmentation types

Type

Context manager

DMET

EWF

Comments

IAO

iao_fragmentation(minao='auto')

No

Yes

Default for EWF

IAO+PAO

iaopao_fragmentation(minao='auto')

Yes

Yes

SAO

sao_fragmentation()

Yes

Yes

Site basis

site_fragmentation()

Yes

Yes

For lattice models only

The minimal reference basis for the IAO construction can be selected with the minao argument. By default a suitable basis set will be chosen automatically 'auto'.

Note

Generally different fragmentations should not be combined, as the resulting fragments are not guaranteed to be orthogonal. The only exception to this are the IAO and IAO+PAO fragmentation, which can be combined as long as no atom is added twice:

with emb.iao_fragmentation() as f:
    f.add_atomic_fragment(0)
with emb.iaopao_fragmentation() as f:
    f.add_atomic_fragment(1)

Adding Fragments

Within the fragmentation context manager, fragments can be added using the methods add_atomic_fragment and add_orbital_fragment.

The add_atomic_fragment method can accept both atom indices and symbols and can further filter specific orbital types with the orbital_filter argument. The capabilities are best demonstrated in example examples/ewf/molecules/12-custom-fragments.py:

 1import pyscf
 2import pyscf.gto
 3import pyscf.scf
 4import pyscf.cc
 5
 6import vayesta
 7import vayesta.ewf
 8
 9mol = pyscf.gto.Mole()
10mol.atom = """
11Se	0.0000	0.0000	0.2807
12O 	0.0000	1.3464	-0.5965
13O 	0.0000	-1.3464	-0.5965
14"""
15mol.basis = "cc-pVDZ"
16mol.output = "pyscf.out"
17mol.build()
18
19# Hartree-Fock
20mf = pyscf.scf.RHF(mol)
21mf.kernel()
22
23# Reference full system CCSD:
24cc = pyscf.cc.CCSD(mf)
25cc.kernel()
26
27# Embedded CCSD
28emb = vayesta.ewf.EWF(mf, bath_options=dict(threshold=1e-6))
29with emb.iao_fragmentation() as f:
30    # Fragment containing the 1s state of O and 1s and 2s states of Se
31    f.add_atomic_fragment(["Se", "O"], orbital_filter=["1s", "Se 2s"])
32    # Atoms can be specified by labels or indices
33    # Fragment containing the 2s state at O and 3s and 4s states of Se
34    f.add_atomic_fragment([0, 1, 2], orbital_filter=["O 2s", "3s", "4s"])
35    # Fragment containing the 2p x- and y-states on the oxygen and the 2p, 3p y- and z- states on selenium
36    # Note that the oxygen does not have 3p IAO states, such that it is not necessary to specify the element for these states
37    f.add_atomic_fragment([0, 1, 2], orbital_filter=["O 2py", "O 2pz", "Se 2p", "3py", "3pz"])
38    # All 4p states on Se and all px states (2px on O, 2-4px on Se)
39    f.add_atomic_fragment(["Se", "O"], orbital_filter=["4p", "px"])
40    # 3d states on Se
41    f.add_atomic_fragment(0, orbital_filter=["3d"])
42emb.kernel()
43
44print("E(HF)=        %+16.8f Ha" % mf.e_tot)
45print("E(CCSD)=      %+16.8f Ha" % cc.e_tot)
46print("E(Emb. CCSD)= %+16.8f Ha" % emb.e_tot)

The add_orbital_fragment method allows selecting orbitals (IAOs, PAOs, or SAOs) to define fragments, either via their indices or labels.