def test(): # bond length in bohr dist = 2.0 # positions of protons posH1 = (0.0, 0.0, -dist / 2.0) posH2 = (0.0, 0.0, +dist / 2.0) atomlist = [(1, posH1), (1, posH2)] # Set resolution of multicenter grid settings.radial_grid_factor = 20 settings.lebedev_order = 23 # energy of continuum orbital E = 1.0 # same functional as used in the calculation of pseudo orbitals xc = XCFunctionals.libXCFunctional(Parameters.pseudo_orbital_x, Parameters.pseudo_orbital_c) dft = BasissetFreeDFT(atomlist, xc) print "initial orbital guess from DFTB calculation" orbitals = dft.getOrbitalGuess() norb = len(orbitals) # all orbitals are doubly occupied nelec = 2 * norb bound_orbitals = dft.getOrbitalGuess() # effective potential rho = density_func(bound_orbitals) veff = effective_potential_func(atomlist, rho, xc, nelec=nelec) # radius that separates inner from outer region r0 = vdw_sphere_radius(atomlist, fac=2.0) # basis sets bs_core = AtomicBasisSet(atomlist, orbital_set="core") bs_valence = AtomicBasisSet(atomlist, orbital_set="valence") bs_continuum = AtomicScatteringBasisSet(atomlist, E, lmax=1) # combine basis functions from all basis sets bfs = bs_core.bfs + bs_valence.bfs + bs_continuum.bfs A, D, S = fvvm_matrix_elements(atomlist, bfs, veff, E, r0) # solve generalized eigenvalue problem # A.C = b.D.C b, C = sla.eigh(A, D) return b, C
def loadContinuum(self): E = self.photo_kinetic_energy k = np.sqrt(2 * E) wavelength = 2.0 * np.pi / k # determine the radius of the sphere where the angular distribution is calculated. It should be # much larger than the extent of the molecule (xmin, xmax), (ymin, ymax), (zmin, zmax) = Cube.get_bbox(self.atomlist, dbuff=0.0) dx, dy, dz = xmax - xmin, ymax - ymin, zmax - zmin Rmax0 = self.settings.getOption("Averaging", "sphere radius Rmax") Rmax = max([dx, dy, dz]) + Rmax0 Npts = max(int(Rmax), 1) * 50 print "Radius of sphere around molecule, Rmax = %s bohr" % Rmax print "Points on radial grid, Npts = %d" % Npts self.bs_free = AtomicScatteringBasisSet(self.atomlist, E, rmin=0.0, rmax=Rmax + 2 * wavelength, Npts=Npts) self.SKT_bf, SKT_ff = load_slako_scattering(self.atomlist, E) if self.settings.getOption("Continuum Orbital", "Ionization transitions") == "inter-atomic": inter_atomic = True else: inter_atomic = False print "inter-atomic transitions: %s" % inter_atomic self.Dipole = ScatteringDipoleMatrix(self.atomlist, self.valorbs, self.SKT_bf, inter_atomic=inter_atomic).real # if self.activate_average.isChecked(): print "ORIENTATION AVERAGING" npts_euler = self.settings.getOption("Averaging", "Euler angle grid points") npts_theta = self.settings.getOption("Averaging", "polar angle grid points") self.orientation_averaging = PAD.OrientationAveraging_small_memory( self.Dipole, self.bs_free, Rmax, E, npts_euler=npts_euler, npts_theta=npts_theta) else: print "NO AVERAGING"
def test_continuum_normalization(): """ check that atomic continuum orbitals are correctly normalized """ # single atom at origin atomlist = [(1, [0, 0, 0])] # energy E=1/2 k^2 for energy in np.linspace(0.1, 0.5, 10): # atomic continuum orbitals basis = AtomicScatteringBasisSet(atomlist, energy, lmax=3) # scaling_factors = continuum_normalization(basis.bfs, energy) print "energy= %s" % energy print "scaling factors" print scaling_factors # Since atomic continuum orbitals are correctly normalized by # default, all scaling factors should be 1 assert la.norm(abs(scaling_factors - 1.0)) < 1.0e-2
def test_mesh_size(): """ check that the integrals <i|(H-E)|j> between scattering orbitals do not depend on the mesh size, i.e. are converged for a the grid sizes used """ import matplotlib.pyplot as plt # H2 atomlist = [(1, [0, 0, -0.7]), (1, [0, 0, +0.7])] E = 0.1 basis = AtomicScatteringBasisSet(atomlist, E) pot0 = AtomicPotentialSuperposition(atomlist, confined=False) for radial_grid_factor in [1, 2, 3, 4]: for lebedev_order in [11, 17, 23]: L = scattering_integrals(atomlist, basis.bfs, pot0, radial_grid_factor=radial_grid_factor, lebedev_order=lebedev_order) import matplotlib.pyplot as plt plt.imshow(L) plt.show()
def averaged_pad_scan(xyz_file, dyson_file, selected_orbitals, npts_euler, npts_theta, nskip, inter_atomic, sphere_radius): molecule_name = os.path.basename(xyz_file).replace(".xyz", "") atomlist = XYZ.read_xyz(xyz_file)[-1] # shift molecule to center of mass print "shift molecule to center of mass" pos = XYZ.atomlist2vector(atomlist) masses = AtomicData.atomlist2masses(atomlist) pos_com = MolCo.shift_to_com(pos, masses) atomlist = XYZ.vector2atomlist(pos_com, atomlist) # compute molecular orbitals with DFTB tddftb = LR_TDDFTB(atomlist) tddftb.setGeometry(atomlist, charge=0) options={"nstates": 1} try: tddftb.getEnergies(**options) except DFTB.Solver.ExcitedStatesNotConverged: pass valorbs, radial_val = load_pseudo_atoms(atomlist) if dyson_file == None: print "tight-binding Kohn-Sham orbitals are taken as Dyson orbitals" H**O, LUMO = tddftb.dftb2.getFrontierOrbitals() bound_orbs = tddftb.dftb2.getKSCoefficients() if selected_orbitals == None: # all orbitals selected_orbitals = range(0,bound_orbs.shape[1]) else: selected_orbitals = eval(selected_orbitals, {}, {"H**O": H**O+1, "LUMO": LUMO+1}) print "Indeces of selected orbitals (counting from 1): %s" % selected_orbitals orbital_names = ["orb_%s" % o for o in selected_orbitals] selected_orbitals = np.array(selected_orbitals, dtype=int)-1 # counting from 0 dyson_orbs = bound_orbs[:,selected_orbitals] ionization_energies = -tddftb.dftb2.getKSEnergies()[selected_orbitals] else: print "coeffients for Dyson orbitals are read from '%s'" % dyson_file orbital_names, ionization_energies, dyson_orbs = load_dyson_orbitals(dyson_file) ionization_energies = np.array(ionization_energies) / AtomicData.hartree_to_eV print "" print "*******************************************" print "* PHOTOELECTRON ANGULAR DISTRIBUTIONS *" print "*******************************************" print "" # determine the radius of the sphere where the angular distribution is calculated. It should be # much larger than the extent of the molecule (xmin,xmax),(ymin,ymax),(zmin,zmax) = Cube.get_bbox(atomlist, dbuff=0.0) dx,dy,dz = xmax-xmin,ymax-ymin,zmax-zmin Rmax = max([dx,dy,dz]) + sphere_radius Npts = max(int(Rmax),1) * 50 print "Radius of sphere around molecule, Rmax = %s bohr" % Rmax print "Points on radial grid, Npts = %d" % Npts nr_dyson_orbs = len(orbital_names) # compute PADs for all selected orbitals for iorb in range(0, nr_dyson_orbs): print "computing photoangular distribution for orbital %s" % orbital_names[iorb] data_file = "betas_" + molecule_name + "_" + orbital_names[iorb] + ".dat" pad_data = [] print " SCAN" nskip = max(1, nskip) # save table fh = open(data_file, "w") print " Writing table with betas to %s" % data_file print>>fh, "# ionization from orbital %s IE = %6.3f eV" % (orbital_names[iorb], ionization_energies[iorb]*AtomicData.hartree_to_eV) print>>fh, "# inter_atomic: %s npts_euler: %s npts_theta: %s rmax: %s" % (inter_atomic, npts_euler, npts_theta, Rmax) print>>fh, "# PKE/eV sigma beta1 beta2 beta3 beta4" for i,E in enumerate(slako_tables_scattering.energies): if i % nskip != 0: continue print " PKE = %6.6f Hartree (%4.4f eV)" % (E, E*AtomicData.hartree_to_eV) k = np.sqrt(2*E) wavelength = 2.0 * np.pi/k bs_free = AtomicScatteringBasisSet(atomlist, E, rmin=0.0, rmax=Rmax+2*wavelength, Npts=Npts) SKT_bf, SKT_ff = load_slako_scattering(atomlist, E) Dipole = ScatteringDipoleMatrix(atomlist, valorbs, SKT_bf, inter_atomic=inter_atomic).real orientation_averaging = PAD.OrientationAveraging_small_memory(Dipole, bs_free, Rmax, E, npts_euler=npts_euler, npts_theta=npts_theta) pad,betasE = orientation_averaging.averaged_pad(dyson_orbs[:,iorb]) pad_data.append( [E*AtomicData.hartree_to_eV] + list(betasE) ) # save PAD for this energy print>>fh, "%10.6f %10.6e %+10.6e %+10.6f %+10.6e %+10.6e" % tuple(pad_data[-1]) fh.flush() fh.close()
def test_hmi_continuum(): """ check that the continuum wavefunction of H2+ really are solutions of Schroedinger's equation, i.e. have (H-E)\phi = 0 everywhere starting from an LCAO guess for the continuum orbital, try to find the exact solution by adding orbital corrections iteratively """ # First we compute the exact wavefunction of the hydrogen molecular ion. from DFTB.Scattering import HMI # The bond length and charges cannot be changed, since the # separation constants were calculated only for the H2+ ion at R=2! R = 2.0 Za = 1.0 Zb = 1.0 # energy of continuum orbital E = 0.5 ## sigma (m=0) orbital m = 0 n = 0 trig = 'cos' # separation constant Lsep = HMI.SeparationConstants(R, Za, Zb) Lsep.load_separation_constants() Lfunc = Lsep.L_interpolated(m, n) c2 = 0.5 * E * R**2 mL, nL, L = Lfunc(c2) parity = (-1)**(mL + nL) phi_exact = HMI.create_wavefunction(mL, L, R * (Za + Zb), 0.0, R, c2, parity, trig) # Old implementation of H2+ wavefunctions, the wavefunction # looks indistinguishable from the exact wavefunction, but the # non-zero residue shows that is contains large errors. from DFTB.Scattering.hydrogen_molecular_ion import DimerWavefunctions wfn = DimerWavefunctions(R, Za, Zb, plot=False) delta, (Rfunc, Sfunc, Pfunc), wavefunction_exact = wfn.getContinuumOrbital( m, n, trig, E) def phi_exact_DW(x, y, z): return wavefunction_exact((x, y, z), None) # Set resolution of multicenter grid settings.radial_grid_factor = 10 settings.lebedev_order = 41 # Next we compute the wavefunction using the basis set free method atomlist = [(int(Za), (0.0, 0.0, -R / 2.0)), (int(Zb), (0.0, 0.0, +R / 2.0))] # no other electrons, only nuclear potential def potential(x, y, z): nuc = 0.0 * x for Zi, posi in atomlist: ri = np.sqrt((x - posi[0])**2 + (y - posi[1])**2 + (z - posi[2])**2) nuc += -Zi / ri return nuc # electron-electron interaction def potential_ee(x, y, z): return 0.0 * x # Set resolution of multicenter grid settings.radial_grid_factor = 10 settings.lebedev_order = 41 # residual of exact wavefunction (should be zero) residual_exact = residual_func(atomlist, phi_exact, potential, E) residual_ee_exact = residual_ee_func(atomlist, phi_exact, potential_ee, E) residual_exact_DW = residual_func(atomlist, phi_exact_DW, potential, E) # Laplacian laplacian_exact = laplacian_func(atomlist, phi_exact) import matplotlib.pyplot as plt plt.clf() r = np.linspace(-15.0, 15.0, 5000) x = 0 * r y = 0 * r z = r # plot exact wavefunction plt.plot(r, phi_exact(x, y, z), label="$\phi$ exact") # phi_exact_xyz = phi_exact(x, y, z) phi_exact_DW_xyz = phi_exact_DW(x, y, z) scale = phi_exact_xyz.max() / phi_exact_DW_xyz.max() plt.plot(r, scale * phi_exact_DW_xyz, label="$\phi$ exact (DimerWavefunction)") # and residual plt.plot(r, residual_exact(x, y, z), label=r"$(H-E)\phi$ (exact, old)") plt.plot(r, residual_exact_DW(x, y, z), ls="-.", label=r"$(H-E)\phi$ (exact, DimerWavefunction, old)") plt.plot(r, residual_ee_exact(x, y, z), ls="--", label=r"$(H-E)\phi$ (exact, new)") # kinetic energy plt.plot(r, -0.5 * laplacian_exact(x, y, z), ls="--", label=r"$-\frac{1}{2}\nabla^2 \phi$") # potential energy plt.plot(r, (potential(x, y, z) - E) * phi_exact(x, y, z), ls="--", label=r"$(V-E)\phi$") ## The initial guess for the \sigma continuum orbital ## is a regular Coulomb function centered on the midpoint ## between the two protons. #phi0 = regular_coulomb_func(E, +2, 0, 0, 0.0, center=(0.0, 0.0, 0.0)) """ ## The initial guess for the \sigma continuum orbital is ## the sum of two hydrogen s continuum orbitals bs = AtomicScatteringBasisSet(atomlist, E, lmax=0) phi0 = add_two_functions(atomlist, bs.bfs[0], bs.bfs[1], 1.0/np.sqrt(2.0), 1.0/np.sqrt(2.0)) """ """ ## start with exact wavefunction phi0 = phi_exact """ # The initial guess for the \sigma continuum orbital is # a hydrogen continuum orbital in the center bs = AtomicScatteringBasisSet([(1, (0.0, 0.0, 0.0))], E, lmax=0) phi0 = bs.bfs[0] plt.plot(r, phi0(x, y, z), ls="-.", label="guess $\phi_0$") plt.legend() plt.show() #phi = improve_continuum_orbital(atomlist, phi0, potential_ee, E, thresh=1.0e-6) phi = relax_continuum_orbital(atomlist, phi0, potential_ee, E, thresh=1.0e-6) import matplotlib.pyplot as plt plt.clf() r = np.linspace(-15.0, 15.0, 5000) x = 0 * r y = 0 * r z = r phi_exact_xyz = phi_exact(x, y, z) phi_xyz = phi(x, y, z) # scale numerical phi such that the maxima agree scale = phi_exact_xyz.max() / phi_xyz.max() phi_xyz *= scale print("scaling factor s = %s" % scale) plt.plot(r, phi_exact_xyz, label="exact") plt.plot(r, phi_xyz, label="numerical") plt.legend() plt.show()
def test_h2_continuum_orbital(): """ The sigma continuum orbital is approximated as a linear combination of two s continuum orbitals and is then improved iteratively """ # bond length in bohr dist = 2.0 # positions of protons posH1 = (0.0, 0.0, -dist / 2.0) posH2 = (0.0, 0.0, +dist / 2.0) atomlist = [(1, posH1), (1, posH2)] # Set resolution of multicenter grid settings.radial_grid_factor = 20 settings.lebedev_order = 23 # energy of continuum orbital E = 1.0 # same functional as used in the calculation of pseudo orbitals xc = XCFunctionals.libXCFunctional(Parameters.pseudo_orbital_x, Parameters.pseudo_orbital_c) dft = BasissetFreeDFT(atomlist, xc) print("initial orbital guess from DFTB calculation") orbitals = dft.getOrbitalGuess() norb = len(orbitals) # all orbitals are doubly occupied nelec = 2 * norb bound_orbitals = dft.getOrbitalGuess() # electron density (closed shell) rho = density_func(bound_orbitals) # effective Kohn-Sham potential veff = effective_potential_func(atomlist, rho, xc, nelec=nelec, nuclear=True) # effective Kohn-Sham potential without nuclear attraction # (only electron-electron interaction) veff_ee = effective_potential_func(atomlist, rho, xc, nelec=nelec, nuclear=False) ps = AtomicPotentialSet(atomlist) lmax = 0 bs = AtomicScatteringBasisSet(atomlist, E, lmax=lmax) #test_AO_basis(atomlist, bs, ps, E) R = residual2_matrix(atomlist, veff, ps, bs) S = continuum_overlap(bs.bfs, E) print("continuum overlap") print(S) print("residual^2 matrix") print(R) eigvals, eigvecs = sla.eigh(R, S) print(eigvals) print("eigenvector belonging to lowest eigenvalue") print(eigvecs[:, 0]) # LCAO continuum orbitals continuum_orbitals = orbital_transformation(atomlist, bs.bfs, eigvecs) # improve continuum orbital by adding a correction term # # phi = phi0 + dphi # # The orbital correction dphi is the solution of the inhomogeneous # Schroedinger equation # # (H-E)dphi = -(H-E)phi0 # print("orbital correction...") phi0 = continuum_orbitals[0] phi = improve_continuum_orbital(atomlist, phi0, veff_ee, E)
Cube.function_to_cubefile( atomlist, wavefunction2, filename="/tmp/h2+_continuum_orbital_%d_%d_%s.cube" % (m, n, str(E).replace(".", "p")), dbuff=15.0, ppb=2.5) from DFTB.Scattering import PAD PAD.asymptotic_density(wavefunction2, 20.0, E) plt.show() # build LCAO of two atomic s-waves with PKE=5 eV from DFTB.Scattering.SlakoScattering import AtomicScatteringBasisSet bs = AtomicScatteringBasisSet(atomlist, E) lcao_continuum = np.array([ +1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, -1.0, 0.0, 0.0, 0.0, 0.0,
def variational_kohn(atomlist, energy, lmax=2, rmax=300.0, npts_r=60, radial_grid_factor=3, lebedev_order=23): """ find the scattering states for a molecular Hamiltonian H = T + sum_k Vk where the molecular potential consists of a superposition of atomic effective potentials Vk Asymptotically the scattering orbital can be classified by its angular momentum l,m. The scattering orbital is assumed to have the form: |Psi_{l,m}> = cos(delta_{l,m}) |S_{l,m}> + sin(delta_{l,m}) |C_{l,m}> + sum_i C_{l,m,i} |Ai,li,mi> where |S_{l,m}> and |C_{l,m}> are sine- and cosine-like free solutions of (T-E)|Psi> = 0 with angular quantum numbers l,m and |Ai,li,mi> are atomic continuum orbitals, which are solutions of (T+Vi-E)|Psi> = 0. delta_{l,m} is the phase shift The coefficients C and the phase shift delta for each partial wave are determined from a variational principle. Parameters ---------- atomlist : list of tuples (Z,[x,y,z]) with atom number and nuclear position energy : energy E=1/2 k^2 (in Hartree) Optional -------- lmax : highest angular momentum of atomic continuum orbitals 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] radial_grid_factor: factor by which the number of radial grid points is increased for integration on the interval [0,+inf] lebedev_order : order of Lebedev grid for angular integral Returns ------- continuum_orbitals : list of continuum orbitals |Psi_{l,m}>, which are instances of LinearCombinationWavefunction phase_shifts : list of phase shift delta_{l,m} lms : list of tuples (l,m) with the asymptotic angular momentum quantum numbers for each continuum orbital """ # molecular potential V = sum_k Vk pot0 = AtomicPotentialSuperposition(atomlist, confined=False) # basis of atomic scattering states # The atomic basis functions |i> = |A,l,m> belonging to atom A # are solution of # (T + V_A - E)|A,l,m> = 0 # basis = AtomicScatteringBasisSet(atomlist, energy, lmax=lmax) # sets of linearly independent solutions for each partial wave l,m # (T - E)|S_{l,m}> = 0 # (T - E)|C_{l,m}> = 0 basis0 = AsymptoticSolutionsBasisSet(atomlist, energy, lmax=lmax) # The basis functions can thus be grouped into three categories # * atomic scattering functions centered on each atom : |A,l,m> # for each atom A # and for l=0,...,lmax m=-lmax,...,0,...,lmax # * sine-like spherical waves at the center of mass : |S_{l,m}> # for l=0,...,lmax m=-lmax,...,0,...,lmax # * cosine-like spherical waves at the center of mass : |C_{l,m}> # for l=0,...,lmax m=-lmax,...,0,...,lmax # count basis functions # * on atomic centers nb = len(basis.bfs) # * sine-like basis functions ns = len(basis0.bfsS) # * cosine-like basis functions nc = len(basis0.bfsC) # Now we need to compute the matrix elements of the operator L=H-E # between all basis functions # |A,l,m> |S_{l,m}> |C_{l,m}> bfs = basis.bfs + basis0.bfsS + basis0.bfsC L = scattering_integrals(atomlist, bfs, pot0, radial_grid_factor=radial_grid_factor, lebedev_order=lebedev_order) # extract blocks from L-matrix, which is not Hermitian # # Lbb Lbs Lbc # L = Lsb Lss Lsc # Lcb Lcs Lcc # # Lbb[i,j] = <Ai,li,mi|(H-E)|Aj,lj,mj> Lbb = L[:nb, :nb] # Lbs[i,j] = <Ai,li,mi|(H-E)|S_{lj,mj}> Lbs = L[:nb, nb:nb + ns] # Lbc[i,j] = <Ai,li,mi|(H-E)|C_{lj,mj}> Lbc = L[:nb, nb + ns:nb + ns + nc] # Lss[i,j] = <S_{li,mi}|(H-E)|S_{lj,mj}> Lss = L[nb:nb + ns, nb:nb + ns] # Lsb[i,j] = <S_{li,mi}|(H-E)|Aj,lj,mj> Lsb = L[nb:nb + ns, :nb] # Lsc[i,j] = <S_{li,mi}|(H-E)|C_{lj,mj}> Lsc = L[nb:nb + ns, nb + ns:nb + ns + nc] # Lcs[i,j] = <C_{li,mi}|(H-E)|S_{lj,mj}> Lcs = L[nb + ns:nb + ns + nc, nb:nb + ns] # Lcb[i,j] = <C_{li,mi}|(H-E)|Aj,lj,mj> Lcb = L[nb + ns:nb + ns + nc, :nb] # Lcc[i,j] = <C_{li,mi}|(H-E)|C_{lj,mj}> Lcc = L[nb + ns:nb + ns + nc, nb + ns:nb + ns + nc] # Determine coefficients cS and cC in the expansion # |uS_{l,m}> = |S_{l,m}> + sum_k cS_{l,m;Ak,lk,mk} |Ak,lk,mk> # |uC_{l,m}> = |C_{l,m}> + sum_k cC_{l,m;Ak,lk,mk} |Ak,lk,mk> # from the conditions # (1) <Ai,li,mi|(H-E)|uS_{l,m}> = 0 # (2) <Ai,li,mi|(H-E)|uC_{l,m}> = 0 cS = -la.solve(Lbb, Lbs) # solution of eqn. (1) Lbs + Lbb.cS = 0 cC = -la.solve(Lbb, Lbc) # solution of eqn. (2) Lbc + Lbb.cC = 0 # The partial waves are assumed to have the form # |Psi_{l,m}> = |uS_{l,m}> + t_{l,m} |uC_{l,m}> # The tangent of the phase shift t_{l,m} = tan(delta_{l,m}) is # determined from the conditions # (3) <uS_{l,m}|(H-E)|Psi_{l,m}> = 0 # (4) <uC_{l,m}|(H-E)|Psi_{l,m}> = 0 # This leads to the following 2 x 2 system of equations # (Mss Msc) (1) # ( ) * ( ) = 0 # (Mcs Mcc) (t) # Both equations cannot be satisfied at the same time, so we # have to choose one of them # # coefficients are real, so we do not really need complex conjugation cSt = cS.conjugate().transpose() cCt = cC.conjugate().transpose() Mss = Lss + np.dot(cSt, Lbs) + np.dot(Lsb, cS) + np.dot( cSt, np.dot(Lbb, cS)) Msc = Lsc + np.dot(cSt, Lbc) + np.dot(Lsb, cC) + np.dot( cSt, np.dot(Lbb, cC)) Mcs = Lcs + np.dot(cCt, Lbs) + np.dot(Lcb, cS) + np.dot( cCt, np.dot(Lbb, cS)) Mcc = Lcc + np.dot(cCt, Lbc) + np.dot(Lcb, cC) + np.dot( cCt, np.dot(Lbb, cC)) continuum_orbitals = [] phase_shifts = [] lms = [] for i in range(0, ns): l, m = basis0.bfsS[i].l, basis0.bfsS[i].m # solve Mss + t Msc = 0 t = -Mss[i, i] / Msc[i, i] # the other possibility would be to solve Mcs + t Mcc = 0 # t = -Mcs[i,i]/Mcc[i,i] # compute phase shift delta = np.arctan(t) # Because a global phase is irrelevant, the phase shift is only # determined module pi. sin(pi+delta) = -sin(delta) while delta < 0.0: delta += np.pi # and normalization factor # sin(kr+delta) = sin(kr) cos(delta) + sin(delta) cos(kr) # = cos(delta) [ sin(kr) + tan(delta) cos(kr)] # normalized # |Psi_{l,m}> = cos(delta) |uS_{l,m}> + sin(delta) |uC_{l,m}> # = cos(delta) |S_{l,m}> + sin(delta) |C_{l,m}> # + sum_k [ cos(delta) cS_{l,m;Ak,lk,mk} + sin(delta) cC_{l,m;Ak,lk,mk} ] |Ak,lk,mk> # basis functions and coefficients bfs = [basis0.bfsS[i], basis0.bfsC[i]] + basis.bfs coeffs = np.zeros(2 + nb) coeffs[0] = np.cos(delta) coeffs[1] = np.sin(delta) coeffs[2:] = np.cos(delta) * cS[:, i] + np.sin(delta) * cC[:, i] # create continuum wavefunction psi = LinearCombinationWavefunction(bfs, coeffs) continuum_orbitals.append(psi) phase_shifts.append(delta) lms.append((l, m)) # All continuum orbitals are rescaled / normalized, such that asymptotically # their radial parts tend to 1/r sin(k*r + delta) scaling_factors = continuum_normalization(continuum_orbitals, energy, rmax=rmax, npts_r=npts_r, lebedev_order=lebedev_order) for iw, wfn in enumerate(continuum_orbitals): # scale coefficients of linear combination of basis functions by sqrt(Is/I) wfn.coeffs *= scaling_factors[iw] return continuum_orbitals, phase_shifts, lms
def test_lcao_continuum(): import matplotlib.pyplot as plt # bond length in bohr dist = 2.0 # positions of protons posH1 = (0.0, 0.0, -dist / 2.0) posH2 = (0.0, 0.0, +dist / 2.0) atomlist = [(1, posH1), (1, posH2)] # Set resolution of multicenter grid settings.radial_grid_factor = 20 settings.lebedev_order = 23 # energy of continuum orbital E = 1.0 # same functional as used in the calculation of pseudo orbitals xc = XCFunctionals.libXCFunctional(Parameters.pseudo_orbital_x, Parameters.pseudo_orbital_c) dft = BasissetFreeDFT(atomlist, xc) print("initial orbital guess from DFTB calculation") orbitals = dft.getOrbitalGuess() norb = len(orbitals) # all orbitals are doubly occupied nelec = 2 * norb bound_orbitals = dft.getOrbitalGuess() # effective potential rho = density_func(bound_orbitals) veff = effective_potential_func(atomlist, rho, xc, nelec=nelec) ps = AtomicPotentialSet(atomlist) r = np.linspace(-15.0, 15.0, 10000) x = 0.0 * r y = 0.0 * r z = r for lmax in [0, 1, 2, 3]: bs = AtomicScatteringBasisSet(atomlist, E, lmax=lmax) #test_AO_basis(atomlist, bs, ps, E) R = residual2_matrix(atomlist, veff, ps, bs) S = continuum_overlap(bs.bfs, E) print("continuum overlap") print(S) print("residual^2 matrix") print(R) eigvals, eigvecs = sla.eigh(R, S) print(eigvals) print("eigenvector belonging to lowest eigenvalue") print(eigvecs[:, 0]) # LCAO continuum orbitals continuum_orbitals = orbital_transformation(atomlist, bs.bfs, eigvecs) # improve continuum orbital by adding a correction term # # phi = phi0 + dphi # # The orbital correction dphi is the solution of the inhomogeneous # Schroedinger equation # # (H-E)dphi = -(H-E)phi0 # print("orbital correction...") phi0 = continuum_orbitals[0] phi = improve_continuum_orbital(atomlist, phi0, veff, E) exit(-1) residual_0 = residual_func(atomlist, phi0, veff, E) def source(x, y, z): return -residual_0(x, y, z) delta_phi = inhomogeneous_schroedinger(atomlist, veff, source, E) residual_d = residual_func(atomlist, delta_phi, veff, E) a, b = variational_mixture_continuum(atomlist, phi0, delta_phi, veff, E) phi = add_two_functions(atomlist, phi0, delta_phi, a, b) residual = residual_func(atomlist, phi, veff, E) plt.plot(r, 1.0 / np.sqrt(2.0) * bs.bfs[0](x, y, z), label=r"AO") plt.plot(r, phi0(x, y, z), label=r"$\phi_0$") plt.plot(r, delta_phi(x, y, z), label=r"$\Delta \phi$") plt.plot(r, phi(x, y, z), label=r"$\phi_0 + \Delta \phi$") plt.legend() plt.show() """ dphi = delta_phi(x,y,z) imin = np.argmin(abs(r-1.0)) dphi[abs(r) < 1.0] = dphi[imin] - (dphi[abs(r) < 1.0] - dphi[imin]) plt.plot(r, dphi, label=r"$\Delta \phi$") """ plt.plot(r, residual_0(x, y, z), label=r"$(H-E) \phi_0$") plt.plot(r, residual_d(x, y, z), label=r"$(H-E)\Delta \phi$") plt.plot(r, residual(x, y, z), label=r"$(H-E)(a \phi_0 + b \Delta \phi)$") plt.plot(r, a * residual_0(x, y, z) + b * residual_d(x, y, z), ls="-.", label=r"$(H-E)(a \phi_0 + b \Delta \phi)$ (separate)") plt.legend() plt.show() averaged_angular_distribution(atomlist, bound_orbitals, continuum_orbitals, E) # save continuum MOs to cubefiles for i, phi in enumerate(continuum_orbitals): def func(grid, dV): x, y, z = grid return phi(x, y, z) Cube.function_to_cubefile( atomlist, func, filename="/tmp/cmo_lmax_%2.2d_orb%4.4d.cube" % (lmax, i), ppb=5.0) # for i, phi in enumerate(continuum_orbitals): residual = residual_func(atomlist, phi, veff, E) delta_e = energy_correction(atomlist, residual, phi, method="Becke") print(" orbital %d energy <%d|H-E|%d> = %e" % (i, i, i, delta_e)) l, = plt.plot(r, phi(x, y, z), label=r"$\phi_{%d}$ ($l_{max}$ = %d)" % (i, lmax)) plt.plot(r, residual(x, y, z), ls="-.", label=r"$(H-E)\phi_{%d}$" % i, color=l.get_color()) plt.legend() plt.show()
def test_hmi_lcao_continuum(): """ """ # First we compute the exact wavefunction of the hydrogen molecular ion. from DFTB.Scattering import HMI R = 2.0 Za = 1.0 Zb = 1.0 E = 0.5 ## sigma (m=0) orbital m = 0 n = 0 L = m + n trig = 'cos' # separation constant Lsep = HMI.SeparationConstants(R, Za, Zb) Lsep.load_separation_constants() Lfunc = Lsep.L_interpolated(m, n) c2 = 0.5 * E * R**2 mL, nL, L = Lfunc(c2) parity = (-1)**(mL + nL) phi_exact = HMI.create_wavefunction(mL, L, R * (Za + Zb), 0.0, R, c2, parity, trig) """ from DFTB.Scattering.hydrogen_molecular_ion import DimerWavefunctions wfn = DimerWavefunctions(R,Za,Zb, plot=False) delta, (Rfunc,Sfunc,Pfunc),wavefunction_exact = wfn.getContinuumOrbital(m,n,trig,E) def phi_exact(x,y,z): return wavefunction_exact((x,y,z), None) """ # Next we compute the wavefunction using the basis set free method atomlist = [(int(Za), (0.0, 0.0, -R / 2.0)), (int(Zb), (0.0, 0.0, +R / 2.0))] # no other electrons, only nuclear potential def potential(x, y, z): nuc = 0.0 * x for Zi, posi in atomlist: ri = np.sqrt((x - posi[0])**2 + (y - posi[1])**2 + (z - posi[2])**2) nuc += -Zi / ri return nuc # Set resolution of multicenter grid settings.radial_grid_factor = 40 settings.lebedev_order = 21 # residual of exact wavefunction (should be zero) residual_exact = residual_func(atomlist, phi_exact, potential, E) # Laplacian laplacian_exact = laplacian_func(atomlist, phi_exact) import matplotlib.pyplot as plt plt.clf() r = np.linspace(-15.0, 15.0, 5000) x = 0 * r y = 0 * r z = r # plot exact wavefunction plt.plot(r, phi_exact(x, y, z), label="$\phi$ exact") # and residual plt.plot(r, residual_exact(x, y, z), label=r"$(H-E)\phi$ (exact)") # kinetic energy plt.plot(r, -0.5 * laplacian_exact(x, y, z), ls="--", label=r"$-\frac{1}{2}\nabla^2 \phi$") # potential energy plt.plot(r, (potential(x, y, z) - E) * phi_exact(x, y, z), ls="--", label=r"$(V-E)\phi$") # LCAO continuum orbitals ps = AtomicPotentialSet(atomlist) lmax = 4 bs = AtomicScatteringBasisSet(atomlist, E, lmax=lmax) R = residual2_matrix(atomlist, potential, ps, bs) S = continuum_overlap(bs.bfs, E) print("continuum overlap") print(S) print("residual^2 matrix") print(R) eigvals, eigvecs = sla.eigh(R, S) print(eigvals) print("eigenvector belonging to lowest eigenvalue") print(eigvecs[:, 0]) continuum_orbitals = orbital_transformation(atomlist, bs.bfs, eigvecs) plt.cla() plt.plot(r, phi_exact(x, y, z), ls="--", label="$\phi$ (exact)") for i in range(0, len(continuum_orbitals)): plt.plot(r, continuum_orbitals[i](x, y, z)) plt.legend() plt.show()