def continuum_normalization(continuum_orbitals, energy, rmax=300.0, npts_r=60, lebedev_order=23): """ continuum orbitals are normalized such that asymptotically their radial part behaves as S(r) = 1/r sin(k*r + delta) To find the factor by which the unnormalized continuum orbitals have to be multiplied, we integrate the modulus squared over one wavelength 2*pi/k starting at a large radius R: /pi /2 pi /R+2pi/k 2 2 I(wfn) = | sin(th) dth | dph | dr r |wfn(r,th,ph)| /0 /0 /R and compare with the corresponding integral for the asymptotic form S(r) /R+2pi/k 2 I(S) = | dr sin (k*r+delta) = pi / k /R The wavefunction is then normalized as wfn(normalized) = sqrt(I(S)/I(wfn)) * wfn(unnormalized) Parameters ---------- continuum_orbitals: list of unnormalized continuum orbitals, instances of LinearCombinationWavefunction or AtomicBasisFunction energy : energy of continuum orbital (in Hartree), E=1/2 k^2 Optional -------- rmax : radius (in bohr) at which the integration states npts_r : number of points in quadrature rule for radial integration on interval [rmax,rmax+2pi/k] lebedev_order : order of Lebedev grid for angular integral Returns ------- scaling_factors : list of scaling factors sqrt(I(S)/I(wfn)) """ # wavenumber k = np.sqrt(2 * energy) # angular grid Lmax, (th, ph, angular_weights) = select_angular_grid(lebedev_order) sc = np.sin(th) * np.cos(ph) ss = np.sin(th) * np.sin(ph) c = np.cos(th) # radial grid # sample points and weights for Gauss-Legendre quadrature on the interval [-1,1] leggauss_pts, leggauss_weights = legendre.leggauss(npts_r) # For Gaussian quadrature the integral [-1,1] has to be changed into [rmax,rmax+2pi/k]. # new endpoints of interval a = rmax b = rmax + 2.0 * np.pi / k # transformed sampling points and weights r = 0.5 * (b - a) * leggauss_pts + 0.5 * (a + b) radial_weights = 0.5 * (b - a) * leggauss_weights # cartesian coordinates of grid x = outerN(r, sc).flatten() y = outerN(r, ss).flatten() z = outerN(r, c).flatten() r2 = x * x + y * y + z * z weights = outerN(radial_weights, 4.0 * np.pi * angular_weights).flatten() # desired integral I(S) for asymptotic orbital Is = np.pi / k # For each continuum orbital, we compute I(wfn) and the scaling factor scaling_factors = np.zeros(len(continuum_orbitals)) for iw, wfn in enumerate(continuum_orbitals): # integral I(wfn) I = np.sum(r2 * abs(wfn.amp(x, y, z))**2 * weights) # scaling factor scaling_factors[iw] = np.sqrt(Is / I) return scaling_factors
def continuum_overlap(bfs, E, rlarge=50000.0): """ Continuum orbitals cannot be normalized, but a quantity equivalent to the overlap matrix may be defined by integrating the product of two continuum wavefunctions over one wavelength at a large distance R away from the origin /R+lambda /2pi /pi 2 * S = | dr | dph | sin(th) dth r phi (r,th,ph) phi (r,th,ph) ij /R / /0 i j where lambda = 2pi/k = 2 pi / sqrt(2*E) is the wavelength Parameters ---------- bfs : list of callables, continuum basis functions with energy E E : float, energy E in Hartree Optional -------- rlarge : float, radius R in bohr where the integration starts """ # number of basis functions Nbfs = len(bfs) # 'overlap' matrix S = np.zeros((Nbfs, Nbfs)) # angular grid Lmax, (th, ph, angular_weights) = select_angular_grid(settings.lebedev_order) Nang = len(th) sc = np.sin(th) * np.cos(ph) ss = np.sin(th) * np.sin(ph) c = np.cos(th) # wavelength k = np.sqrt(2 * E) lamb = 2.0 * np.pi / k print("k = %s" % k) # radial grid on the interval [rlarge,rlarge+lambda] Nr = 1000 r = np.linspace(rlarge, rlarge + lamb, Nr) # r^2 dr radial_weights = r**2 * np.ediff1d(r, to_end=r[-1] - r[-2]) # cartesian coordinates of grid x = outerN(r, sc).flatten() y = outerN(r, ss).flatten() z = outerN(r, c).flatten() weights = outerN(radial_weights, 4.0 * np.pi * angular_weights).flatten() # Npts = Nr * Nang for i in range(0, Nbfs): phi_i = bfs[i] for j in range(0, Nbfs): phi_j = bfs[j] # If the continuum orbitals are normalized such that the # radial functions have the form # r->oo # phi(r) --------> 1/r sin(kr + delta_l) # # then the diagonal elements of the overlap matrix # should be # S[i,i] = pi/k S[i, j] = np.sum(phi_i(x, y, z) * phi_j(x, y, z) * weights) print("overlap (unnormalized)") print(S) # For printing we normalize the overlap matrix # so that the diagonal elements are equal to 1. print("overlap (normalized)") print((S * k / np.pi)) return S
def fvvm_matrix_elements(atomlist, bfs, potential, E, r0): """ compute matrix elements of eqns. (8) and (9) in Ref. [1] which are needed for the finite-volume variational method. The integration is limited to a sphere of radius r0 around the origin. Parameters ========== atomlist : list of tuples (Z,[x,y,z]) with atomic numbers and positions bfs : list of callables, `n` real basis functions potential : callable, potential(x,y,z) evaluates effective potential E : float, energy of continuum orbital r0 : float, radius of sphere which defines the integration volume Returns ======= A : n x n matrix D : n x n matrix S : n x n matrix with overlap in finite volume """ # number of basis functions nbfs = len(bfs) print " %d basis functions" % nbfs # volume integral # __ __ # A = <\/chi |\/chi > + 2 <chi |(V-E)|chi > # ij i j V i j V A = np.zeros((nbfs, nbfs)) # # surface integral # D = <chi |chi > # ij i j surface of V D = np.zeros((nbfs, nbfs)) # # overlap inside volume # # S = <chi |chi > # ij i j V S = np.zeros((nbfs, nbfs)) # angular grid for surface integral Lmax, (th, ph, angular_weights) = select_angular_grid(settings.lebedev_order) # cartesian coordinates of Lebedev points on sphere of # radius r0 xS = r0 * np.sin(th) * np.cos(ph) yS = r0 * np.sin(th) * np.sin(ph) zS = r0 * np.cos(th) for i in range(0, nbfs): # basis function i chi_i = bfs[i] for j in range(i, nbfs): # basis function j print " integrals between basis functions i= %d and j= %d" % ( i + 1, j + 1) chi_j = bfs[j] # volume integral A_ij grad_prod_func = gradient_product(atomlist, chi_i, chi_j) def integrand(x, y, z): # (grad chi_i).(grad chi_j) integ = grad_prod_func(x, y, z) # 2 * chi_i * chi_j (V-E) integ += 2 * chi_i(x, y, z) * chi_j( x, y, z) * (potential(x, y, z) - E) # Outside the integration volume the integrand # is set to zero. r = np.sqrt(x * x + y * y + z * z) integ[r0 < r] = 0.0 return integ A[i, j] = integral(atomlist, integrand) # A is symmetric A[j, i] = A[i, j] # surface integral D_ij D[i, j] = 4.0 * np.pi * r0**2 * np.sum( angular_weights * chi_i(xS, yS, zS) * chi_j(xS, yS, zS)) # D is symmetric D[j, i] = D[i, j] # overlap integral S_ij inside the volume def integrand(x, y, z): integ = chi_i(x, y, z) * chi_j(x, y, z) # Outside the integration volume the integrand # is set to zero. r = np.sqrt(x * x + y * y + z * z) integ[r0 < r] = 0.0 return integ S[i, j] = integral(atomlist, integrand) # S is symmetric S[j, i] = S[i, j] return A, D, S