Esempio n. 1
0
    def relocalize_states (self, loc2bas, fragments, oneRDM_loc, natorb=False, canonicalize=False):
        '''Do Boys localization on a subspace and assign resulting states to the various fragments using projection operators.
           Optionally diagonalize either the fock or the density matrix inside each subspace. Canonicalize overrides natorb'''

        fock_loc = self.loc_rhf_fock_bis (oneRDM_loc)
        ao2bas = boys.Boys (self.mol, np.dot (self.ao2loc, loc2bas)).kernel ()
        loc2bas = reduce (np.dot, [self.ao2loc.conjugate ().T, self.ao_ovlp, ao2bas])
        
        weights = np.asarray ([np.einsum ('ip,ip->p', loc2bas[f.frag_orb_list,:].conjugate (), loc2bas[f.frag_orb_list,:]) for f in fragments])
        frag_assignments = np.argmax (weights, axis=0)

        loc2bas_assigned = []        
        for idx, frag in enumerate (fragments):
            pick_orbs = (frag_assignments == idx)
            norbs = np.count_nonzero (pick_orbs)
            print ("{} states found for fragment {}".format (norbs, frag.frag_name))
            loc2pick = loc2bas[:,pick_orbs]
            if canonicalize and norbs:
                f = represent_operator_in_basis (fock_loc, loc2pick)
                evals, evecs = matrix_eigen_control_options (f, sort_vecs=1, only_nonzero_vals=False)
                loc2pick = np.dot (loc2pick, evecs)
            elif natorb and norbs:
                f = represent_operator_in_basis (oneRDM_loc, loc2pick)
                evals, evecs = matrix_eigen_control_options (f, sort_vecs=-1, only_nonzero_vals=False)
                loc2pick = np.dot (loc2pick, evecs)
            loc2bas_assigned.append (loc2pick)
        return loc2bas_assigned
Esempio n. 2
0
    def get_trial_nos (self, aobasis=False, loc2wmas=None, oneRDM_loc=None, fock=None, jmol_shift=False, try_symmetrize=True):
        if oneRDM_loc is None: oneRDM_loc = self.oneRDM_loc
        if fock is None:
            fock = self.activeFOCK
        elif isinstance (fock, str) and fock == 'calculate':
            fock = self.loc_rhf_fock_bis (oneRDM_loc)
        if loc2wmas is None: loc2wmas = [np.zeros ((self.norbs_tot, 0), dtype=self.ao2loc.dtype)]
        elif isinstance (loc2wmas, np.ndarray):
            if loc2wmas.ndim == 2: loc2wmas = loc2wmas[None,:,:]
            loc2wmas = [loc2amo for loc2amo in loc2wmas]
        occ_wmas = [np.zeros (0) for ix in loc2wmas]
        symm_wmas = [np.zeros (0) for ix in loc2wmas]
        for ix, loc2amo in enumerate (loc2wmas):
            occ_wmas[ix], loc2wmas[ix], symm_wmas[ix] = matrix_eigen_control_options (oneRDM_loc, symmetry=self.loc2symm, subspace=loc2amo,
                sort_vecs=-1, only_nonzero_vals=False, strong_symm=self.enforce_symmetry)
        occ_wmas = np.concatenate (occ_wmas)
        symm_wmas = np.concatenate (symm_wmas)
        loc2wmas = np.concatenate (loc2wmas, axis=-1)
        nelec_wmas = int (round (compute_nelec_in_subspace (oneRDM_loc, loc2wmas)))

        loc2wmcs = get_complementary_states (loc2wmas, symmetry=self.loc2symm, enforce_symmetry=self.enforce_symmetry)
        norbs_wmas = loc2wmas.shape[1]
        norbs_wmcs = loc2wmcs.shape[1]
        ene_wmcs, loc2wmcs, symm_wmcs = matrix_eigen_control_options (fock, symmetry=self.loc2symm, subspace=loc2wmcs, sort_vecs=1, only_nonzero_vals=False, strong_symm=self.enforce_symmetry)
            
        assert ((self.nelec_tot - nelec_wmas) % 2 == 0), 'Non-even number of unactive electrons {}'.format (self.nelec_tot - nelec_wmas)
        norbs_core = (self.nelec_tot - nelec_wmas) // 2
        norbs_virt = norbs_wmcs - norbs_core
        loc2wmis = loc2wmcs[:,:norbs_core]
        symm_wmis = symm_wmcs[:norbs_core]
        loc2wmxs = loc2wmcs[:,norbs_core:]
        symm_wmxs = symm_wmcs[norbs_core:]
        
        if self.mol.symmetry:
            symm_wmis = {self.mol.irrep_name[x]: np.count_nonzero (symm_wmis==x) for x in np.unique (symm_wmis)}
            err = measure_subspace_blockbreaking (loc2wmis, self.loc2symm)
            print ("Trial wave function inactive-orbital irreps = {}, err = {}".format (symm_wmis, err))
            symm_wmas = {self.mol.irrep_name[x]: np.count_nonzero (symm_wmas==x) for x in np.unique (symm_wmas)} 
            err = measure_subspace_blockbreaking (loc2wmas, self.loc2symm)
            print ("Trial wave function active-orbital irreps = {}, err = {}".format (symm_wmas, err))
            symm_wmxs = {self.mol.irrep_name[x]: np.count_nonzero (symm_wmxs==x) for x in np.unique (symm_wmxs)}
            err = measure_subspace_blockbreaking (loc2wmxs, self.loc2symm)
            print ("Trial wave function external-orbital irreps = {}, err = {}".format (symm_wmxs, err))

        loc2no = np.concatenate ((loc2wmcs[:,:norbs_core], loc2wmas, loc2wmcs[:,norbs_core:]), axis=1)
        occ_no = np.concatenate ((2*np.ones (norbs_core), occ_wmas, np.zeros (norbs_virt)))
        ene_no = np.concatenate ((ene_wmcs[:norbs_core], np.zeros (norbs_wmas), ene_wmcs[norbs_core:]))
        assert (len (occ_no) == len (ene_no) and loc2no.shape[1] == len (occ_no)), '{} {} {}'.format (loc2no.shape, len (ene_no), len (occ_no))
        norbs_occ = norbs_core + norbs_wmas
        if jmol_shift:
            print ("Shifting natural-orbital energies so that jmol puts them in the correct order:")
            if ene_no[norbs_core-1] > 0: ene_no[:norbs_core] -= ene_no[norbs_core-1] + 1e-6
            if ene_no[norbs_occ] < 0: ene_no[norbs_occ:] -= ene_no[norbs_occ] - 1e-6
            assert (np.all (np.diff (ene_no) >=0)), ene_no
        if aobasis:
            return self.ao2loc @ loc2no, ene_no, occ_no
        return loc2no, ene_no, occ_no
