Differences
This shows you the differences between two versions of the page.
Next revision | Previous revision | ||
documentation:tutorials:nio_crystal_field:nixs_l23 [2016/10/08 21:25] – created Maurits W. Haverkort | documentation:tutorials:nio_crystal_field:nixs_l23 [2019/02/21 08:25] (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 :-). | ||
+ | ### | ||
+ | |||
+ | ### | ||
+ | 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 set the output of the program to a minimum | ||
+ | Verbosity(0) | ||
+ | |||
+ | -- we need a 2p and 3d shell | ||
+ | NF=16 | ||
+ | NB=0 | ||
+ | IndexDn_2p={0, | ||
+ | IndexUp_2p={1, | ||
+ | IndexDn_3d={6, | ||
+ | IndexUp_3d={7, | ||
+ | |||
+ | OppSx | ||
+ | OppSy | ||
+ | OppSz | ||
+ | OppSsqr =NewOperator(" | ||
+ | OppSplus=NewOperator(" | ||
+ | OppSmin =NewOperator(" | ||
+ | |||
+ | OppLx | ||
+ | OppLy | ||
+ | OppLz | ||
+ | OppLsqr =NewOperator(" | ||
+ | OppLplus=NewOperator(" | ||
+ | OppLmin =NewOperator(" | ||
+ | |||
+ | OppJx | ||
+ | OppJy | ||
+ | OppJz | ||
+ | OppJsqr =NewOperator(" | ||
+ | OppJplus=NewOperator(" | ||
+ | OppJmin =NewOperator(" | ||
+ | |||
+ | Oppldots=NewOperator(" | ||
+ | |||
+ | -- define the coulomb operator | ||
+ | -- we here define the part depending on F0 seperately from the part depending on F2 | ||
+ | -- when summing we can put in the numerical values of the slater integrals | ||
+ | OppF0 =NewOperator(" | ||
+ | OppF2 =NewOperator(" | ||
+ | OppF4 =NewOperator(" | ||
+ | |||
+ | Akm = PotentialExpandedOnClm(" | ||
+ | OpptenDq = NewOperator(" | ||
+ | |||
+ | Akm = PotentialExpandedOnClm(" | ||
+ | OppNeg = NewOperator(" | ||
+ | Akm = PotentialExpandedOnClm(" | ||
+ | OppNt2g = NewOperator(" | ||
+ | |||
+ | Oppcldots= NewOperator(" | ||
+ | OppUpdF0 = NewOperator(" | ||
+ | OppUpdF2 = NewOperator(" | ||
+ | OppUpdG1 = NewOperator(" | ||
+ | OppUpdG3 = NewOperator(" | ||
+ | |||
+ | -- in crystal field theory U drops out of the equation | ||
+ | U | ||
+ | F2dd = 11.142 | ||
+ | F4dd = 6.874 | ||
+ | F0dd = U+(F2dd+F4dd)*2/ | ||
+ | -- in crystal field theory U drops out of the equation | ||
+ | Upd | ||
+ | F2pd = 6.667 | ||
+ | G1pd = 4.922 | ||
+ | G3pd = 2.796 | ||
+ | F0pd = Upd + G1pd*1/15 + G3pd*3/70 | ||
+ | tenDq | ||
+ | zeta_3d = 0.081 | ||
+ | zeta_2p = 11.498 | ||
+ | Bz = 0.000001 | ||
+ | |||
+ | Hamiltonian = F0dd*OppF0 + F2dd*OppF2 + F4dd*OppF4 + tenDq*OpptenDq + zeta_3d*Oppldots + Bz*(2*OppSz + OppLz) | ||
+ | |||
+ | XASHamiltonian = Hamiltonian + zeta_2p * Oppcldots + F2pd * OppUpdF2 + G1pd * OppUpdG1 + G3pd * OppUpdG3 | ||
+ | |||
+ | -- we now can create the lowest Npsi eigenstates: | ||
+ | Npsi=3 | ||
+ | -- in order to make sure we have a filling of 2 electrons we need to define some restrictions | ||
+ | StartRestrictions = {NF, NB, {" | ||
+ | |||
+ | psiList = Eigensystem(Hamiltonian, | ||
+ | |||
+ | oppList={Hamiltonian, | ||
+ | |||
+ | print(" | ||
+ | for key,psi in pairs(psiList) do | ||
+ | expvalue = psi * oppList * psi | ||
+ | for k,v in pairs(expvalue) do | ||
+ | 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: | ||
+ | ### | ||
+ | |{{: | ||
+ | ^ $2p$ to $3d$ excitations as one would measure using non-resonant inelastic x-ray scattering. | ||
+ | |||
+ | |||
+ | ### | ||
+ | The output to standard out is: | ||
+ | <file Quanty_Output nIXS_L23.out> | ||
+ | -2.444 | ||
+ | -2.444 | ||
+ | -2.444 | ||
+ | for q= 9 per a0 ( 17.007535121086 per A) The ratio of k=1 and k=3 transition strength is: | ||
+ | </ | ||
+ | ### | ||
+ | |||
+ | |||
+ | ===== Table of contents ===== | ||
+ | {{indexmenu> |