def get_1p_or_1h(self): occa, occb, vira, virb = self._get_slices() self.e_mp2 = 0 hs = [] for a, b in [(0, 1), (1, 0)]: if a == 0: occa, occb, vira, virb = self._get_slices() else: occb, occa, virb, vira = self._get_slices() eri_aa_ovov = self.eri[a, a][occa, vira, occa, vira] eri_ab_ovov = self.eri[a, b][occa, vira, occb, virb] eo_a = self.hf.e[a][occa] eo_b = self.hf.e[b][occb] ev_a = self.hf.e[a][vira] ev_b = self.hf.e[b][virb] t2_aa = eri_aa_ovov.copy() t2_aa /= util.dirsum('i,a,j,b->iajb', eo_a, -ev_a, eo_a, -ev_a) t2a_aa = t2_aa - t2_aa.swapaxes(0, 2).copy() t2_ab = eri_ab_ovov.copy() t2_ab /= util.dirsum('i,a,j,b->iajb', eo_a, -ev_a, eo_b, -ev_b) self.e_mp2 += util.einsum('iajb,iajb->', t2a_aa, eri_aa_ovov) * 0.5 self.e_mp2 += util.einsum('iajb,iajb->', t2_ab, eri_ab_ovov) * 0.5 h = np.diag(eo_a) h += util.einsum('a,iakb,jakb->ij', ev_a, t2a_aa, t2a_aa) h += util.einsum('a,iakb,jakb->ij', ev_a, t2_ab, t2_ab) h += util.einsum('a,ibka,jbka->ij', ev_b, t2_ab, t2_ab) h -= util.einsum('k,iakb,jakb->ij', eo_a, t2a_aa, t2a_aa) * 0.5 h -= util.einsum('k,iakb,jakb->ij', eo_b, t2_ab, t2_ab) * 0.5 h -= util.einsum('k,iakb,jakb->ij', eo_b, t2_ab, t2_ab) * 0.5 h -= util.einsum('i,iakb,jakb->ij', eo_a, t2a_aa, t2a_aa) * 0.25 h -= util.einsum('i,iakb,jakb->ij', eo_a, t2_ab, t2_ab) * 0.25 h -= util.einsum('i,iakb,jakb->ij', eo_a, t2_ab, t2_ab) * 0.25 h -= util.einsum('j,iakb,jakb->ij', eo_a, t2a_aa, t2a_aa) * 0.25 h -= util.einsum('j,iakb,jakb->ij', eo_a, t2_ab, t2_ab) * 0.25 h -= util.einsum('j,iakb,jakb->ij', eo_a, t2_ab, t2_ab) * 0.25 h += util.einsum('iakb,jakb->ij', t2a_aa, eri_aa_ovov) * 0.5 h -= util.einsum('iakb,jbka->ij', t2a_aa, eri_aa_ovov) * 0.5 h += util.einsum('iakb,jakb->ij', t2_ab, eri_ab_ovov) h += util.einsum('jakb,iakb->ij', t2a_aa, eri_aa_ovov) * 0.5 h -= util.einsum('jakb,kaib->ij', t2a_aa, eri_aa_ovov) * 0.5 h += util.einsum('jakb,iakb->ij', t2_ab, eri_ab_ovov) hs.append(h) self.h_1p_or_1h = tuple(hs)
def build_se(self): if self.rpa is None: self.solve_casida() e_rpa, v_rpa, xpy = self.rpa naux_gf = self.gf.naux chempot = self.chempot c = self.gf.v[:self.nphys] co = c[:, self.gf.e < chempot] cv = c[:, self.gf.e >= chempot] xyia = util.mo2qo(self.eri, c, co, cv).reshape(self.nphys, naux_gf, -1) omega = util.einsum('xyk,ks->xys', xyia, xpy) e_gf = self.gf.e e_rpa_s = np.outer(np.sign(e_gf - chempot), e_rpa) e = util.dirsum('i,ij->ij', e_gf, e_rpa_s).flatten() v = omega.reshape((self.nphys, -1)) self.se = aux.Aux(e, v, chempot=self.chempot) log.write('naux (se,build) = %d\n' % self.se.naux, self.verbose) return self.se
def gradient(x, se, fock, nelec, occupancy=2.0, buf=None, return_val=False): ''' Gradient function for the minimization with respect to error. Parameters ---------- x : float objective parameter, i.e. chemical potential on the auxiliary space buf : ndarray, optional array to store the extended Fock matrix occupancy : float, optional orbital occupancy, i.e. 2 for RHF and 1 for UHF, default 2 Returns ------- fx : float objective function value, i.e. error in number of electrons dx : float gradient value ''' w, v, chempot, error = diag_fock_ext(se, fock, nelec, chempot=x, buf=buf, occupancy=occupancy) nocc = np.sum(w < chempot) phys = slice(None, se.nphys) aux = slice(se.nphys, None) occ = slice(None, nocc) vir = slice(nocc, None) h_1 = -np.dot(v[aux, vir].conj().T, v[aux, occ]) z_ai = -h_1 / util.dirsum('i,a->ai', w[occ], -w[vir]) c_occ = np.dot(v[phys, vir], z_ai) rdm1 = np.dot(v[phys, occ], c_occ.conj().T) * 4 ne = np.trace(rdm1).real d = 2 * error * ne if return_val: return error**2, d else: return d
def get_1p_or_1h(self): occ, vir = self._get_slices() eri_ovov = self.eri[occ, vir, occ, vir] eo = self.hf.e[occ] ev = self.hf.e[vir] t2 = eri_ovov.copy() t2 /= util.dirsum('i,a,j,b->iajb', eo, -ev, eo, -ev) t2a = t2 - t2.swapaxes(0, 2).copy() os = self.options['os_factor'] ss = self.options['ss_factor'] self.e_mp2 = util.einsum('iajb,iajb->', t2, eri_ovov) * (os + ss) self.e_mp2 -= util.einsum('iajb,ibja->', t2, eri_ovov) * ss h = np.diag(eo) h += util.einsum('a,iakb,jakb->ij', ev, t2a, t2a) h += util.einsum('a,iakb,jakb->ij', ev, t2, t2) h += util.einsum('a,ibka,jbka->ij', ev, t2, t2) h -= util.einsum('k,iakb,jakb->ij', eo, t2a, t2a) * 0.5 h -= util.einsum('k,iakb,jakb->ij', eo, t2, t2) * 0.5 h -= util.einsum('k,iakb,jakb->ij', eo, t2, t2) * 0.5 h -= util.einsum('i,iakb,jakb->ij', eo, t2a, t2a) * 0.25 h -= util.einsum('i,iakb,jakb->ij', eo, t2, t2) * 0.25 h -= util.einsum('i,iakb,jakb->ij', eo, t2, t2) * 0.25 h -= util.einsum('j,iakb,jakb->ij', eo, t2a, t2a) * 0.25 h -= util.einsum('j,iakb,jakb->ij', eo, t2, t2) * 0.25 h -= util.einsum('j,iakb,jakb->ij', eo, t2, t2) * 0.25 h += util.einsum('iakb,jakb->ij', t2a, eri_ovov) * 0.5 h -= util.einsum('iakb,jbka->ij', t2a, eri_ovov) * 0.5 h += util.einsum('iakb,jakb->ij', t2, eri_ovov) h += util.einsum('jakb,iakb->ij', t2a, eri_ovov) * 0.5 h -= util.einsum('jakb,kaib->ij', t2a, eri_ovov) * 0.5 h += util.einsum('jakb,iakb->ij', t2, eri_ovov) self.h_1p_or_1h = h
def build_dfrmp2_part(eo, ev, ixq, qja, wtol=1e-12, ss_factor=1.0, os_factor=1.0): ''' Builds a set of auxiliaries representing all (i,j,a) or (a,b,i) diagrams for a restricted reference. Parameters ---------- eo : (o) ndarray occupied (virtual) energies ev : (v) ndarray virtual (occupied) energies ixq : (o,n,q) ndarray density-fitted two-electron integrals indexed as occupied, physical, auxiliary (physical, virtual, auxiliary) qja : (q,o,v) ndarray density-fitted two-electron integrals indexed as auxiliary, occupied, virtual (auxiliary, virtual, occupied) wtol : float, optional threshold for an eigenvalue to be considered zero ss_factor : float, optional same spin factor, default 1.0 os_factor : float, optional opposite spin factor, default 1.0 Returns ------- e : (m) ndarray auxiliary energies v : (n,m) ndarray auxiliary couplings ''' nphys = ixq.shape[1] ndf, nocc, nvir = qja.shape npoles = nocc * nocc * nvir e = np.zeros((npoles), dtype=types.float64) v = np.zeros((nphys, npoles), dtype=ixq.dtype) pos_factor = np.sqrt(0.5 * os_factor) neg_factor = np.sqrt(0.5 * os_factor + ss_factor) dia_factor = np.sqrt(os_factor) ixq = ixq.reshape((nocc * nphys, ndf)) qja = qja.reshape((ndf, nocc * nvir)) n0 = 0 for i in range(nocc): nja = i * nvir am = slice(n0, n0 + nja) bm = slice(n0 + nja, n0 + nja * 2) cm = slice(n0 + nja * 2, n0 + nja * 2 + nvir) xq = ixq[i * nphys:(i + 1) * nphys] qa = qja[:, i * nvir:(i + 1) * nvir] xja = np.dot(ixq[:i * nphys].conj(), qa) xja = util.reshape_internal(xja, (i, nphys, nvir), (0, 1), (nphys, i * nvir)) xia = np.dot(xq.conj(), qja[:, :i * nvir]).reshape((nphys, -1)) xa = np.dot(xq.conj(), qa) e[am] = e[bm] = eo[i] + util.dirsum('i,a->ia', eo[:i], -ev).ravel() e[cm] = 2 * eo[i] - ev v[:, am] = neg_factor * (xja - xia) v[:, bm] = pos_factor * (xja + xia) v[:, cm] = dia_factor * xa n0 += nja * 2 + nvir mask = np.absolute(util.complex_sum(v * v, axis=0)) >= wtol e = e[mask] v = v[:, mask] assert e.shape[0] == v.shape[1] return e, v
def build_dfrmp2_part_direct(eo, ev, ixq, qja, wtol=1e-12, ss_factor=1.0, os_factor=1.0): ''' Builds a set of auxiliaries representing all (i,j,a) or (a,b,i) diagrams for a restricted reference. Uses a generator which iterates over blocks. Parameters ---------- eo : (o) ndarray occupied (virtual) energies ev : (v) ndarray virtual (occupied) energies ixq : (o,n,q) ndarray density-fitted two-electron integrals indexed as occupied, physical, auxiliary (physical, virtual, auxiliary) qja : (q,o,v) ndarray density-fitted two-electron integrals indexed as auxiliary, occupied, virtual (auxiliary, virtual, occupied) wtol : float, optional threshold for an eigenvalue to be considered zero ss_factor : float, optional same spin factor, default 1.0 os_factor : float, optional opposite spin factor, deafult 1.0 Yields ------ e : (m) ndarray auxiliary energies v : (n,m) ndarray auxiliary couplings ''' nphys = ixq.shape[1] ndf, nocc, nvir = qja.shape npoles = nocc * nocc * nvir pos_factor = np.sqrt(0.5 * os_factor) neg_factor = np.sqrt(0.5 * os_factor + ss_factor) dia_factor = np.sqrt(os_factor) ixq = ixq.reshape((nocc * nphys, ndf)) qja = qja.reshape((ndf, nocc * nvir)) for i in range(nocc): nja = i * nvir xq = ixq[i * nphys:(i + 1) * nphys] qa = qja[:, i * nvir:(i + 1) * nvir] xja = np.dot(ixq[:i * nphys].conj(), qa) xja = util.reshape_internal(xja, (i, nphys, nvir), (0, 1), (nphys, i * nvir)) xia = np.dot(xq.conj(), qja[:, :i * nvir]).reshape((nphys, -1)) xa = np.dot(xq.conj(), qa) ea = eb = eo[i] + util.dirsum('i,a->ia', eo[:i], -ev).ravel() ec = 2 * eo[i] - ev va = neg_factor * (xja - xia) vb = pos_factor * (xja + xia) vc = dia_factor * xa if len(ea): yield ea, va if len(eb): yield eb, vb if len(ec): yield ec, vc