Esempio n. 3
0
def orthonormalize_a_basis (overlapping_basis, ovlp=1, num_zero_atol=params.num_zero_ltol, symmetry=None, enforce_symmetry=False):
    if (is_basis_orthonormal (overlapping_basis)):
        return overlapping_basis
    c2b = np.asarray (overlapping_basis)
    cOc = np.asarray (ovlp)
    if enforce_symmetry:
        c2n = np.zeros ((overlapping_basis.shape[0], 0), dtype=overlapping_basis.dtype)
        for c2s in symmetry:
            s2c = c2s.conjugate ().T
            s2b = s2c @ c2b
            sOs = s2c @ cOc @ c2s if cOc.shape == ((c2b.shape[0], c2b.shape[0])) else (s2c * cOc) @ c2s
            s2n = orthonormalize_a_basis (s2b, ovlp=sOs, num_zero_atol=num_zero_atol, symmetry=None, enforce_symmetry=False)
            c2n = np.append (c2n, c2s @ s2n, axis=1)
        return (c2n)

    b2c = c2b.conjugate ().T
    bOb = b2c @ cOc @ c2b if cOc.shape == ((c2b.shape[0], c2b.shape[0])) else (b2c * cOc) @ c2b
    assert (not is_matrix_zero (bOb)), "overlap matrix is zero! problem with basis?"
    assert (np.allclose (bOb, bOb.conjugate ().T)), "overlap matrix not hermitian! problem with basis?"
    assert (np.abs (np.trace (bOb)) > num_zero_atol), "overlap matrix zero or negative trace! problem with basis?"
     
    evals, evecs = matrix_eigen_control_options (bOb, sort_vecs=-1, only_nonzero_vals=True, num_zero_atol=num_zero_atol)
    if len (evals) == 0:
        return np.zeros ((c2b.shape[0], 0), dtype=c2b.dtype)
    p2x = np.asarray (evecs)
    c2x = c2b @ p2x 
    assert (not np.any (evals < 0)), "overlap matrix has negative eigenvalues! problem with basis?"

    # I want c2n = c2x * x2n
    # x2n defined such that n2c * c2n = I
    # n2x * x2c * c2x * x2n = n2x * evals_xx * x2n = I
    # therefore
    # x2n = evals_xx^{-1/2}
    x2n = np.asarray (np.diag (np.reciprocal (np.sqrt (evals))))
    c2n = c2x @ x2n
    n2c = c2n.conjugate ().T
    nOn = n2c @ cOc @ c2n if cOc.shape == ((c2b.shape[0], c2b.shape[0])) else (n2c * cOc) @ c2n
    if not is_basis_orthonormal (c2n):
        # Assuming numerical problem due to massive degeneracy; remove constant from diagonal to improve solver?
        assert (np.all (np.isclose (np.diag (nOn), 1))), np.diag (nOn) - 1
        nOn[np.diag_indices_from (nOn)] -= 1
        evals, evecs = matrix_eigen_control_options (nOn, sort_vecs=-1, only_nonzero_vals=False)
        n2x = np.asarray (evecs)
        c2x = c2n @ n2x
        x2n = np.asarray (np.diag (np.reciprocal (np.sqrt (evals + 1))))
        c2n = c2x @ x2n
        n2c = c2n.conjugate ().T
        nOn = n2c @ cOc @ c2n if cOc.shape == ((c2b.shape[0], c2b.shape[0])) else (n2c * cOc) @ c2n
        assert (is_basis_orthonormal (c2n)), "failed to orthonormalize basis even after two tries somehow\n" + str (
            prettyprint_ndarray (nOn)) + "\n" + str (np.linalg.norm (nOn - np.eye (c2n.shape[1]))) + "\n" + str (evals)

    return np.asarray (c2n)
Esempio n. 4
0
def get_unfrozen_states(oneRDMfroz_loc):
    _, loc2froz = matrix_eigen_control_options(oneRDMfroz_loc,
                                               only_nonzero_vals=True)
    if loc2froz.shape[1] == loc2froz.shape[0]:
        raise RuntimeError(
            "No unfrozen states: eigenbasis of oneRDMfroz_loc is complete!")
    return get_complementary_states(loc2froz)
Esempio n. 5
0
def eigen_weaksymm(the_matrix,
                   the_blocks,
                   subspace=None,
                   sort_vecs=1,
                   only_nonzero_vals=False,
                   atol=params.num_zero_atol,
                   rtol=params.num_zero_rtol):
    if the_blocks is None: the_blocks = [np.eye(the_matrix.shape[0])]
    if subspace is None: subspace = np.eye(the_matrix.shape[0])
    subspace_matrix = represent_operator_in_basis(the_matrix, subspace)
    evals, evecs = matrix_eigen_control_options(
        subspace_matrix,
        symm_blocks=None,
        sort_vecs=sort_vecs,
        only_nonzero_vals=only_nonzero_vals,
        num_zero_atol=atol)
    evecs = subspace @ evecs
    idx_unchk = np.ones(len(evals), dtype=np.bool_)
    while np.count_nonzero(idx_unchk > 0):
        chk_1st_eval = evals[idx_unchk][0]
        idx_degen = np.isclose(evals, chk_1st_eval, rtol=rtol, atol=atol)
        if np.count_nonzero(idx_degen) > 1:
            evecs[:, idx_degen] = align_states(evecs[:, idx_degen],
                                               the_blocks,
                                               atol=atol,
                                               rtol=rtol)
        idx_unchk[idx_degen] = False
    return evals, evecs, assign_blocks_weakly(evecs, the_blocks)
Esempio n. 6
0
def get_states_from_projector (the_projector, num_zero_atol=params.num_zero_atol):
    proj_cc = np.asarray (the_projector)
    assert (np.allclose (proj_cc, proj_cc.H)), "projector must be hermitian\n" + str (np.linalg.norm (proj_cc - proj_cc.conjugate ().T))
    assert (is_matrix_idempotent (proj_cc)), "projector must be idempotent\n" + str (np.linalg.norm ((proj_cc @ proj_cc) - proj_cc))
    evals, evecs = matrix_eigen_control_options (proj_cc, sort_vecs=-1, only_nonzero_vals=True, num_zero_atol=num_zero_atol)
    idx = np.isclose (evals, 1)
    return evecs[:,idx]
