def gradient_product(atomlist, f1, f2): """ create function for evaluating the product of the gradients of two functions __ __ (\/ f1).(\/ f2) Parameters ========== atomlist : list of tuples (Zat,[xi,yi,zi]) with atomic numbers and nuclear coordinates f1,f2 : callable, f1(x,y,z) and f2(x,y,z) should evaluate the functions at the grid points specified by x = [x0,x1,...,xn], y = [y0,y1,...yn] and z = [z0,z1,...,zn] Returns ======= grad_prod_func : callable, evaluates (grad f1).(grad f2)(x,y,z) """ # Translate molecular geometry into the format understood # by the multicenter integration code atomic_numbers, atomic_coordinates = atomlist2arrays(atomlist) grad_prod_func = multicenter_gradient_product( f1, f2, atomic_coordinates, atomic_numbers, cusps_separate=True, radial_grid_factor=settings.radial_grid_factor, lebedev_order=settings.lebedev_order) return grad_prod_func
def electronic_dipole(atomlist, bfA, bfB): """ electric dipole between two basis functions (a|e*r|b) Parameters ========== atomlist : list of tuples (Z,[x,y,z]) with atomic numbers and positions bfA, bfB : instances of AtomicBasisFunction or callables, e.g. bfA(x,y,z) etc. Returns ======= Dab : numpy array with components [Dx,Dy,Dz] of dipole matrix elements (a|e*x|b), (a|e*y|b), (a|e*z|b) """ # Bring geometry data into a form understood by the module MolecularIntegrals atomic_numbers, atomic_coordinates = atomlist2arrays(atomlist) # 1. define integrand, xyz=0,1,2 selects the x,y or z-component def dipole_density(x,y,z, xyz=0): er = [x,y,z] return bfA(x,y,z) * er[xyz] * bfB(x,y,z) # # 2. integrate density on a multicenter grid Dab = np.zeros(3, dtype=complex) for xyz in [0,1,2]: Dab[xyz] = multicenter_integration(lambda x,y,z: dipole_density(x,y,z, xyz=xyz), atomic_coordinates, atomic_numbers, radial_grid_factor=settings.radial_grid_factor, lebedev_order=settings.lebedev_order) return Dab.real
def overlap(atomlist, bfA, bfB): """ overlap between two basis functions (a|b) Parameters ========== atomlist : list of tuples (Z,[x,y,z]) with atomic numbers and positions bfA, bfB : instances of AtomicBasisFunction or callables, e.g. bfA(x,y,z) etc. Returns ======= Sab : float, overlap integral """ # Bring geometry data into a form understood by the module MolecularIntegrals atomic_numbers, atomic_coordinates = atomlist2arrays(atomlist) # Now we compute the integrals numerically on a multicenter grid. # 1. define integrand s = a b def s_integrand(x,y,z): return bfA(x,y,z).conjugate() * bfB(x,y,z) # # 2. integrate density on a multicenter grid Sab = multicenter_integration(s_integrand, atomic_coordinates, atomic_numbers, radial_grid_factor=settings.radial_grid_factor, lebedev_order=settings.lebedev_order) return Sab
def integral(atomlist, f): """ compute the integral / | f(r) dV / Parameters ========== atomlist : list of tuples (Z,[x,y,z]) with atomic numbers and positions f : callable, f(x,y,z) Returns ======= integ : float, volume integral """ # Bring geometry data into a form understood by the module MolecularIntegrals atomic_numbers, atomic_coordinates = atomlist2arrays(atomlist) # compute integral on a multicenter grid integ = multicenter_integration(f, atomic_coordinates, atomic_numbers, radial_grid_factor=settings.radial_grid_factor, lebedev_order=settings.lebedev_order) return integ
def nuclear_repulsion(atomlist): """ repulsion energy between nuclei Z(i) Z(j) V_rep = sum sum ------------- i j>i |R(i)-R(j)| Parameters ========== atomlist : list of tuples (Z,[x,y,z]) for each atom, molecular geometry Returns ======= Vrep : float, repulsion energy """ atomic_numbers, atomic_coordinates = atomlist2arrays(atomlist) Nat = len(atomic_numbers) Vrep = 0.0 for i in range(0, Nat): Zi = atomic_numbers[i] Ri = atomic_coordinates[:,i] for j in range(i+1,Nat): Zj = atomic_numbers[j] Rj = atomic_coordinates[:,j] Vrep += Zi*Zj / la.norm(Ri-Rj) return Vrep
def linear_combination(atomlist, orbitals, coeffs): norb = len(orbitals) assert len(orbitals) == len(coeffs) atomic_numbers, atomic_coordinates = atomlist2arrays(atomlist) lcao_orb = multicenter_operation(orbitals, lambda orbs: sum( [coeffs[j] * orbs[j] for j in range(0, norb)] ), atomic_coordinates, atomic_numbers, radial_grid_factor=settings.radial_grid_factor, lebedev_order=settings.lebedev_order) return lcao_orb
def nuclear(atomlist, bfA, bfB): """ matrix element of nuclear attraction (-Zk) (a| sum_k -------- |b) |r-Rk| Parameters ========== atomlist : list of tuples (Z,[x,y,z]) with atomic numbers and positions bfA, bfB : instances of AtomicBasisFunction or callables, e.g. bfA(x,y,z) etc. Returns ======= Nab : float, integral of nuclear attraction """ # Bring data into a form understood by the module MolecularIntegrals atomic_numbers, atomic_coordinates = atomlist2arrays(atomlist) # Now we compute the integrals numerically on a multicenter grid. # # 1. define integrand n = sum_k (-Zk)/|r-Rk| * a(r) * b(r) def nuc_integrand(x, y, z): # nuclear attraction potential nuc = 0.0 * x for Zk, Rk in atomlist: # position of k-th nucleus X, Y, Z = Rk # Zk is the atomic number of k-th nucleus nuc -= Zk / np.sqrt((x - X)**2 + (y - Y)**2 + (z - Z)**2) # product of bra and ket wavefunctions rhoAB = bfA(x, y, z) * bfB(x, y, z) nuc *= rhoAB return nuc # # 2. integrate nuclear attraction energy density on multicenter grid Nab = multicenter_integration( nuc_integrand, atomic_coordinates, atomic_numbers, radial_grid_factor=settings.radial_grid_factor, lebedev_order=settings.lebedev_order) return Nab
def kinetic(atomlist, bfA, bfB): """ matrix element of kinetic energy __2 (a| - 1/2 \/ |b) Parameters ========== atomlist : list of tuples (Z,[x,y,z]) with atomic numbers and positions bfA, bfB : instances of AtomicBasisFunction or callables, e.g. bfA(x,y,z) etc. Returns ======= Tab : float, kinetic energy integral """ # Bring data into a form understood by the module MolecularIntegrals atomic_numbers, atomic_coordinates = atomlist2arrays(atomlist) # Now we compute the integrals numerically on a multicenter grid. # __2 # 1. find Laplacian \/ b lapB = multicenter_laplacian( bfB, atomic_coordinates, atomic_numbers, radial_grid_factor=settings.radial_grid_factor, lebedev_order=settings.lebedev_order) # __2 # 2. define integrand t = -1/2 a \/ b def t_integrand(x, y, z): return -0.5 * bfA(x, y, z) * lapB(x, y, z) # # 3. integrate kinetic energy density on multicenter grid Tab = multicenter_integration( t_integrand, atomic_coordinates, atomic_numbers, radial_grid_factor=settings.radial_grid_factor, lebedev_order=settings.lebedev_order) return Tab
# define density function def rho(x, y, z): r = np.sqrt(x * x + y * y + z * z) n = np.exp(-alpha * r) # normalization constant such that # / # N | exp(-alpha*r) dV = 1 # / nrm = alpha**3 / (8.0 * np.pi) return nrm * n print "(grad rho)^2 ..." # GGA xc-functionals need the gradient squared of the density # __ 2 # sigma = (\/ rho) atomic_numbers, atomic_coordinates = atomlist2arrays(atomlist) sigma = multicenter_gradient2( rho, atomic_coordinates, atomic_numbers, radial_grid_factor=settings.radial_grid_factor, lebedev_order=settings.lebedev_order) s = sigma(r, 0 * r, 0 * r) n = rho(r, 0 * r, 0 * r) # dimensional quantity x = |grad rho|/rho^{4/3} x = np.sqrt(s) / n**(4.0 / 3.0) # van Leeuven and Baerends correction to exchange-correction # potential with asymptotic -1/r behaviour. The functional
def lithium_cation_continuum(l, m, k): """ compute continuum orbital in the electrostatic potential of the Li^+ core Parameters ---------- l,m : angular quantum numbers of asymptotic solution e.g. l=0,m=0 s-orbital l=1,m=+1 px-orbital k : length of wavevector in a.u., the energy of the continuum orbital is E=1/2 k^2 """ # Li^+ atom atomlist = [(3, (0.0, 0.0, 0.0))] charge = +1 # choose resolution of multicenter grids for bound orbitals settings.radial_grid_factor = 20 # controls size of radial grid settings.lebedev_order = 25 # controls size of angular grid # 1s core orbitals for Li+^ atom RDFT = BasissetFreeDFT(atomlist, None, charge=charge) # bound_orbitals = RDFT.getOrbitalGuess() Etot, bound_orbitals, orbital_energies = RDFT.solveKohnSham() # choose resolution of multicenter grids for continuum orbitals settings.radial_grid_factor = 120 # controls size of radial grid settings.lebedev_order = 41 # controls size of angular grid # show number of radial and angular points in multicenter grid print_grid_summary(atomlist, settings.lebedev_order, settings.radial_grid_factor) print "electron density..." # electron density of two electrons in the 1s core orbital rho = density_func(bound_orbitals) print "effective potential..." # potential energy for Li nucleus and 2 core electrons potential = effective_potential_func(atomlist, rho, None, nelec=2) def v0(x, y, z): r = np.sqrt(x * x + y * y + z * z) return -1.0 / r def v1(x, y, z): return potential(x, y, z) - v0(x, y, z) # The continuum orbital is specified by its energy and asymptotic # angular momentum (E,l,m) # energy of continuum orbital E = 0.5 * k**2 # angular quantum numbers of asymptotic solution assert abs(m) <= l print " " print "Asymptotic continuum wavefunction" print "=================================" print " energy E= %e Hartree ( %e eV )" % ( E, E * AtomicData.hartree_to_eV) print " wavevector k= %e a.u." % k print " angular moment l= %d m= %+d" % (l, m) print " " # asymptotically correct solution for V0 = -1/r (hydrogen) Cf = regular_coulomb_func(E, charge, l, m, 0.0) phi0 = Cf # right-hand side of inhomogeneous Schroedinger equation def source(x, y, z): return -v1(x, y, z) * phi0(x, y, z) # # solve (H0 + V1 - E) dphi = - V1 phi0 # for orbital correction dphi atomic_numbers, atomic_coordinates = atomlist2arrays(atomlist) print "Schroedinger equation..." dphi = multicenter_inhomogeneous_schroedinger( potential, source, E, atomic_coordinates, atomic_numbers, radial_grid_factor=settings.radial_grid_factor, lebedev_order=settings.lebedev_order) # Combine asymptotically correct solution with correction # phi = phi0 + dphi phi = add_two_functions(atomlist, phi0, dphi, 1.0, 1.0) # residual for phi0 R0 = (H-E)phi0 residual0 = residual_func(atomlist, phi0, potential, E) # residual for final solution R = (H-E)phi residual = residual_func(atomlist, phi, potential, E) # spherical average of residual function residual_avg = spherical_average_residual_func(atomlist[0], residual) # The phase shift is determined by matching the radial wavefunction # to a shifted and scaled Coulomb function at a number of radial # sampling points drawn from the interval [rmin, rmax]. # On the one hand rmin < rmax should be chosen large enough, # so that the continuum orbital approaches its asymptotic form, # on the other hand rmax should be small enough that the accuracy # of the solution due to the sparse r-grid is still high enough. # A compromise has to be struck depending on the size of the radial grid. # The matching points are spread over several periods, # but not more than 30 bohr. wavelength = 2.0 * np.pi / k print "wavelength = %e" % wavelength rmin = 70.0 rmax = rmin + max(10 * wavelength, 30.0) Npts = 100 # determine phase shift and scaling factor by a least square # fit the the regular Coulomb function scale, delta = phaseshift_lstsq(atomlist, phi, E, charge, l, m, rmin, rmax, Npts) print "scale factor (relative to Coulomb wave) = %s" % scale print "phase shift (relative to Coulomb wave) = %e " % delta # normalize wavefunction, so that 1/scale phi(x,y,z) approaches # asymptotically a phase-shifted Coulomb wave phi_norm = multicenter_operation( [phi], lambda fs: fs[0] / scale, atomic_coordinates, atomic_numbers, radial_grid_factor=settings.radial_grid_factor, lebedev_order=settings.lebedev_order) # The continuum orbital should be orthogonal to the bound # orbitals belonging to the same Hamiltonian. I think this # should come out correctly by default. print " " print " Overlaps between bound orbitals and continuum orbital" print " =====================================================" for ib, bound_orbital in enumerate(bound_orbitals): olap_bc = overlap(atomlist, bound_orbital, phi_norm) print " <bound %d| continuum> = %e" % (ib + 1, olap_bc) print "" # shifted regular coulomb function Cf_shift = regular_coulomb_func(E, charge, l, m, delta) # save radial wavefunctions and spherically averaged residual # radial part of Coulomb wave without phase shift phi0_rad = radial_wave_func(atomlist, phi0, l, m) # radial part of shifted Coulomb wave Cf_shift_rad = radial_wave_func(atomlist, Cf_shift, l, m) # radial part of scattering solution phi_norm_rad = radial_wave_func(atomlist, phi_norm, l, m) print "" print "# RADIAL_WAVEFUNCTIONS" print "# Asymptotic wavefunction:" print "# charge Z= %+d" % charge print "# energy E= %e k= %e" % (E, k) print "# angular momentum l= %d m= %+d" % (l, m) print "# phase shift delta= %e rad" % delta print "# " print "# R/bohr Coulomb Coulomb radial wavefunction spherical avg. residual" print "# shifted R_{l,m}(r) <|(H-E)phi|^2>" import sys # write table to console r = np.linspace(1.0e-3, 100, 1000) data = np.vstack((r, phi0_rad(r), Cf_shift_rad(r), phi_rad(r), residual_avg(r))).transpose() np.savetxt(sys.stdout, data, fmt=" %+e ") print "# END" print ""
def electron_repulsion_integral_rho(atomlist, rhoAB, rhoCD, rho, xc): """ compute the electron repulsion integral (ab|{1/r12+f_xc[rho]}|cd) Parameters ========== atomlist : list of tuples (Z,[x,y,z]) with atomic numbers and positions rhoAB, rhoCD : callables, rhoAB(x,y,z) computes the product a(x,y,z)*b(x,y,z) and rhoCD(x,y,z) computes the product c(x,y,z)*d(x,y,z) rho : callable, rho(x,y,z) computes the total electron density xc : instance of libXCFunctional, only LDA functional will work properly as the density gradient is not available Returns ======= Iabcd : float, sum of Hartree and exchange correlation part of electron integral """ # bring data into a form understood by the module MolecularIntegrals atomic_numbers, atomic_coordinates = atomlist2arrays(atomlist) # Now we compute the integrals numerically on a multicenter grid. # # compute electrostatic Hartree term # (ab|1/r12|cd) # 1. solve the Poisson equation to get the electrostatic potential # Vcd(r) due to the charge distribution c(r)*d(r) Vcd = multicenter_poisson(rhoCD, atomic_coordinates, atomic_numbers, radial_grid_factor=settings.radial_grid_factor, lebedev_order=settings.lebedev_order) # # 2. integrate a(r)*b(r)*Vcd(r) def Iabcd_hartree_integrand(x, y, z): return rhoAB(x, y, z) * Vcd(x, y, z) # Coulomb integral Iabcd_hartree = multicenter_integration( Iabcd_hartree_integrand, atomic_coordinates, atomic_numbers, radial_grid_factor=settings.radial_grid_factor, lebedev_order=settings.lebedev_order) # # compute contribution from exchange-correlation functional # (ab|f_xc[rho]|cd) def Iabcd_fxc_integrand(x, y, z): return rhoAB(x, y, z) * xc.fxc(rho(x, y, z)) * rhoCD(x, y, z) Iabcd_xc = multicenter_integration( Iabcd_fxc_integrand, atomic_coordinates, atomic_numbers, radial_grid_factor=settings.radial_grid_factor, lebedev_order=settings.lebedev_order) Iabcd = Iabcd_hartree + Iabcd_xc # check that density integrates to the correct number of electrons total_elec_charge = multicenter_integration( rho, atomic_coordinates, atomic_numbers, radial_grid_factor=settings.radial_grid_factor, lebedev_order=settings.lebedev_order) total_nuc_charge = sum([Zi for (Zi, posi) in atomlist]) #print "total electronic charge : %e" % total_elec_charge #print "total nuclear charge : %e" % total_nuc_charge assert abs(total_elec_charge - total_nuc_charge) < 1.0e-3 #print "Hartree contribution (ab|1/r12|cd) = %+e" % Iabcd_hartree #print "XC-contribution (ab|f_xc[rho0]|cd) = %+e" % Iabcd_xc return Iabcd