Beispiel #1
0
 def evaluate(self, r, rlY_lm):
     timer.start('tsoe eval')
     x_MM = self.zeros()
     G_LLL = gaunt(self.lmaxgaunt)
     rlY_L = rlY_lm.toarray(self.lmaxspline)
     for x_mm, oe in self.slice(x_MM):
         oe.evaluate(r, rlY_L, G_LLL, x_mm)
     timer.stop('tsoe eval')
     return x_MM
Beispiel #2
0
 def derivative(self, r, Rhat, rlY_lm, drlYdR_lmc):
     x_cMM = self.zeros((3, ))
     G_LLL = gaunt(self.lmaxgaunt)
     rlY_L = rlY_lm.toarray(self.lmaxspline)
     drlYdR_Lc = drlYdR_lmc.toarray(self.lmaxspline)
     # print(drlYdR_lmc.shape)
     for x_cmm, oe in self.slice(x_cMM):
         oe.derivative(r, Rhat, rlY_L, G_LLL, drlYdR_Lc, x_cmm)
     return x_cMM
Beispiel #3
0
    def get_magnetic_integrals(self, rgd, phi_jg, phit_jg):
        """Calculate PAW-correction matrix elements of r x nabla.

        ::

          /  _       _          _     ~   _      ~   _
          | dr [phi (r) O  phi (r) - phi (r) O  phi (r)]
          /        1     x    2         1     x    2

                       d      d
          where O  = y -- - z --
                 x     dz     dy

        and similar for y and z."""

        G_LLL = gaunt(max(self.l_j))
        Y_LLv = nabla(max(self.l_j))

        r_g = rgd.r_g
        dr_g = rgd.dr_g
        rnabla_iiv = np.zeros((self.ni, self.ni, 3))
        i1 = 0
        for j1 in range(self.nj):
            l1 = self.l_j[j1]
            nm1 = 2 * l1 + 1
            i2 = 0
            for j2 in range(self.nj):
                l2 = self.l_j[j2]
                nm2 = 2 * l2 + 1
                f1f2or = np.dot(
                    phi_jg[j1] * phi_jg[j2] - phit_jg[j1] * phit_jg[j2],
                    r_g**2 * dr_g)
                for v in range(3):
                    v1 = (v + 1) % 3
                    v2 = (v + 2) % 3
                    # term from radial wfs does not contribute
                    # term from spherical harmonics derivatives
                    G = np.zeros((nm1, nm2))
                    for l3 in range(abs(l1 - 1), l1 + 2):
                        for m3 in range(0, (2 * l3 + 1)):
                            L3 = l3**2 + m3
                            try:
                                G += np.outer(
                                    G_LLL[L3, l1**2:l1**2 + nm1, 1 + v1],
                                    Y_LLv[L3, l2**2:l2**2 + nm2, v2])
                                G -= np.outer(
                                    G_LLL[L3, l1**2:l1**2 + nm1, 1 + v2],
                                    Y_LLv[L3, l2**2:l2**2 + nm2, v1])
                            except IndexError:
                                pass  # L3 might be too large, ignore
                    rnabla_iiv[i1:i1 + nm1, i2:i2 + nm2, v] += f1f2or * G
                i2 += nm2
            i1 += nm1
        return (4 * pi / 3) * rnabla_iiv
Beispiel #4
0
 def calculate_T_Lqp(self, lcut, nq, _np, nj, jlL_i):
     G_LLL = gaunt(max(self.l_j))
     Lcut = (2 * lcut + 1)**2
     T_Lqp = np.zeros((Lcut, nq, _np))
     p = 0
     i1 = 0
     for j1, l1, L1 in jlL_i:
         for j2, l2, L2 in jlL_i[i1:]:
             if j1 < j2:
                 q = j2 + j1 * nj - j1 * (j1 + 1) // 2
             else:
                 q = j1 + j2 * nj - j2 * (j2 + 1) // 2
             T_Lqp[:, q, p] = G_LLL[L1, L2, :Lcut]
             p += 1
         i1 += 1
     return T_Lqp