Esempio n. 7
0
def count_linind_states (the_states, ovlp=1, num_zero_atol=params.num_zero_atol):
    c2b = np.asarray (the_states)
    b2c = c2b.conjugate ().T
    cOc = np.asarray (ovlp)
    nbas = c2b.shape[0]
    nstates = c2b.shape[1]
    bOb = b2c @ cOc @ c2b if cOc.shape == ((nbas, nbas)) else b2c @ c2b
    if is_matrix_zero (bOb) or np.abs (np.trace (bOb)) <= num_zero_atol: return 0
    evals = matrix_eigen_control_options (bOb, only_nonzero_vals=True)[0]
    return len (evals)
Esempio n. 8
0
def symmetrize_basis(the_basis,
                     the_blocks,
                     sorting_metric=None,
                     sort_vecs=1,
                     do_eigh_metric=True,
                     check_metric_block_adapted=True,
                     atol=params.num_zero_atol,
                     rtol=params.num_zero_rtol):
    atol_scl = atol * the_basis.shape[0]
    rtol_scl = rtol * the_basis.shape[0]
    if the_blocks is None: the_blocks = [np.eye(the_basis.shape[0])]
    assert (is_subspace_block_adapted(the_basis, the_blocks, tol=atol)
            ), 'Basis space must be block-adapted before blockifying states'
    symmetrized_basis = align_states(the_basis, the_blocks)
    labels = assign_blocks(symmetrized_basis, the_blocks)
    assert (is_basis_orthonormal(symmetrized_basis, atol=atol,
                                 rtol=rtol)), "? labels = {}".format(labels)

    if sorting_metric is None:
        return symmetrized_basis, labels
    else:
        if sorting_metric.shape[0] == the_basis.shape[0]:
            metric_symm = represent_operator_in_basis(sorting_metric,
                                                      symmetrized_basis)
        else:
            assert (
                sorting_metric.shape[0] == the_basis.shape[1]
            ), 'The sorting metric must be in either the row or column basis of the orbital matrix that is being symmetrized'
            metric_symm = represent_operator_in_basis(
                sorting_metric,
                the_basis.conjugate().T @ symmetrized_basis)
        if check_metric_block_adapted:
            assert (is_operator_block_adapted(metric_symm, labels, tol=atol))
        metric_evals, evecs, labels = matrix_eigen_control_options(
            metric_symm,
            symm_blocks=labels,
            sort_vecs=sort_vecs,
            only_nonzero_vals=False,
            num_zero_atol=atol)
        symmetrized_basis = symmetrized_basis @ evecs
        return symmetrized_basis, labels, metric_evals
Esempio n. 9
0
    def setup_wm_core_scf(self, fragments, calcname):

        self.restore_wm_full_scf()
        oneRDMcorr_loc = sum((frag.oneRDMas_loc for frag in fragments))
        if np.all(np.isclose(oneRDMcorr_loc, 0)):
            print("Null correlated 1-RDM; default settings for wm wvfn")
            self.activeFOCK = represent_operator_in_basis(
                self.fullFOCKao, self.ao2loc)
            self.activeJKidem = self.activeFOCK - self.activeOEI
            self.activeJKcorr = np.zeros((self.norbs_tot, self.norbs_tot))
            self.oneRDMcorr_loc = oneRDMcorr_loc
            self.loc2idem = np.eye(self.norbs_tot)
            self.nelec_idem = self.nelec_tot
            return

        loc2corr = np.concatenate([frag.loc2amo for frag in fragments], axis=1)
        loc2idem = get_complementary_states(loc2corr)
        evecs = matrix_eigen_control_options(represent_operator_in_basis(
            self.loc_oei(), loc2idem),
                                             sort_vecs=1,
                                             only_nonzero_vals=False)[1]
        loc2idem = np.dot(loc2idem, evecs)

        # I want to alter the outputs of self.loc_oei (), self.loc_rhf_fock (), and the get_wm_1RDM_etc () functions.
        # self.loc_oei ()      = P_idem * (activeOEI + JKcorr) * P_idem
        # self.loc_rhf_fock () = P_idem * (activeOEI + JKcorr + JKidem) * P_idem
        # The get_wm_1RDM_etc () functions will need to add oneRDMcorr_loc to their final return value
        # The chemical potential is so that identically zero eigenvalues from the projection into the idem space don't get confused
        # with numerically-zero eigenvalues in the idem space: all occupied orbitals must have negative energy

        # Make true output 1RDM from fragments to use as guess for wm mcscf calculation
        oneRDMguess_loc = np.zeros_like(oneRDMcorr_loc)
        for f in itertools.product(fragments, fragments):
            loc2frag = [i.loc2frag for i in f]
            oneRDMguess_loc += sum(
                (0.5 * project_operator_into_subspace(i.oneRDM_loc, *loc2frag)
                 for i in f))

        nelec_corr = np.trace(oneRDMcorr_loc)
        if is_close_to_integer(nelec_corr,
                               100 * params.num_zero_atol) == False:
            raise ValueError(
                "nelec_corr not an integer! {}".format(nelec_corr))
        nelec_idem = int(round(self.nelec_tot - nelec_corr))
        JKcorr = self.loc_rhf_jk_bis(oneRDMcorr_loc)
        oneRDMidem_loc = self.get_wm_1RDM_from_scf_on_OEI(
            self.loc_oei() + JKcorr,
            nelec=nelec_idem,
            loc2wrk=loc2idem,
            oneRDMguess_loc=oneRDMguess_loc,
            output=calcname + '_trial_wvfn.log')
        JKidem = self.loc_rhf_jk_bis(oneRDMidem_loc)
        print("trace of oneRDMcorr_loc = {}".format(np.trace(oneRDMcorr_loc)))
        print("trace of oneRDMidem_loc = {}".format(np.trace(oneRDMidem_loc)))
        print("trace of oneRDM_loc in corr basis = {}".format(
            np.trace(
                represent_operator_in_basis(
                    oneRDMcorr_loc + oneRDMidem_loc,
                    orthonormalize_a_basis(loc2corr)))))
        svals = get_overlapping_states(loc2idem, loc2corr)[2]
        print("trace of <idem|corr|idem> = {}".format(np.sum(svals * svals)))
        print(loc2corr.shape)
        print(loc2idem.shape)

        ########################################################################################################
        self.activeFOCK = self.activeOEI + JKidem + JKcorr
        self.activeJKidem = JKidem
        self.activeJKcorr = JKcorr
        self.oneRDMcorr_loc = oneRDMcorr_loc
        self.loc2idem = loc2idem
        self.nelec_idem = nelec_idem
        ########################################################################################################

        # Analysis: 1RDM and total energy
        print("Analyzing LASSCF trial wave function")
        oei = self.activeOEI + (JKcorr + JKidem) / 2
        fock = self.activeFOCK
        oneRDM = oneRDMidem_loc + oneRDMcorr_loc
        E = self.activeCONST + np.tensordot(oei, oneRDM, axes=2)
        for frag in fragments:
            if frag.norbs_as > 0:
                if frag.E2_cum == 0 and np.amax(np.abs(
                        frag.twoCDMimp_amo)) > 0:
                    V = self.dmet_tei(frag.loc2amo)
                    L = frag.twoCDMimp_amo
                    frag.E2_cum = np.tensordot(V, L, axes=4) / 2
                E += frag.E2_cum
        print("LASSCF trial wave function total energy: {:.6f}".format(E))
        self.oneRDM_loc = oneRDM
        self.e_tot = E

        # Molden
        fock_idem = represent_operator_in_basis(fock, loc2idem)
        oneRDM_corr = represent_operator_in_basis(oneRDM, loc2corr)
        idem_evecs = matrix_eigen_control_options(fock_idem,
                                                  sort_vecs=1,
                                                  only_nonzero_vals=False)[1]
        corr_evecs = matrix_eigen_control_options(oneRDM_corr,
                                                  sort_vecs=-1,
                                                  only_nonzero_vals=False)[1]
        loc2molden = np.append(np.dot(loc2idem, idem_evecs),
                               np.dot(loc2corr, corr_evecs),
                               axis=1)
        wm_ene = np.einsum('ip,ij,jp->p', loc2molden, fock, loc2molden)
        wm_ene[-loc2corr.shape[1]:] = 0
        wm_occ = np.einsum('ip,ij,jp->p', loc2molden, oneRDM, loc2molden)
        ao2molden = np.dot(self.ao2loc, loc2molden)
        molden.from_mo(self.mol,
                       calcname + '_trial_wvfn.molden',
                       ao2molden,
                       occ=wm_occ,
                       ene=wm_ene)
