Differences
This shows you the differences between two versions of the page.
| documentation:tutorials:nio_ligand_field:rixs_l23m45 [2016/10/09 15:39] – created Maurits W. Haverkort | documentation:tutorials:nio_ligand_field:rixs_l23m45 [2025/11/20 03:29] (current) – external edit 127.0.0.1 | ||
|---|---|---|---|
| Line 1: | Line 1: | ||
| + | {{indexmenu_n> | ||
| + | ====== RIXS $L_{2, | ||
| + | ### | ||
| + | Using the function ResonantSpectra we can calculate inelastic x-ray scattering. (or other second order processess) Here we show the example of $d$-$d$ excitations in NiO using RIXS. Due to an ever increasing resolution this method gains rapidly in popularity and impact. See for example the combined work of Ghiringhelli \textit{et al.} For NiO see \cite{Ghiringhelli: | ||
| + | ### | ||
| + | |||
| + | ### | ||
| + | A script to calculate these spectra is: | ||
| + | <code Quanty RIXS_L23M45.Quanty> | ||
| + | -- this example calculates the resonant inelastic x-ray scattering in NiO we look at the | ||
| + | -- Ni L23M45 edge, i.e. we make an excitation from 2p to 3d (L23) and decay from the | ||
| + | -- 3d shell back to the 2p shell (final " | ||
| + | -- measure d-d excitations and or magnons. | ||
| + | |||
| + | -- We use the definitions of all operators and basis orbitals as defined in the file | ||
| + | -- 44_50_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.120 | ||
| + | |||
| + | 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 | ||
| + | |||
| + | |||
| + | -- Not the main task of this example, but in order to know where the resonance sits it is | ||
| + | -- nice to calculate the x-ray absorption spectra | ||
| + | XASSpectra = CreateSpectra(XASHamiltonian, | ||
| + | XASSpectra.Print({{" | ||
| + | |||
| + | -- We also calculate the fluorescence yield spectra which are equal to the integral over the | ||
| + | -- RIXS spectra | ||
| + | FYSpectra = CreateFluorescenceYield(XASHamiltonian, | ||
| + | FYSpectra.Print({{" | ||
| + | |||
| + | -- and we calculate the RIXS spectra. Note that in order to calculate RIXS you need to | ||
| + | -- specify two Hamiltonians. | ||
| + | -- These calculations can be lengthy and one should think that for each incoming energy | ||
| + | -- (E1) one needs to do a new calculation. In this case there are thus 120 calculations | ||
| + | RIXSSpectra = CreateResonantSpectra(XASHamiltonian, | ||
| + | RIXSSpectra.Print({{" | ||
| + | |||
| + | print(" | ||
| + | |||
| + | -- Once finished we can make some nice plots. The spectra are saved to disk in ASCII format | ||
| + | -- but I like to add gnuplot scripts so you can look at pictures imeidately | ||
| + | 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 out ' | ||
| + | set terminal postscript portrait enhanced color " | ||
| + | |||
| + | unset colorbox | ||
| + | |||
| + | energyshift=857.6 | ||
| + | |||
| + | set multiplot | ||
| + | set size 0.5,0.55 | ||
| + | set origin 0,0 | ||
| + | |||
| + | set ylabel " | ||
| + | set xlabel " | ||
| + | |||
| + | set yrange [852:860] | ||
| + | set xrange [-0.5:7.5] | ||
| + | |||
| + | plot "< | ||
| + | |||
| + | set origin 0.5,0 | ||
| + | |||
| + | set yrange [869:877] | ||
| + | set xrange [-0.5:7.5] | ||
| + | |||
| + | plot "< | ||
| + | |||
| + | unset multiplot | ||
| + | |||
| + | set out ' | ||
| + | set terminal postscript portrait enhanced color " | ||
| + | |||
| + | set multiplot | ||
| + | set size 0.25,1.0 | ||
| + | set origin 0,0 | ||
| + | |||
| + | set ylabel "E (eV)" font " | ||
| + | set xlabel " | ||
| + | set yrange [energyshift-10: | ||
| + | set xrange [-0.3:0] | ||
| + | plot " | ||
| + | " | ||
| + | |||
| + | set size 0.8,1.0 | ||
| + | set origin 0.2,0.0 | ||
| + | |||
| + | set xlabel " | ||
| + | unset ylabel | ||
| + | unset ytics | ||
| + | set xrange [-0.5:7.5] | ||
| + | |||
| + | ofset = 0.25 | ||
| + | scale=3 | ||
| + | |||
| + | plot for [i=0:120] " | ||
| + | |||
| + | unset multiplot | ||
| + | ]] | ||
| + | |||
| + | -- 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(" | ||
| + | os.execute(" | ||
| + | </ | ||
| + | ### | ||
| + | |||
| + | It produces two plots. | ||
| + | |||
| + | | {{ : | ||
| + | ^ Resonant inelastic x-ray scattering spectra for different incoming and outgoing photon energies. ^ | ||
| + | |||
| + | The first shows on the left the XAS spectra in black and the integrated RIXS spectra (FY) in blue. The resonant enhanced $d$-$d$ excitations are shown in the right plot. | ||
| + | |||
| + | | {{ : | ||
| + | ^ Resonant inelastic x-ray scattering spectra for different incoming and outgoing photon energies. ^ | ||
| + | |||
| + | The second plot shows the same RIXS spectra, but now as an intensity map. | ||
| + | |||
| + | ### | ||
| + | The text output of the script is: | ||
| + | <file Quanty_Output RIXS_L23M45.out> | ||
| + | # < | ||
| + | 1 | ||
| + | 2 | ||
| + | 3 | ||
| + | Start of LanczosTriDiagonalizeKrylovMC | ||
| + | Start of LanczosTriDiagonalizeKrylovMC | ||
| + | Finished calculating the spectra now start plotting. | ||
| + | This might take more time than the calculation | ||
| + | </ | ||
| + | ### | ||
| + | |||
| + | ===== Table of contents ===== | ||
| + | {{indexmenu> | ||