Beispiel #5
0
    def get_derivative_integrals(self, rgd, phi_jg, phit_jg):
        """Calculate PAW-correction matrix elements of nabla.

        ::

          /  _       _  d       _     ~   _  d   ~   _
          | dr [phi (r) -- phi (r) - phi (r) -- phi (r)]
          /        1    dx    2         1    dx    2

        and similar for y and z."""

        G_LLL = gaunt(max(self.l_j))
        Y_LLv = nabla(max(self.l_j))

        r_g = rgd.r_g
        dr_g = rgd.dr_g
        nabla_iiv = np.empty((self.ni, self.ni, 3))
        i1 = 0
        for j1 in range(self.nj):
            l1 = self.l_j[j1]
            nm1 = 2 * l1 + 1
            i2 = 0
            for j2 in range(self.nj):
                l2 = self.l_j[j2]
                nm2 = 2 * l2 + 1
                f1f2or = np.dot(
                    phi_jg[j1] * phi_jg[j2] - phit_jg[j1] * phit_jg[j2],
                    r_g * dr_g)
                dphidr_g = np.empty_like(phi_jg[j2])
                rgd.derivative(phi_jg[j2], dphidr_g)
                dphitdr_g = np.empty_like(phit_jg[j2])
                rgd.derivative(phit_jg[j2], dphitdr_g)
                f1df2dr = np.dot(
                    phi_jg[j1] * dphidr_g - phit_jg[j1] * dphitdr_g,
                    r_g**2 * dr_g)
                for v in range(3):
                    Lv = 1 + (v + 2) % 3
                    nabla_iiv[i1:i1 + nm1, i2:i2 + nm2, v] = (
                        (4 * pi / 3)**0.5 * (f1df2dr - l2 * f1f2or) *
                        G_LLL[Lv, l2**2:l2**2 + nm2, l1**2:l1**2 + nm1].T +
                        f1f2or *
                        Y_LLv[l1**2:l1**2 + nm1, l2**2:l2**2 + nm2, v])
                i2 += nm2
            i1 += nm1
        return nabla_iiv
Beispiel #6
0
    def calculate_exx(self, s=None):
        if s is None:
            self.exx = sum(self.calculate_exx(s)
                           for s in range(self.nspins)) / self.nspins
            return self.exx

        states = []
        lmax = 0
        for ch in self.channels:
            l = ch.l
            for n, phi_g in enumerate(ch.phi_ng):
                f = ch.f_n[n]
                if f > 0 and ch.s == s:
                    states.append((l, f * self.nspins / 2.0 / (2 * l + 1),
                                   phi_g))
                    if l > lmax:
                        lmax = l

        G_LLL = gaunt(lmax)

        exx = 0.0
        j1 = 0
        for l1, f1, phi1_g in states:
            f = 1.0
            for l2, f2, phi2_g in states[j1:]:
                n_g = phi1_g * phi2_g
                for l in range((l1 + l2) % 2, l1 + l2 + 1, 2):
                    G = (G_LLL[l1**2:(l1 + 1)**2,
                               l2**2:(l2 + 1)**2,
                               l**2:(l + 1)**2]**2).sum()
                    vr_g = self.rgd.poisson(n_g, l)
                    e = f * self.rgd.integrate(vr_g * n_g, -1) / 4 / pi
                    exx -= e * G * f1 * f2
                f = 2.0
            j1 += 1

        return exx