Esempio n. 10
0
def solve (frag, guess_1RDM, chempot_imp):

    # Augment OEI with the chemical potential
    OEI = frag.impham_OEI_C - chempot_imp

    # Do I need to get the full RHF solution?
    guess_orbs_av = len (frag.imp_cache) == 2 or frag.norbs_as > 0 

    # Get the RHF solution
    mol = gto.Mole()
    abs_2MS = int (round (2 * abs (frag.target_MS)))
    abs_2S = int (round (2 * abs (frag.target_S)))
    sign_MS = int (np.sign (frag.target_MS)) or 1
    mol.spin = abs_2MS
    mol.verbose = 0 
    if frag.mol_stdout is None:
        mol.output = frag.mol_output
        mol.verbose = 0 if frag.mol_output is None else lib.logger.DEBUG
    mol.atom.append(('H', (0, 0, 0)))
    mol.nelectron = frag.nelec_imp
    if frag.enforce_symmetry:
        mol.groupname  = frag.symmetry
        mol.symm_orb   = get_subspace_symmetry_blocks (frag.loc2imp, frag.loc2symm)
        mol.irrep_name = frag.ir_names
        mol.irrep_id   = frag.ir_ids
    mol.max_memory = frag.ints.max_memory
    mol.build ()
    if frag.mol_stdout is None:
        frag.mol_stdout = mol.stdout
    else:
        mol.stdout = frag.mol_stdout
        mol.verbose = 0 if frag.mol_output is None else lib.logger.DEBUG
    if frag.enforce_symmetry: mol.symmetry = True
    #mol.incore_anyway = True
    mf = scf.RHF(mol)
    mf.get_hcore = lambda *args: OEI
    mf.get_ovlp = lambda *args: np.eye(frag.norbs_imp)
    mf.energy_nuc = lambda *args: frag.impham_CONST
    if frag.impham_CDERI is not None:
        mf = mf.density_fit ()
        mf.with_df._cderi = frag.impham_CDERI
    else:
        mf._eri = ao2mo.restore(8, frag.impham_TEI, frag.norbs_imp)
    mf = fix_my_RHF_for_nonsinglet_env (mf, frag.impham_OEI_S)
    mf.__dict__.update (frag.mf_attr)
    if guess_orbs_av: mf.max_cycle = 2
    mf.scf (guess_1RDM)
    if (not mf.converged) and (not guess_orbs_av):
        if np.any (np.abs (frag.impham_OEI_S) > 1e-8) and mol.spin != 0:
            raise NotImplementedError('Gradient and Hessian fixes for nonsinglet environment of Newton-descent ROHF algorithm')
        print ("CASSCF RHF-step not converged on fixed-point iteration; initiating newton solver")
        mf = mf.newton ()
        mf.kernel ()

    # Instability check and repeat
    if not guess_orbs_av:
        for i in range (frag.num_mf_stab_checks):
            if np.any (np.abs (frag.impham_OEI_S) > 1e-8) and mol.spin != 0:
                raise NotImplementedError('ROHF stability-check fixes for nonsinglet environment')
            mf.mo_coeff = mf.stability ()[0]
            guess_1RDM = mf.make_rdm1 ()
            mf = scf.RHF(mol)
            mf.get_hcore = lambda *args: OEI
            mf.get_ovlp = lambda *args: np.eye(frag.norbs_imp)
            mf._eri = ao2mo.restore(8, frag.impham_TEI, frag.norbs_imp)
            mf = fix_my_RHF_for_nonsinglet_env (mf, frag.impham_OEI_S)
            mf.scf (guess_1RDM)
            if not mf.converged:
                mf = mf.newton ()
                mf.kernel ()

    E_RHF = mf.e_tot
    print ("CASSCF RHF-step energy: {}".format (E_RHF))

    # Get the CASSCF solution
    CASe = frag.active_space[0]
    CASorb = frag.active_space[1] 
    checkCAS =  (CASe <= frag.nelec_imp) and (CASorb <= frag.norbs_imp)
    if (checkCAS == False):
        CASe = frag.nelec_imp
        CASorb = frag.norbs_imp
    if (abs_2MS > abs_2S):
        CASe = ((CASe + sign_MS * abs_2S) // 2, (CASe - sign_MS * abs_2S) // 2)
    else:
        CASe = ((CASe + sign_MS * abs_2MS) // 2, (CASe - sign_MS * abs_2MS) // 2)
    if frag.impham_CDERI is not None:
        mc = mcscf.DFCASSCF(mf, CASorb, CASe)
    else:
        mc = mcscf.CASSCF(mf, CASorb, CASe)
    smult = abs_2S + 1 if frag.target_S is not None else (frag.nelec_imp % 2) + 1
    mc.fcisolver = csf_solver (mf.mol, smult, symm=frag.enforce_symmetry)
    if frag.enforce_symmetry: mc.fcisolver.wfnsym = frag.wfnsym
    mc.max_cycle_macro = 50 if frag.imp_maxiter is None else frag.imp_maxiter
    mc.conv_tol = min (1e-9, frag.conv_tol_grad**2)  
    mc.ah_start_tol = mc.conv_tol / 10
    mc.ah_conv_tol = mc.conv_tol / 10
    mc.__dict__.update (frag.corr_attr)
    mc = fix_my_CASSCF_for_nonsinglet_env (mc, frag.impham_OEI_S)
    norbs_amo = mc.ncas
    norbs_cmo = mc.ncore
    norbs_imo = frag.norbs_imp - norbs_amo
    nelec_amo = sum (mc.nelecas)
    norbs_occ = norbs_amo + norbs_cmo
    #mc.natorb = True

    # Guess orbitals
    ci0 = None
    dm_imp = frag.get_oneRDM_imp ()
    fock_imp = mf.get_fock (dm=dm_imp)
    if len (frag.imp_cache) == 2:
        imp2mo, ci0 = frag.imp_cache
        print ("Taking molecular orbitals and ci vector from cache")
    elif frag.norbs_as > 0:
        nelec_imp_guess = int (round (np.trace (frag.oneRDMas_loc)))
        norbs_cmo_guess = (frag.nelec_imp - nelec_imp_guess) // 2
        print ("Projecting stored amos (frag.loc2amo; spanning {} electrons) onto the impurity basis and filling the remainder with default guess".format (nelec_imp_guess))
        imp2mo, my_occ = project_amo_manually (frag.loc2imp, frag.loc2amo, fock_imp, norbs_cmo_guess, dm=frag.oneRDMas_loc)
    elif frag.loc2amo_guess is not None:
        print ("Projecting stored amos (frag.loc2amo_guess) onto the impurity basis (no amo dm available)")
        imp2mo, my_occ = project_amo_manually (frag.loc2imp, frag.loc2amo_guess, fock_imp, norbs_cmo, dm=None)
        frag.loc2amo_guess = None
    else:
        dm_imp = np.asarray (mf.make_rdm1 ())
        while dm_imp.ndim > 2:
            dm_imp = dm_imp.sum (0)
        imp2mo = mf.mo_coeff
        fock_imp = mf.get_fock (dm=dm_imp)
        fock_mo = represent_operator_in_basis (fock_imp, imp2mo)
        _, evecs = matrix_eigen_control_options (fock_mo, sort_vecs=1)
        imp2mo = imp2mo @ evecs
        my_occ = ((dm_imp @ imp2mo) * imp2mo).sum (0)
        print ("No stored amos; using mean-field canonical MOs as initial guess")
    # Guess orbital processing
    if callable (frag.cas_guess_callback):
        mo = reduce (np.dot, (frag.ints.ao2loc, frag.loc2imp, imp2mo))
        mo = frag.cas_guess_callback (frag.ints.mol, mc, mo)
        imp2mo = reduce (np.dot, (frag.imp2loc, frag.ints.ao2loc.conjugate ().T, frag.ints.ao_ovlp, mo))
        frag.cas_guess_callback = None

    # Guess CI vector
    if len (frag.imp_cache) != 2 and frag.ci_as is not None:
        loc2amo_guess = np.dot (frag.loc2imp, imp2mo[:,norbs_cmo:norbs_occ])
        metric = np.arange (CASorb) + 1
        gOc = np.dot (loc2amo_guess.conjugate ().T, (frag.ci_as_orb * metric[None,:]))
        umat_g, svals, umat_c = matrix_svd_control_options (gOc, sort_vecs=1, only_nonzero_vals=True)
        if (svals.size == norbs_amo):
            print ("Loading ci guess despite shifted impurity orbitals; singular value error sum: {}".format (np.sum (svals - metric)))
            imp2mo[:,norbs_cmo:norbs_occ] = np.dot (imp2mo[:,norbs_cmo:norbs_occ], umat_g)
            ci0 = transform_ci_for_orbital_rotation (frag.ci_as, CASorb, CASe, umat_c)
        else:
            print ("Discarding stored ci guess because orbitals are too different (missing {} nonzero svals)".format (norbs_amo-svals.size))

    # Symmetry align if possible
    imp2unac = frag.align_imporbs_symm (np.append (imp2mo[:,:norbs_cmo], imp2mo[:,norbs_occ:], axis=1), sorting_metric=fock_imp,
        sort_vecs=1, orbital_type='guess unactive', mol=mol)[0]
    imp2mo[:,:norbs_cmo] = imp2unac[:,:norbs_cmo]
    imp2mo[:,norbs_occ:] = imp2unac[:,norbs_cmo:]
    #imp2mo[:,:norbs_cmo] = frag.align_imporbs_symm (imp2mo[:,:norbs_cmo], sorting_metric=fock_imp, sort_vecs=1, orbital_type='guess inactive', mol=mol)[0]
    imp2mo[:,norbs_cmo:norbs_occ], umat = frag.align_imporbs_symm (imp2mo[:,norbs_cmo:norbs_occ], sorting_metric=fock_imp,
        sort_vecs=1, orbital_type='guess active', mol=mol)
    #imp2mo[:,norbs_occ:] = frag.align_imporbs_symm (imp2mo[:,norbs_occ:], sorting_metric=fock_imp, sort_vecs=1, orbital_type='guess external', mol=mol)[0]
    if frag.enforce_symmetry:
        imp2mo = cleanup_subspace_symmetry (imp2mo, mol.symm_orb)
        err_symm = measure_subspace_blockbreaking (imp2mo, mol.symm_orb)
        err_orth = measure_basis_nonorthonormality (imp2mo)
        print ("Initial symmetry error after cleanup = {}".format (err_symm))
        print ("Initial orthonormality error after cleanup = {}".format (err_orth))
    if ci0 is not None: ci0 = transform_ci_for_orbital_rotation (ci0, CASorb, CASe, umat)
        

    # Guess orbital printing
    if frag.mfmo_printed == False and frag.ints.mol.verbose:
        ao2mfmo = reduce (np.dot, [frag.ints.ao2loc, frag.loc2imp, imp2mo])
        print ("Writing {} {} orbital molden".format (frag.frag_name, 'CAS guess'))
        molden.from_mo (frag.ints.mol, frag.filehead + frag.frag_name + '_mfmorb.molden', ao2mfmo, occ=my_occ)
        frag.mfmo_printed = True
    elif len (frag.active_orb_list) > 0: # This is done AFTER everything else so that the _mfmorb.molden always has consistent ordering
        print('Applying caslst: {}'.format (frag.active_orb_list))
        imp2mo = mc.sort_mo(frag.active_orb_list, mo_coeff=imp2mo)
        frag.active_orb_list = []
    if len (frag.frozen_orb_list) > 0:
        mc.frozen = copy.copy (frag.frozen_orb_list)
        print ("Applying frozen-orbital list (this macroiteration only): {}".format (frag.frozen_orb_list))
        frag.frozen_orb_list = []

    if frag.enforce_symmetry: imp2mo = lib.tag_array (imp2mo, orbsym=label_orb_symm (mol, mol.irrep_id, mol.symm_orb, imp2mo, s=mf.get_ovlp (), check=False))

    t_start = time.time()
    E_CASSCF = mc.kernel(imp2mo, ci0)[0]
    if (not mc.converged) and np.all (np.abs (frag.impham_OEI_S) < 1e-8):
        mc = mc.newton ()
        E_CASSCF = mc.kernel(mc.mo_coeff, mc.ci)[0]
    if not mc.converged:
        print ('Assuming ci vector is poisoned; discarding...')
        imp2mo = mc.mo_coeff.copy ()
        mc = mcscf.CASSCF(mf, CASorb, CASe)
        smult = abs_2S + 1 if frag.target_S is not None else (frag.nelec_imp % 2) + 1
        mc.fcisolver = csf_solver (mf.mol, smult)
        E_CASSCF = mc.kernel(imp2mo)[0]
        if not mc.converged:
            if np.any (np.abs (frag.impham_OEI_S) > 1e-8):
                raise NotImplementedError('Gradient and Hessian fixes for nonsinglet environment of Newton-descent CASSCF algorithm')
            mc = mc.newton ()
            E_CASSCF = mc.kernel(mc.mo_coeff, mc.ci)[0]
    assert (mc.converged)

    '''
    mc.conv_tol = 1e-12
    mc.ah_start_tol = 1e-10
    mc.ah_conv_tol = 1e-12
    E_CASSCF = mc.kernel(mc.mo_coeff, mc.ci)[0]
    if not mc.converged:
        mc = mc.newton ()
        E_CASSCF = mc.kernel(mc.mo_coeff, mc.ci)[0]
    #assert (mc.converged)
    '''
    
    # Get twoRDM + oneRDM. cs: MC-SCF core, as: MC-SCF active space
    # I'm going to need to keep some representation of the active-space orbitals

    # Symmetry align if possible
    oneRDM_amo, twoRDM_amo = mc.fcisolver.make_rdm12 (mc.ci, mc.ncas, mc.nelecas)
    fock_imp = mc.get_fock ()
    mc.mo_coeff[:,:norbs_cmo] = frag.align_imporbs_symm (mc.mo_coeff[:,:norbs_cmo], sorting_metric=fock_imp, sort_vecs=1, orbital_type='optimized inactive', mol=mol)[0]
    mc.mo_coeff[:,norbs_cmo:norbs_occ], umat = frag.align_imporbs_symm (mc.mo_coeff[:,norbs_cmo:norbs_occ],
        sorting_metric=oneRDM_amo, sort_vecs=-1, orbital_type='optimized active', mol=mol)
    mc.mo_coeff[:,norbs_occ:] = frag.align_imporbs_symm (mc.mo_coeff[:,norbs_occ:], sorting_metric=fock_imp, sort_vecs=1, orbital_type='optimized external', mol=mol)[0]
    if frag.enforce_symmetry:
        amo2imp = mc.mo_coeff[:,norbs_cmo:norbs_occ].conjugate ().T
        mc.mo_coeff = cleanup_subspace_symmetry (mc.mo_coeff, mol.symm_orb)
        umat = umat @ (amo2imp @ mc.mo_coeff[:,norbs_cmo:norbs_occ])
        err_symm = measure_subspace_blockbreaking (mc.mo_coeff, mol.symm_orb)
        err_orth = measure_basis_nonorthonormality (mc.mo_coeff)
        print ("Final symmetry error after cleanup = {}".format (err_symm))
        print ("Final orthonormality error after cleanup = {}".format (err_orth))
    mc.ci = transform_ci_for_orbital_rotation (mc.ci, CASorb, CASe, umat)

    # Cache stuff
    imp2mo = mc.mo_coeff #mc.cas_natorb()[0]
    loc2mo = np.dot (frag.loc2imp, imp2mo)
    imp2amo = imp2mo[:,norbs_cmo:norbs_occ]
    loc2amo = loc2mo[:,norbs_cmo:norbs_occ]
    frag.imp_cache = [mc.mo_coeff, mc.ci]
    frag.ci_as = mc.ci
    frag.ci_as_orb = loc2amo.copy ()
    t_end = time.time()

    # oneRDM
    oneRDM_imp = mc.make_rdm1 ()

    # twoCDM
    oneRDM_amo, twoRDM_amo = mc.fcisolver.make_rdm12 (mc.ci, mc.ncas, mc.nelecas)
    oneRDMs_amo = np.stack (mc.fcisolver.make_rdm1s (mc.ci, mc.ncas, mc.nelecas), axis=0)
    oneSDM_amo = oneRDMs_amo[0] - oneRDMs_amo[1] if frag.target_MS >= 0 else oneRDMs_amo[1] - oneRDMs_amo[0]
    oneSDM_imp = represent_operator_in_basis (oneSDM_amo, imp2amo.conjugate ().T)
    print ("Norm of spin density: {}".format (linalg.norm (oneSDM_amo)))
    # Note that I do _not_ do the *real* cumulant decomposition; I do one assuming oneSDM_amo = 0.
    # This is fine as long as I keep it consistent, since it is only in the orbital gradients for this impurity that
    # the spin density matters. But it has to stay consistent!
    twoCDM_amo = get_2CDM_from_2RDM (twoRDM_amo, oneRDM_amo)
    twoCDM_imp = represent_operator_in_basis (twoCDM_amo, imp2amo.conjugate ().T)
    print('Impurity CASSCF energy (incl chempot): {}; spin multiplicity: {}; time to solve: {}'.format (E_CASSCF, spin_square (mc)[1], t_end - t_start))

    # Active-space RDM data
    frag.oneRDMas_loc  = symmetrize_tensor (represent_operator_in_basis (oneRDM_amo, loc2amo.conjugate ().T))
    frag.oneSDMas_loc  = symmetrize_tensor (represent_operator_in_basis (oneSDM_amo, loc2amo.conjugate ().T))
    frag.twoCDMimp_amo = twoCDM_amo
    frag.loc2mo  = loc2mo
    frag.loc2amo = loc2amo
    frag.E2_cum  = np.tensordot (ao2mo.restore (1, mc.get_h2eff (), mc.ncas), twoCDM_amo, axes=4) / 2
    frag.E2_cum += (mf.get_k (dm=oneSDM_imp) * oneSDM_imp).sum () / 4
    # The second line compensates for my incorrect cumulant decomposition. Anything to avoid changing the checkpoint files...

    # General impurity data
    frag.oneRDM_loc = frag.oneRDMfroz_loc + symmetrize_tensor (represent_operator_in_basis (oneRDM_imp, frag.imp2loc))
    frag.oneSDM_loc = frag.oneSDMfroz_loc + frag.oneSDMas_loc
    frag.twoCDM_imp = None # Experiment: this tensor is huge. Do I actually need to keep it? In principle, of course not.
    frag.E_imp      = E_CASSCF + np.einsum ('ab,ab->', chempot_imp, oneRDM_imp)

    return None
Esempio n. 11
0
def project_amo_manually (loc2imp, loc2gamo, fock_mf, norbs_cmo, dm=None):
    norbs_amo = loc2gamo.shape[1]
    amo2imp = np.dot (loc2gamo.conjugate ().T, loc2imp)
    ovlp = np.dot (amo2imp, amo2imp.conjugate ().T)
    '''
    print ("Do impurity orbitals span guess amos?")
    print (prettyprint (ovlp, fmt='{:5.2f}'))
    if dm is not None:
        print ("Density matrix?")
        print (prettyprint (represent_operator_in_basis (dm, loc2gamo), fmt='{:5.2f}'))
    '''
    proj = np.dot (amo2imp.conjugate ().T, amo2imp)
    evals, evecs = matrix_eigen_control_options (proj, sort_vecs=-1, only_nonzero_vals=False)
    imp2amo = np.copy (evecs[:,:norbs_amo])
    imp2imo = np.copy (evecs[:,norbs_amo:])
    fock_imo = represent_operator_in_basis (fock_mf, imp2imo)
    _, evecs = matrix_eigen_control_options (fock_imo, sort_vecs=1, only_nonzero_vals=False)
    imp2imo = np.dot (imp2imo, evecs)
    imp2cmo = imp2imo[:,:norbs_cmo]
    imp2vmo = imp2imo[:,norbs_cmo:]
    # Sort amo in order to apply stored ci vector
    imp2gamo = np.dot (loc2imp.conjugate ().T, loc2gamo)
    amoOgamo = np.dot (imp2amo.conjugate ().T, imp2gamo)
    #print ("Overlap matrix between guess-active and active:")
    #print (prettyprint (amoOgamo, fmt='{:5.2f}'))
    Pgamo1_amo = np.einsum ('ik,jk->ijk', amoOgamo, amoOgamo.conjugate ())
    imp2ramo = np.zeros_like (imp2amo)
    ramo_evals = np.zeros (imp2ramo.shape[1], dtype=imp2ramo.dtype)
    while (Pgamo1_amo.shape[0] > 0):
        max_eval = 0
        argmax_eval = -1
        argmax_evecs = None
        for idx in range (Pgamo1_amo.shape[2]):
            evals, evecs = matrix_eigen_control_options (Pgamo1_amo[:,:,idx], sort_vecs=-1, only_nonzero_vals=False)
            if evals[0] > max_eval:
                max_eval = evals[0]
                max_evecs = evecs
                argmax_eval = idx
        #print ("With {} amos to go, assigned highest eigenvalue ({}) to {}".format (Pgamo1_amo.shape[0], max_eval, argmax_eval))
        ramo_evals[argmax_eval] = max_eval
        imp2ramo[:,argmax_eval] = np.einsum ('ij,j->i', imp2amo, max_evecs[:,0])
        imp2amo = np.dot (imp2amo, max_evecs[:,1:])
        amoOgamo = np.dot (imp2amo.conjugate ().T, imp2gamo)
        Pgamo1_amo = np.einsum ('ik,jk->ijk', amoOgamo, amoOgamo.conjugate ())
    imp2amo = imp2ramo
    print ("Fidelity of projection of guess active orbitals onto impurity space:\n{}".format (ramo_evals))
    amoOgamo = np.dot (imp2amo.conjugate ().T, imp2gamo)
    idx_signflip = np.diag (amoOgamo) < 0
    imp2amo[:,idx_signflip] *= -1
    amoOgamo = np.dot (imp2amo.conjugate ().T, imp2gamo)
    ''' 
    print ("Overlap matrix between guess-active and active:")
    print (prettyprint (amoOgamo, fmt='{:5.2f}'))
    O = np.dot (imp2amo.conjugate ().T, imp2amo) - np.eye (imp2amo.shape[1]) 
    print ("Overlap error between active and active: {}".format (linalg.norm (O)))
    O = np.dot (imp2amo.conjugate ().T, imp2cmo)    
    print ("Overlap error between active and occupied: {}".format (linalg.norm (O)))
    O = np.dot (imp2amo.conjugate ().T, imp2vmo)    
    print ("Overlap error between active and virtual: {}".format (linalg.norm (O)))
    '''
    my_occ = np.zeros (loc2imp.shape[1], dtype=np.float64)
    my_occ[:norbs_cmo] = 2
    my_occ[norbs_cmo:][:imp2amo.shape[1]] = 1
    if dm is not None:
        loc2amo = np.dot (loc2imp, imp2amo)
        evals, evecs = matrix_eigen_control_options (represent_operator_in_basis (dm, loc2amo), sort_vecs=-1, only_nonzero_vals=False)
        imp2amo = np.dot (imp2amo, evecs)
        print ("Guess density matrix eigenvalues for guess amo: {}".format (evals))
        my_occ[norbs_cmo:][:imp2amo.shape[1]] = evals
    imp2mo = np.concatenate ([imp2cmo, imp2amo, imp2vmo], axis=1)
    return imp2mo, my_occ
Esempio n. 12
0
def get_overlapping_states (bra_basis, ket_basis, across_operator=None, inner_symmetry=None, outer_symmetry=(None, None), enforce_symmetry=False,
        max_nrvecs=0, max_nlvecs=0, num_zero_atol=params.num_zero_atol, only_nonzero_vals=True, full_matrices=False):
    c2p = np.asarray (bra_basis)
    c2q = np.asarray (ket_basis)
    cOc = 1 if across_operator is None else np.asarray (across_operator)
    assert (c2p.shape[0] == c2q.shape[0]), "you need to give the two spaces in the same basis"
    assert (c2p.shape[1] <= c2p.shape[0]), "you need to give the first state in a complete basis (c2p). Did you accidentally transpose it?"
    assert (c2q.shape[1] <= c2q.shape[0]), "you need to give the second state in a complete basis (c2q). Did you accidentally transpose it?"
    assert (max_nlvecs <= c2p.shape[1]), "you can't ask for more left states than are in your left space"
    assert (max_nrvecs <= c2q.shape[1]), "you can't ask for more right states than are in your right space"
    if np.any (across_operator):
        assert (c2p.shape[0] == cOc.shape[0] and c2p.shape[0] == cOc.shape[1]), "when specifying an across_operator, it's dimensions need to be the same as the external basis"
    get_labels = (not (inner_symmetry is None)) or (not (outer_symmetry[0] is None)) or (not (outer_symmetry[1] is None))

    try:
        rets = matrix_svd_control_options (cOc, lspace=c2p, rspace=c2q, full_matrices=full_matrices,
            symmetry=inner_symmetry,
            lspace_symmetry=outer_symmetry[0],
            rspace_symmetry=outer_symmetry[1],
            strong_symm=enforce_symmetry,
            sort_vecs=-1, only_nonzero_vals=only_nonzero_vals, num_zero_atol=num_zero_atol)
        c2l, svals, c2r = rets[:3]
        if get_labels: llab, rlab = rets[3:]
    except linalg.LinAlgError as e:
        print ("LinAlgError in SVD! Analyzing...")
        if isinstance (cOc, np.ndarray):
            print ("Shape of across_operator: {}".format (cOc.shape))
            print ("Any NANs in across_operator? {}".format (np.count_nonzero (np.isnan (cOc))))
            print ("Any INFs in across_operator? {}".format (np.count_nonzero (np.isinf (cOc))))
            print ("min/max across_operator: {}/{}".format (np.amin (cOc), np.amax (cOc)))
        print ("Shape of bra_basis: {}".format (c2p.shape))
        print ("Any NANs in bra_basis? {}".format (np.count_nonzero (np.isnan (c2p))))
        print ("Any INFs in bra_basis? {}".format (np.count_nonzero (np.isinf (c2p))))
        print ("min/max bra_basis: {}/{}".format (np.amin (c2p), np.amax (c2p)))
        print ("Shape of ket_basis: {}".format (c2p.shape))
        print ("Any NANs in ket_basis? {}".format (np.count_nonzero (np.isnan (c2p))))
        print ("Any INFs in ket_basis? {}".format (np.count_nonzero (np.isinf (c2p))))
        print ("min/max ket_basis: {}/{}".format (np.amin (c2p), np.amax (c2p)))
        proj_l = c2p @ c2p.conjugate ().T
        if isinstance (cOc, np.ndarray):
            proj_l = cOc @ proj_l @ cOc
        r_symmetry = inner_symmetry if outer_symmetry[1] is None else outer_symmetry[1]
        rets = matrix_eigen_control_options (proj_l, subspace=c2q, symmetry=r_symmetry, strong_symm=enforce_symmetry, sort_vecs=-1,
            only_nonzero_vals=False, num_zero_atol=num_zero_atol)
        evals_r, c2r = rets[:2]
        if get_labels: rlab = rets[2]
        proj_r = c2q @ c2q.conjugate ().T
        if isinstance (cOc, np.ndarray):
            proj_r = cOc @ proj_r @ cOc 
        l_symmetry = inner_symmetry if outer_symmetry[0] is None else outer_symmetry[0]
        rets = matrix_eigen_control_options (proj_r, subspace=c2p, symmetry=l_symmetry, strong_symm=enforce_symmetry, sort_vecs=-1,
            only_nonzero_vals=False, num_zero_atol=num_zero_atol)
        evals_l, c2l = rets[:2]
        if get_labels: llab = rets[2]
        print ("These pairs of eigenvalues should be equal and all positive:")
        for el, er in zip (evals_l, evals_r):
            print (el, er)
        mlen = min (len (evals_l), len (evals_r))
        if len (evals_l) > mlen: print ("More left-hand eigenvalues: {}".format (evals_l[mlen:]))
        if len (evals_r) > mlen: print ("More left-hand eigenvalues: {}".format (evals_r[mlen:]))
        raise (e)

    # Truncate the basis if requested
    max_nlvecs = max_nlvecs or c2l.shape[1]
    max_nrvecs = max_nrvecs or c2r.shape[1]

    # But you can't truncate it smaller than it already is
    max_nlvecs = min (max_nlvecs, c2l.shape[1])
    max_nrvecs = min (max_nrvecs, c2r.shape[1])
    c2l = c2l[:,:max_nlvecs]
    c2r = c2r[:,:max_nrvecs]

    if get_labels: return c2l, c2r, svals, llab, rlab
    return c2l, c2r, svals