Differences
This shows you the differences between two versions of the page.
Next revision | Previous revision | ||
documentation:tutorials:nio_ligand_field:nixs_l23 [2016/10/09 15:51] – created Maurits W. Haverkort | documentation:tutorials:nio_ligand_field:nixs_l23 [2018/03/20 11:08] (current) – Maurits W. Haverkort | ||
---|---|---|---|
Line 1: | Line 1: | ||
+ | {{indexmenu_n> | ||
+ | ====== nIXS $L_{2,3}$ ====== | ||
+ | ### | ||
+ | Besides low energy transitions nIXS can be used as a core level spectroscopy technique. One then measures resonances with non-resonant inelastic x-ray scattering :-). | ||
+ | ### | ||
+ | |||
+ | ### | ||
+ | This tutorial uses Radial wave functions in order to calculate the length of q dependence. You can download them in a zip file here {{ : | ||
+ | ### | ||
+ | |||
+ | ### | ||
+ | The input script: | ||
+ | <code Quanty NIXS_L23.Quanty> | ||
+ | -- using inelastic x-ray scattering one can not only measure low energy excitations, | ||
+ | -- but equally well core to core transitions. This allows one to probe for example | ||
+ | -- 3p to 3d transitions using octupole operators. | ||
+ | |||
+ | -- We use the definitions of all operators and basis orbitals as defined in the file | ||
+ | -- include and can afterwards directly continue by creating the Hamiltonian | ||
+ | -- and calculating the spectra | ||
+ | |||
+ | dofile(" | ||
+ | |||
+ | -- The parameters and scheme needed is the same as the one used for XAS | ||
+ | |||
+ | -- We follow the energy definitions as introduced in the group of G.A. Sawatzky (Groningen) | ||
+ | -- J. Zaanen, G.A. Sawatzky, and J.W. Allen PRL 55, 418 (1985) | ||
+ | -- for parameters of specific materials see | ||
+ | -- A.E. Bockquet et al. PRB 55, 1161 (1996) | ||
+ | -- After some initial discussion the energies U and Delta refer to the center of a configuration | ||
+ | -- The L^10 d^n | ||
+ | -- The L^9 d^n+1 configuration has an energy Delta | ||
+ | -- The L^8 d^n+2 configuration has an energy 2*Delta+Udd | ||
+ | -- | ||
+ | -- If we relate this to the onsite energy of the L and d orbitals we find | ||
+ | -- 10 eL + n ed + n(n-1) | ||
+ | -- 9 eL + (n+1) ed + (n+1)n | ||
+ | -- 8 eL + (n+2) ed + (n+1)(n+2) U/2 == 2*Delta+U | ||
+ | -- 3 equations with 2 unknowns, but with interdependence yield: | ||
+ | -- ed = (10*Delta-nd*(19+nd)*U/ | ||
+ | -- eL = nd*((1+nd)*Udd/ | ||
+ | -- | ||
+ | -- For the final state we/they defined | ||
+ | -- The 2p^5 L^10 d^n+1 configuration has an energy 0 | ||
+ | -- The 2p^5 L^9 d^n+2 configuration has an energy Delta + Udd - Upd | ||
+ | -- The 2p^5 L^8 d^n+3 configuration has an energy 2*Delta + 3*Udd - 2*Upd | ||
+ | -- | ||
+ | -- If we relate this to the onsite energy of the p and d orbitals we find | ||
+ | -- 6 ep + 10 eL + n ed + n(n-1) | ||
+ | -- 6 ep + 9 eL + (n+1) ed + (n+1)n | ||
+ | -- 6 ep + 8 eL + (n+2) ed + (n+1)(n+2) Udd/2 + 6 (n+2) Upd == 2*Delta+Udd | ||
+ | -- 5 ep + 10 eL + (n+1) ed + (n+1)(n) | ||
+ | -- 5 ep + 9 eL + (n+2) ed + (n+2)(n+1) Udd/2 + 5 (n+2) Upd == Delta+Udd-Upd | ||
+ | -- 5 ep + 8 eL + (n+3) ed + (n+3)(n+2) Udd/2 + 5 (n+3) Upd == 2*Delta+3*Udd-2*Upd | ||
+ | -- 6 equations with 3 unknowns, but with interdependence yield: | ||
+ | -- epfinal = (10*Delta + (1+nd)*(nd*Udd/ | ||
+ | -- edfinal = (10*Delta - nd*(31+nd)*Udd/ | ||
+ | -- eLfinal = ((1+nd)*(nd*Udd/ | ||
+ | -- | ||
+ | -- | ||
+ | -- | ||
+ | -- note that ed-ep = Delta - nd * U and not Delta | ||
+ | -- note furthermore that ep and ed here are defined for the onsite energy if the system had | ||
+ | -- locally nd electrons in the d-shell. In DFT or Hartree Fock the d occupation is in the end not | ||
+ | -- nd and thus the onsite energy of the Kohn-Sham orbitals is not equal to ep and ed in model | ||
+ | -- calculations. | ||
+ | -- | ||
+ | -- note furthermore that ep and eL actually should be different for most systems. We happily ignore this fact | ||
+ | -- | ||
+ | -- We normally take U and Delta as experimentally determined parameters | ||
+ | |||
+ | -- number of electrons (formal valence) | ||
+ | nd = 8 | ||
+ | -- parameters from experiment (core level PES) | ||
+ | Udd | ||
+ | Upd | ||
+ | Delta | ||
+ | -- parameters obtained from DFT (PRB 85, 165113 (2012)) | ||
+ | F2dd = 11.14 | ||
+ | F4dd = 6.87 | ||
+ | F2pd = 6.67 | ||
+ | G1pd = 4.92 | ||
+ | G3pd = 2.80 | ||
+ | tenDq | ||
+ | tenDqL | ||
+ | Veg | ||
+ | Vt2g = 1.21 | ||
+ | zeta_3d = 0.081 | ||
+ | zeta_2p = 11.51 | ||
+ | Bz = 0.000001 | ||
+ | H112 = 0 | ||
+ | |||
+ | ed = (10*Delta-nd*(19+nd)*Udd/ | ||
+ | eL = nd*((1+nd)*Udd/ | ||
+ | |||
+ | epfinal = (10*Delta + (1+nd)*(nd*Udd/ | ||
+ | edfinal = (10*Delta - nd*(31+nd)*Udd/ | ||
+ | eLfinal = ((1+nd)*(nd*Udd/ | ||
+ | |||
+ | F0dd = Udd + (F2dd+F4dd) * 2/63 | ||
+ | F0pd = Upd + (1/15)*G1pd + (3/70)*G3pd | ||
+ | |||
+ | Hamiltonian = F0dd*OppF0_3d + F2dd*OppF2_3d + F4dd*OppF4_3d + zeta_3d*Oppldots_3d + Bz*(2*OppSz_3d + OppLz_3d) + H112 * (OppSx_3d+OppSy_3d+2*OppSz_3d)/ | ||
+ | | ||
+ | XASHamiltonian = F0dd*OppF0_3d + F2dd*OppF2_3d + F4dd*OppF4_3d + zeta_3d*Oppldots_3d + Bz*(2*OppSz_3d + OppLz_3d)+ H112 * (OppSx_3d+OppSy_3d+2*OppSz_3d)/ | ||
+ | |||
+ | -- we now can create the lowest Npsi eigenstates: | ||
+ | Npsi=3 | ||
+ | -- in order to make sure we have a filling of 8 electrons we need to define some restrictions | ||
+ | StartRestrictions = {NF, NB, {" | ||
+ | |||
+ | psiList = Eigensystem(Hamiltonian, | ||
+ | oppList={Hamiltonian, | ||
+ | |||
+ | -- print of some expectation values | ||
+ | print(" | ||
+ | for i = 1,#psiList do | ||
+ | io.write(string.format(" | ||
+ | for j = 1,#oppList do | ||
+ | expectationvalue = Chop(psiList[i]*oppList[j]*psiList[i]) | ||
+ | io.write(string.format(" | ||
+ | end | ||
+ | io.write(" | ||
+ | end | ||
+ | |||
+ | |||
+ | -- in order to calculate nIXS we need to determine the intensity ratio for the different multipole intensities | ||
+ | -- ( see PRL 99, 257401 (2007) for the formalism ) | ||
+ | -- in short the A^2 interaction is expanded on spherical harmonics and Bessel functions | ||
+ | -- The 3d Wannier functions are expanded on spherical harmonics and a radial wave function | ||
+ | -- For the radial wave-function we calculate <R(r) | j_k(q r) | R(r)> | ||
+ | -- which defines the transition strength for the multipole of order k | ||
+ | |||
+ | -- The radial functions here are calculated for a Ni 2+ atom and stored in the folder NiO_Radial | ||
+ | -- more sophisticated methods can be used | ||
+ | |||
+ | -- read the radial wave functions | ||
+ | -- order of functions | ||
+ | -- r 1S 2S 2P 3S 3P 3D | ||
+ | file = io.open( " | ||
+ | Rnl = {} | ||
+ | for line in file: | ||
+ | RnlLine={} | ||
+ | for i in string.gmatch(line, | ||
+ | table.insert(RnlLine, | ||
+ | end | ||
+ | table.insert(Rnl, | ||
+ | end | ||
+ | |||
+ | -- some constants | ||
+ | a0 = 0.52917721092 | ||
+ | Rydberg = 13.60569253 | ||
+ | Hartree = 2*Rydberg | ||
+ | |||
+ | -- pd transitions from 2p (index 4 in Rnl) to 3d (index 7 in Rnl) | ||
+ | -- <R(r) | j_k(q r) | R(r)> | ||
+ | function RjRpd (q) | ||
+ | Rj1R = 0 | ||
+ | Rj3R = 0 | ||
+ | dr = Rnl[3][1]-Rnl[2][1] | ||
+ | r0 = Rnl[2][1]-2*dr | ||
+ | for ir = 2, #Rnl, 1 do | ||
+ | r = r0 + ir * dr | ||
+ | Rj1R = Rj1R + Rnl[ir][4] * math.SphericalBesselJ(1, | ||
+ | Rj3R = Rj3R + Rnl[ir][4] * math.SphericalBesselJ(3, | ||
+ | end | ||
+ | return Rj1R, Rj3R | ||
+ | end | ||
+ | |||
+ | -- the angular part is given as C(theta_q, phi_q)^* C(theta_r, phi_r) | ||
+ | -- which is a potential expanded on spherical harmonics | ||
+ | function ExpandOnClm(k, | ||
+ | ret={} | ||
+ | for m=-k, k, 1 do | ||
+ | table.insert(ret, | ||
+ | end | ||
+ | return ret | ||
+ | end | ||
+ | |||
+ | -- define nIXS transition operators | ||
+ | function TnIXS_pd(q, theta, phi) | ||
+ | Rj1R, Rj3R = RjRpd(q) | ||
+ | k=1 | ||
+ | A1 = ExpandOnClm(k, | ||
+ | T1 = NewOperator(" | ||
+ | k=3 | ||
+ | A3 = ExpandOnClm(k, | ||
+ | T3 = NewOperator(" | ||
+ | T = T1+T3 | ||
+ | T.Chop() | ||
+ | return T | ||
+ | end | ||
+ | |||
+ | -- q in units per a0 (if you want in units per A take 5*a0 to have a q of 5 per A) | ||
+ | q=9.0 | ||
+ | |||
+ | print(" | ||
+ | |||
+ | -- define some transition operators | ||
+ | qtheta=0 | ||
+ | qphi=0 | ||
+ | Tq001 = TnIXS_pd(q, | ||
+ | |||
+ | qtheta=Pi/2 | ||
+ | qphi=Pi/4 | ||
+ | Tq110 = TnIXS_pd(q, | ||
+ | |||
+ | qtheta=math.acos(math.sqrt(1/ | ||
+ | qphi=Pi/4 | ||
+ | Tq111 = TnIXS_pd(q, | ||
+ | |||
+ | qtheta=math.acos(math.sqrt(9/ | ||
+ | qphi=math.acos(math.sqrt(1/ | ||
+ | Tq123 = TnIXS_pd(q, | ||
+ | |||
+ | -- calculate the spectra | ||
+ | nIXSSpectra = CreateSpectra(XASHamiltonian, | ||
+ | |||
+ | -- print the spectra to a file | ||
+ | nIXSSpectra.Print({{" | ||
+ | |||
+ | -- a gnuplot script to make the plots | ||
+ | gnuplotInput = [[ | ||
+ | set autoscale | ||
+ | set xtic auto | ||
+ | set ytic auto | ||
+ | set style line 1 lt 1 lw 1 lc rgb "# | ||
+ | set style line 2 lt 1 lw 1 lc rgb "# | ||
+ | set style line 3 lt 1 lw 1 lc rgb "# | ||
+ | set style line 4 lt 1 lw 1 lc rgb "# | ||
+ | set style line 5 lt 1 lw 3 lc rgb "# | ||
+ | |||
+ | set xlabel "E (eV)" font " | ||
+ | set ylabel " | ||
+ | |||
+ | set out ' | ||
+ | set size 1.0, 0.3 | ||
+ | set terminal postscript portrait enhanced color " | ||
+ | |||
+ | energyshift=857.6 | ||
+ | |||
+ | plot " | ||
+ | " | ||
+ | " | ||
+ | " | ||
+ | |||
+ | ]] | ||
+ | |||
+ | -- write the gnuplot script to a file | ||
+ | file = io.open(" | ||
+ | file: | ||
+ | file: | ||
+ | |||
+ | -- call gnuplot to execute the script | ||
+ | os.execute(" | ||
+ | -- transform to pdf and eps | ||
+ | os.execute(" | ||
+ | </ | ||
+ | ### | ||
+ | |||
+ | The spectrum produced: | ||
+ | | {{: | ||
+ | ^ nIXS for NiO looking at a $2p$ to $3d$ excitation ^ | ||
+ | |||
+ | ### | ||
+ | The output is to standard out: | ||
+ | <file Quanty_Output> | ||
+ | # < | ||
+ | 1 | ||
+ | 2 | ||
+ | 3 | ||
+ | for q= 9 per a0 ( 17.007535121086 per A) The ratio of k=1 and k=3 transition strength is: | ||
+ | </ | ||
+ | \lstinputlisting[style=output]{../ | ||
+ | ### | ||
+ | |||
+ | |||
+ | ===== Table of contents ===== | ||
+ | {{indexmenu> |