Beispiel #7
0
def constructX(gen, gamma=0):
    """Construct the X_p^a matrix for the given atom.

    The X_p^a matrix describes the valence-core interactions of the
    partial waves.
    """
    # initialize attributes
    uv_j = gen.vu_j  # soft valence states * r:
    lv_j = gen.vl_j  # their repective l quantum numbers
    Nvi = 0
    for l in lv_j:
        Nvi += 2 * l + 1  # total number of valence states (including m)

    # number of core and valence orbitals (j only, i.e. not m-number)
    Njcore = gen.njcore
    Njval = len(lv_j)

    # core states * r:
    uc_j = gen.u_j[:Njcore]
    r, dr, N = gen.r, gen.dr, gen.N
    r2 = r**2

    # potential times radius
    vr = np.zeros(N)

    # initialize X_ii matrix
    X_ii = np.zeros((Nvi, Nvi))

    # make gaunt coeff. list
    lmax = max(gen.l_j[:Njcore] + gen.vl_j)
    G_LLL = gaunt(lmax=lmax)

    # sum over core states
    for jc in range(Njcore):
        lc = gen.l_j[jc]

        # sum over first valence state index
        i1 = 0
        for jv1 in range(Njval):
            lv1 = lv_j[jv1]

            # electron density 1 times radius times length element
            n1c = uv_j[jv1] * uc_j[jc] * dr
            n1c[1:] /= r[1:]

            # sum over second valence state index
            i2 = 0
            for jv2 in range(Njval):
                lv2 = lv_j[jv2]

                # electron density 2
                n2c = uv_j[jv2] * uc_j[jc]
                n2c[1:] /= r2[1:]

                # sum expansion in angular momenta
                for l in range(min(lv1, lv2) + lc + 1):
                    # Int density * potential * r^2 * dr:
                    if gamma == 0:
                        vr = gen.rgd.poisson(n2c, l)
                    else:
                        vr = gen.rgd.yukawa(n2c, l, gamma)
                    nv = np.dot(n1c, vr)

                    # expansion coefficients
                    A_mm = X_ii[i1:i1 + 2 * lv1 + 1, i2:i2 + 2 * lv2 + 1]
                    for mc in range(2 * lc + 1):
                        for m in range(2 * l + 1):
                            G1c = G_LLL[lv1**2:(lv1 + 1)**2, lc**2 + mc,
                                        l**2 + m]
                            G2c = G_LLL[lv2**2:(lv2 + 1)**2, lc**2 + mc,
                                        l**2 + m]
                            A_mm += nv * np.outer(G1c, G2c)
                i2 += 2 * lv2 + 1
            i1 += 2 * lv1 + 1

    # pack X_ii matrix
    X_p = pack2(X_ii)
    return X_p
Beispiel #8
0
def atomic_exact_exchange(atom, type='all'):
    """Returns the exact exchange energy of the atom defined by the
       instantiated AllElectron object 'atom'
    """
    G_LLL = gaunt(lmax=max(atom.l_j))  # Make gaunt coeff. list
    Nj = len(atom.n_j)  # The total number of orbitals

    # determine relevant states for chosen type of exchange contribution
    if type == 'all':
        nstates = mstates = range(Nj)
    else:
        Njcore = core_states(atom.symbol)  # The number of core orbitals
        if type == 'val-val':
            nstates = mstates = range(Njcore, Nj)
        elif type == 'core-core':
            nstates = mstates = range(Njcore)
        elif type == 'val-core':
            nstates = range(Njcore, Nj)
            mstates = range(Njcore)
        else:
            raise RuntimeError('Unknown type of exchange: ', type)

    # Arrays for storing the potential (times radius)
    vr = np.zeros(atom.N)
    vrl = np.zeros(atom.N)

    # do actual calculation of exchange contribution
    Exx = 0.0
    for j1 in nstates:
        # angular momentum of first state
        l1 = atom.l_j[j1]

        for j2 in mstates:
            # angular momentum of second state
            l2 = atom.l_j[j2]

            # joint occupation number
            f12 = 0.5 * (atom.f_j[j1] / (2. * l1 + 1) * atom.f_j[j2] /
                         (2. * l2 + 1))

            # electron density times radius times length element
            nrdr = atom.u_j[j1] * atom.u_j[j2] * atom.dr
            nrdr[1:] /= atom.r[1:]

            # potential times radius
            vr[:] = 0.0

            # L summation
            for l in range(l1 + l2 + 1):
                # get potential for current l-value
                hartree(l, nrdr, atom.r, vrl)

                # take all m1 m2 and m values of Gaunt matrix of the form
                # G(L1,L2,L) where L = {l,m}
                G2 = G_LLL[l1**2:(l1 + 1)**2, l2**2:(l2 + 1)**2,
                           l**2:(l + 1)**2]**2

                # add to total potential
                vr += vrl * np.sum(G2)

            # add to total exchange the contribution from current two states
            Exx += -.5 * f12 * np.dot(vr, nrdr)

    # double energy if mixed contribution
    if type == 'val-core':
        Exx *= 2.

    # return exchange energy
    return Exx
Beispiel #9
0
 def get_gaunt(self, l):
     la = self.la
     lb = self.lb
     G_LLL = gaunt(max(la, lb))
     G_mmm = G_LLL[la**2:(la + 1)**2, lb**2:(lb + 1)**2, l**2:(l + 1)**2]
     return G_mmm
Beispiel #10
0
    def calculate_exx_integrals(self):
        # Find core states:
        core = []
        lmax = 0
        for l, ch in enumerate(self.aea.channels):
            for n, phi_g in enumerate(ch.phi_ng):
                if (l >= len(self.waves_l)
                        or (l < len(self.waves_l)
                            and n + l + 1 not in self.waves_l[l].n_n)):
                    core.append((l, phi_g))
                    if l > lmax:
                        lmax = l

        lmax = max(lmax, len(self.waves_l) - 1)
        G_LLL = gaunt(lmax)

        # Calculate core contribution to EXX energy:
        self.exxcc = 0.0
        j1 = 0
        for l1, phi1_g in core:
            f = 1.0
            for l2, phi2_g in core[j1:]:
                n_g = phi1_g * phi2_g
                for l in range((l1 + l2) % 2, l1 + l2 + 1, 2):
                    G = (G_LLL[l1**2:(l1 + 1)**2, l2**2:(l2 + 1)**2,
                               l**2:(l + 1)**2]**2).sum()
                    vr_g = self.rgd.poisson(n_g, l)
                    e = f * self.rgd.integrate(vr_g * n_g, -1) / 4 / pi
                    self.exxcc -= e * G
                f = 2.0
            j1 += 1

        self.log('EXX (core-core):', self.exxcc, 'Hartree')

        # Calculate core-valence contribution to EXX energy:
        ni = sum(
            len(waves) * (2 * l + 1) for l, waves in enumerate(self.waves_l))

        self.exxcv_ii = np.zeros((ni, ni))

        i1 = 0
        for l1, waves1 in enumerate(self.waves_l):
            for phi1_g in waves1.phi_ng:
                i2 = 0
                for l2, waves2 in enumerate(self.waves_l):
                    for phi2_g in waves2.phi_ng:
                        X_mm = self.exxcv_ii[i1:i1 + 2 * l1 + 1,
                                             i2:i2 + 2 * l2 + 1]
                        if (l1 + l2) % 2 == 0:
                            for lc, phi_g in core:
                                n_g = phi1_g * phi_g
                                for l in range((l1 + lc) % 2,
                                               max(l1, l2) + lc + 1, 2):
                                    vr_g = self.rgd.poisson(phi2_g * phi_g, l)
                                    e = (self.rgd.integrate(vr_g * n_g, -1) /
                                         (4 * pi))
                                    for mc in range(2 * lc + 1):
                                        for m in range(2 * l + 1):
                                            G_L = G_LLL[:, lc**2 + mc,
                                                        l**2 + m]
                                            X_mm += np.outer(
                                                G_L[l1**2:(l1 + 1)**2],
                                                G_L[l2**2:(l2 + 1)**2]) * e
                        i2 += 2 * l2 + 1
                i1 += 2 * l1 + 1
Beispiel #11
0
def two_phi_planewave_integrals(k_Gv, setup=None, Gstart=0, Gend=None,
                                rgd=None, phi_jg=None,
                                phit_jg=None, l_j=None):
    """Calculate PAW-correction matrix elements with planewaves.

    ::
    
      /  _       _   ik.r     _     ~   _   ik.r ~   _
      | dr [phi (r) e    phi (r) - phi (r) e    phi (r)]
      /        1            2         1            2

                     ll   -  /     2                      ~       ~
      = 4pi \sum_lm i  Y (k) | dr r  [ phi (r) phi (r) - phi (r) phi (r) j (kr)
                        lm   /            1       2         1       2     ll

          /
        * | d\Omega Y     Y     Y
          /          l1m1  l2m2  lm

    """

    if Gend is None:
        Gend = len(k_Gv)
        
    if setup is not None:
        rgd = setup.rgd
        l_j = setup.l_j
        # Obtain the phi_j and phit_j
        phi_jg = []
        phit_jg = []
        rcut2 = 2 * max(setup.rcut_j)
        gcut2 = rgd.ceil(rcut2)
        for phi_g, phit_g in zip(setup.data.phi_jg, setup.data.phit_jg):
            phi_g = phi_g.copy()
            phit_g = phit_g.copy()
            phi_g[gcut2:] = phit_g[gcut2:] = 0.
            phi_jg.append(phi_g)
            phit_jg.append(phit_g)
    else:
        assert rgd is not None
        assert phi_jg is not None
        assert l_j is not None

    ng = rgd.N
    r_g = rgd.r_g
    dr_g = rgd.dr_g

    # Construct L (l**2 + m) and j (nl) index
    L_i = []
    j_i = []
    for j, l in enumerate(l_j):
        for m in range(2 * l + 1):
            L_i.append(l**2 + m)
            j_i.append(j)
    ni = len(L_i)
    nj = len(l_j)
    lmax = max(l_j) * 2 + 1

    if setup is not None:
        assert ni == setup.ni and nj == setup.nj

    # Initialize
    npw = k_Gv.shape[0]
    R_jj = np.zeros((nj, nj))
    R_ii = np.zeros((ni, ni))
    phi_Gii = np.zeros((npw, ni, ni), dtype=complex)
    j_lg = np.zeros((lmax, ng))

    # Store (phi_j1 * phi_j2 - phit_j1 * phit_j2 ) for further use
    tmp_jjg = np.zeros((nj, nj, ng))
    for j1 in range(nj):
        for j2 in range(nj):
            tmp_jjg[j1, j2] = (phi_jg[j1] * phi_jg[j2] -
                               phit_jg[j1] * phit_jg[j2])

    G_LLL = gaunt(max(l_j))
    
    # Loop over G vectors
    for iG in range(Gstart, Gend):
        kk = k_Gv[iG]
        k = np.sqrt(np.dot(kk, kk))  # calculate length of q+G

        # Calculating spherical bessel function
        for g, r in enumerate(r_g):
            j_lg[:, g] = sphj(lmax - 1, k * r)[1]

        for li in range(lmax):
            # Radial part
            for j1 in range(nj):
                for j2 in range(nj):
                    R_jj[j1, j2] = np.dot(r_g**2 * dr_g,
                                          tmp_jjg[j1, j2] * j_lg[li])

            for mi in range(2 * li + 1):
                # Angular part
                for i1 in range(ni):
                    L1 = L_i[i1]
                    j1 = j_i[i1]
                    for i2 in range(ni):
                        L2 = L_i[i2]
                        j2 = j_i[i2]
                        R_ii[i1, i2] = G_LLL[L1, L2, li**2 + mi] * R_jj[j1, j2]

                phi_Gii[iG] += R_ii * Y(li**2 + mi,
                                        kk[0] / k,
                                        kk[1] / k,
                                        kk[2] / k) * (-1j)**li
    
    phi_Gii *= 4 * pi

    return phi_Gii.reshape(npw, ni * ni)
Beispiel #12
0
def two_phi_nabla_planewave_integrals(k_Gv, setup=None, Gstart=0, Gend=None,
                                      rgd=None, phi_jg=None,
                                      phit_jg=None, l_j=None):
    """Calculate PAW-correction matrix elements with planewaves and gradient.

    ::
    
      /  _       _   ik.r d       _     ~   _   ik.r d   ~   _
      | dr [phi (r) e     -- phi (r) - phi (r) e     -- phi (r)]
      /        1          dx    2         1          dx    2

    and similar for y and z."""

    if Gend is None:
        Gend = len(k_Gv)
        
    if setup is not None:
        rgd = setup.rgd
        l_j = setup.l_j
        # Obtain the phi_j and phit_j
        phi_jg = []
        phit_jg = []
        rcut2 = 2 * max(setup.rcut_j)
        gcut2 = rgd.ceil(rcut2)
        for phi_g, phit_g in zip(setup.data.phi_jg, setup.data.phit_jg):
            phi_g = phi_g.copy()
            phit_g = phit_g.copy()
            phi_g[gcut2:] = phit_g[gcut2:] = 0.
            phi_jg.append(phi_g)
            phit_jg.append(phit_g)
    else:
        assert rgd is not None
        assert phi_jg is not None
        assert l_j is not None

    ng = rgd.N
    r_g = rgd.r_g
    dr_g = rgd.dr_g

    # Construct L (l**2 + m) and j (nl) index
    L_i = []
    j_i = []
    for j, l in enumerate(l_j):
        for m in range(2 * l + 1):
            L_i.append(l**2 + m)
            j_i.append(j)
    ni = len(L_i)
    nj = len(l_j)
    lmax = max(l_j) * 2 + 1
    
    ljdef = 3
    l2max = max(l_j) * (max(l_j) > ljdef) + ljdef * (max(l_j) <= ljdef)
    G_LLL = gaunt(l2max)
    Y_LLv = nabla(2 * l2max)

    if setup is not None:
        assert ni == setup.ni and nj == setup.nj

    # Initialize
    npw = k_Gv.shape[0]
    R1_jj = np.zeros((nj, nj))
    R2_jj = np.zeros((nj, nj))
    R_vii = np.zeros((3, ni, ni))
    phi_vGii = np.zeros((3, npw, ni, ni), dtype=complex)
    j_lg = np.zeros((lmax, ng))
    
    # Store (phi_j1 * dphidr_j2 - phit_j1 * dphitdr_j2) for further use
    tmp_jjg = np.zeros((nj, nj, ng))
    tmpder_jjg = np.zeros((nj, nj, ng))
    for j1 in range(nj):
        for j2 in range(nj):
            dphidr_g = np.empty_like(phi_jg[j2])
            rgd.derivative(phi_jg[j2], dphidr_g)
            dphitdr_g = np.empty_like(phit_jg[j2])
            rgd.derivative(phit_jg[j2], dphitdr_g)

            tmpder_jjg[j1, j2] = (phi_jg[j1] * dphidr_g -
                                  phit_jg[j1] * dphitdr_g)
            tmp_jjg[j1, j2] = (phi_jg[j1] * phi_jg[j2] -
                               phit_jg[j1] * phit_jg[j2])

    # Loop over G vectors
    for iG in range(Gstart, Gend):
        kk = k_Gv[iG]
        k = np.sqrt(np.dot(kk, kk))  # calculate length of q+G

        # Calculating spherical bessel function
        for g, r in enumerate(r_g):
            j_lg[:, g] = sphj(lmax - 1, k * r)[1]

        for li in range(lmax):
            # Radial part
            for j1 in range(nj):
                for j2 in range(nj):
                    R1_jj[j1, j2] = np.dot(r_g**2 * dr_g,
                                           tmpder_jjg[j1, j2] * j_lg[li])
                    R2_jj[j1, j2] = np.dot(r_g * dr_g,
                                           tmp_jjg[j1, j2] * j_lg[li])

            for mi in range(2 * li + 1):
                if k == 0:
                    Ytmp = Y(li**2 + mi, 1.0, 0, 0)
                else:
                    # Note the spherical bessel gives
                    # zero when k == 0 for li != 0
                    Ytmp = Y(li**2 + mi, kk[0] / k, kk[1] / k, kk[2] / k)

                for v in range(3):
                    Lv = 1 + (v + 2) % 3
                    # Angular part
                    for i1 in range(ni):
                        L1 = L_i[i1]
                        j1 = j_i[i1]
                        for i2 in range(ni):
                            L2 = L_i[i2]
                            j2 = j_i[i2]
                            l2 = l_j[j2]

                            R_vii[v, i1, i2] = ((4 * pi / 3)**0.5 *
                                                np.dot(G_LLL[L1, L2],
                                                       G_LLL[Lv, li**2 + mi])
                                                * (R1_jj[j1, j2] - l2 *
                                                   R2_jj[j1, j2]))

                            R_vii[v, i1, i2] += (R2_jj[j1, j2] *
                                                 np.dot(G_LLL[L1, li**2 + mi],
                                                        Y_LLv[:, L2, v]))

                phi_vGii[:, iG] += (R_vii * Ytmp * (-1j)**li)

    phi_vGii *= 4 * pi

    return phi_vGii.reshape(3, npw, ni * ni)
Beispiel #13
0
    def four_phi_integrals(self):
        """Calculate four-phi integral.

        Calculate the integral over the product of four all electron
        functions in the augmentation sphere, i.e.::

          /
          | d vr  ( phi_i1 phi_i2 phi_i3 phi_i4
          /         - phit_i1 phit_i2 phit_i3 phit_i4 ),

        where phi_i1 is an all electron function and phit_i1 is its
        smooth partner.
        """
        if hasattr(self, 'I4_pp'):
            return self.I4_pp

        # radial grid
        r2dr_g = self.rgd.r_g**2 * self.rgd.dr_g

        phi_jg = self.data.phi_jg
        phit_jg = self.data.phit_jg

        # compute radial parts
        nj = len(self.l_j)
        R_jjjj = np.empty((nj, nj, nj, nj))
        for j1 in range(nj):
            for j2 in range(nj):
                for j3 in range(nj):
                    for j4 in range(nj):
                        R_jjjj[j1, j2, j3, j4] = np.dot(
                            r2dr_g,
                            phi_jg[j1] * phi_jg[j2] * phi_jg[j3] * phi_jg[j4] -
                            phit_jg[j1] * phit_jg[j2] * phit_jg[j3] *
                            phit_jg[j4])

        # prepare for angular parts
        L_i = []
        j_i = []
        for j, l in enumerate(self.l_j):
            for m in range(2 * l + 1):
                L_i.append(l**2 + m)
                j_i.append(j)
        ni = len(L_i)
        # j_i is the list of j values
        # L_i is the list of L (=l**2+m for 0<=m<2*l+1) values
        # https://wiki.fysik.dtu.dk/gpaw/devel/overview.html

        G_LLL = gaunt(max(self.l_j))

        # calculate the integrals
        _np = ni * (ni + 1) // 2  # length for packing
        self.I4_pp = np.empty((_np, _np))
        p1 = 0
        for i1 in range(ni):
            L1 = L_i[i1]
            j1 = j_i[i1]
            for i2 in range(i1, ni):
                L2 = L_i[i2]
                j2 = j_i[i2]
                p2 = 0
                for i3 in range(ni):
                    L3 = L_i[i3]
                    j3 = j_i[i3]
                    for i4 in range(i3, ni):
                        L4 = L_i[i4]
                        j4 = j_i[i4]
                        self.I4_pp[p1, p2] = (
                            np.dot(G_LLL[L1, L2], G_LLL[L3, L4]) *
                            R_jjjj[j1, j2, j3, j4])
                        p2 += 1
                p1 += 1

        # To unpack into I4_iip do:
        # from gpaw.utilities import unpack
        # I4_iip = np.empty((ni, ni, _np)):
        # for p in range(_np):
        #     I4_iip[..., p] = unpack(I4_pp[:, p])

        return self.I4_pp
Beispiel #14
0
    def __init__(
        self,
        w_jg,  # all-lectron partial waves
        wt_jg,  # pseudo partial waves
        nc_g,  # core density
        nct_g,  # smooth core density
        rgd,  # radial grid descriptor
        jl,  # ?
        lmax,  # maximal angular momentum to consider
        e_xc0,  # xc energy of reference atom
        phicorehole_g,  # ?
        fcorehole,  # ?
        tauc_g,  # kinetic core energy array
        tauct_g):  # pseudo kinetic core energy array
        self.e_xc0 = e_xc0
        self.Lmax = (lmax + 1)**2
        self.rgd = rgd
        self.dv_g = rgd.dv_g
        self.Y_nL = Y_nL[:, :self.Lmax]
        self.rnablaY_nLv = rnablaY_nLv[:, :self.Lmax]
        self.ng = ng = len(nc_g)
        self.phi_jg = w_jg
        self.phit_jg = wt_jg

        self.jlL = [(j, l, l**2 + m) for j, l in jl for m in range(2 * l + 1)]
        self.ni = ni = len(self.jlL)
        self.nj = nj = len(jl)
        self.nii = nii = ni * (ni + 1) // 2
        njj = nj * (nj + 1) // 2

        self.tauc_g = tauc_g
        self.tauct_g = tauct_g
        self.tau_npg = None
        self.taut_npg = None

        G_LLL = gaunt(max(l for j, l in jl))
        B_Lqp = np.zeros((self.Lmax, njj, nii))
        p = 0
        i1 = 0
        for j1, l1, L1 in self.jlL:
            for j2, l2, L2 in self.jlL[i1:]:
                if j1 < j2:
                    q = j2 + j1 * nj - j1 * (j1 + 1) // 2
                else:
                    q = j1 + j2 * nj - j2 * (j2 + 1) // 2
                B_Lqp[:, q, p] = G_LLL[L1, L2, :self.Lmax]
                p += 1
            i1 += 1
        self.B_pqL = B_Lqp.T.copy()

        #
        self.n_qg = np.zeros((njj, ng))
        self.nt_qg = np.zeros((njj, ng))
        q = 0
        for j1, l1 in jl:
            for j2, l2 in jl[j1:]:
                self.n_qg[q] = w_jg[j1] * w_jg[j2]
                self.nt_qg[q] = wt_jg[j1] * wt_jg[j2]
                q += 1

        self.nc_g = nc_g
        self.nct_g = nct_g

        if fcorehole != 0.0:
            self.nc_corehole_g = fcorehole * phicorehole_g**2 / (4 * pi)
        else:
            self.nc_corehole_g = None