Differences
This shows you the differences between two versions of the page.
| Next revision | Previous revision | ||
| documentation:tutorials:nio_crystal_field:rixs_l23m1 [2016/10/08 21:19] – created Maurits W. Haverkort | documentation:tutorials:nio_crystal_field:rixs_l23m1 [2025/11/20 03:29] (current) – external edit 127.0.0.1 | ||
|---|---|---|---|
| Line 1: | Line 1: | ||
| + | {{indexmenu_n> | ||
| + | ====== RIXS $L_{2, | ||
| + | ### | ||
| + | Instead of looking at low energy excitations one can use RIXS to look at core-core excitations. A powerful technique comparable to other core level spectroscopies like x-ray absorption and core-level photoemission at once. | ||
| + | ### | ||
| + | |||
| + | ### | ||
| + | Here a script to calculate the Ni $2p$ to $3d$ excitation ($L_{2,3}$) and Ni $3s$ to $2p$ decay ($M_{1}$). | ||
| + | <code Quanty RIXS_L23M1.Quanty> | ||
| + | -- This example calculates core - core RIXS spectra. These type of spectra are | ||
| + | -- highly informative and contain similar information as core level absorption and | ||
| + | -- core level photoemission combined. With generally enhances sensitivity. | ||
| + | |||
| + | -- we focus here on 2p to 3d excitations (L23) and 3s to 2p decay (final hole in the 3s | ||
| + | -- state, the M1 edge) | ||
| + | |||
| + | -- we minimize the output by setting the verbosity to 0 | ||
| + | Verbosity(0) | ||
| + | |||
| + | -- We need the 2p, 3s and 3d shell in the basis | ||
| + | NF=18 | ||
| + | NB=0 | ||
| + | IndexDn_2p={0, | ||
| + | IndexUp_2p={1, | ||
| + | IndexDn_3s={6} | ||
| + | IndexUp_3s={7} | ||
| + | IndexDn_3d={8, | ||
| + | IndexUp_3d={9, | ||
| + | |||
| + | -- just like in the previous example we define several operators acting on the Ni -3d shell | ||
| + | |||
| + | 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(" | ||
| + | |||
| + | -- as in the previous example we define the Coulomb interaction | ||
| + | OppF0 =NewOperator(" | ||
| + | OppF2 =NewOperator(" | ||
| + | OppF4 =NewOperator(" | ||
| + | |||
| + | -- as in the previous example we define the crystal-field operator | ||
| + | Akm = PotentialExpandedOnClm(" | ||
| + | OpptenDq = NewOperator(" | ||
| + | |||
| + | -- and as in the previous example we define operators that count the number of eg and t2g | ||
| + | -- electrons | ||
| + | Akm = PotentialExpandedOnClm(" | ||
| + | OppNeg = NewOperator(" | ||
| + | Akm = PotentialExpandedOnClm(" | ||
| + | OppNt2g = NewOperator(" | ||
| + | |||
| + | -- In order te describe the resonance we need the interaction on the 2p shell (spin-orbit) | ||
| + | Oppcldots= NewOperator(" | ||
| + | |||
| + | -- and the Coulomb interaction between the 2p and 3d shell | ||
| + | OppUpdF0 = NewOperator(" | ||
| + | OppUpdF2 = NewOperator(" | ||
| + | OppUpdG1 = NewOperator(" | ||
| + | OppUpdG3 = NewOperator(" | ||
| + | |||
| + | -- and the Coulomb interaction between the 3s and 3d shell (new here, needed for the final state) | ||
| + | -- note that in the fluoressence yield example we did not need this interaction as the | ||
| + | -- spectra were integrated over all outgoing photon energies | ||
| + | OppUsdF0 = NewOperator(" | ||
| + | OppUsdG2 = NewOperator(" | ||
| + | |||
| + | -- next we define the dipole operator. The dipole operator is given as epsilon.r | ||
| + | -- with epsilon the polarization vector of the light and r the unit position vector | ||
| + | -- We can expand the position vector on (renormalized) spherical harmonics and use | ||
| + | -- the crystal-field operator to create the dipole operator. | ||
| + | |||
| + | -- we both define the dipole operator between the 2p and 3d shell as well as the dipole | ||
| + | -- operator between the 3s and 2p shell | ||
| + | |||
| + | -- x polarized light is defined as x = Cos[phi]Sin[theta] = sqrt(1/2) ( C_1^{(-1)} - C_1^{(1)}) | ||
| + | Akm = {{1, | ||
| + | TXASx = NewOperator(" | ||
| + | T3s2px = NewOperator(" | ||
| + | -- y polarized light is defined as y = Sin[phi]Sin[theta] = sqrt(1/2) I ( C_1^{(-1)} + C_1^{(1)}) | ||
| + | Akm = {{1, | ||
| + | TXASy = NewOperator(" | ||
| + | T3s2py = NewOperator(" | ||
| + | -- z polarized light is defined as z = Cos[theta] = C_1^{(0)} | ||
| + | Akm = {{1,0,1}} | ||
| + | TXASz = NewOperator(" | ||
| + | T3s2pz = NewOperator(" | ||
| + | |||
| + | -- besides linear polarized light one can define circular polarized light as the sum of | ||
| + | -- x and y polarizations with complex prefactors | ||
| + | TXASr = sqrt(1/ | ||
| + | TXASl =-sqrt(1/ | ||
| + | |||
| + | -- we can remove zero's from the dipole operator by chopping it. | ||
| + | TXASr.Chop() | ||
| + | TXASl.Chop() | ||
| + | |||
| + | -- once all operators are defined we can set some parameter values. | ||
| + | |||
| + | -- the value of U drops out of a crystal-field calculation as the total number of electrons | ||
| + | -- is always the same | ||
| + | U | ||
| + | -- F2 and F4 are often referred to in the literature as J_{Hund}. They represent the energy | ||
| + | -- differences between different multiplets. Numerical values can be found in the back of | ||
| + | -- my PhD. thesis for example. http:// | ||
| + | F2dd = 11.142 | ||
| + | F4dd = 6.874 | ||
| + | -- F0 is not the same as U, although they are related. Unimportant in crystal-field theory | ||
| + | -- the difference between U and F0 is so important that I do include it here. (U=0 so F0 is not) | ||
| + | F0dd = U+(F2dd+F4dd)*2/ | ||
| + | -- in crystal field theory U drops out of the equation, also true for the interaction between the | ||
| + | -- Ni 2p and Ni 3d electrons | ||
| + | Upd | ||
| + | -- The Slater integrals between the 2p and 3d shell, again the numerical values can be found | ||
| + | -- in the back of my PhD. thesis. (http:// | ||
| + | F2pd = 6.667 | ||
| + | G1pd = 4.922 | ||
| + | G3pd = 2.796 | ||
| + | -- F0 is not the same as U, although they are related. Unimportant in crystal-field theory | ||
| + | -- the difference between U and F0 is so important that I do include it here. (U=0 so F0 is not) | ||
| + | F0pd = Upd + G1pd*1/15 + G3pd*3/70 | ||
| + | -- We now also need the interaction between the 3s and 3d orbital | ||
| + | Usd | ||
| + | G2sd = 12.560 | ||
| + | -- and again F0 is not equal to U, so include the difference | ||
| + | F0sd = Usd + G2sd/10 | ||
| + | -- tenDq in NiO is 1.1 eV as can be seen in optics or using IXS to measure d-d excitations | ||
| + | tenDq | ||
| + | -- the Ni 3d spin-orbit is small but finite | ||
| + | zeta_3d = 0.081 | ||
| + | -- the Ni 2p spin-orbit is very large and should not be scaled as theory is quite accurate here | ||
| + | zeta_2p = 11.498 | ||
| + | -- we can add a small magnetic field, just to get nice expectation values. (units in eV... ) | ||
| + | Bz = 0.000001 | ||
| + | |||
| + | -- the total Hamiltonian is the sum of the different operators multiplied with their prefactor | ||
| + | Hamiltonian = F0dd*OppF0 + F2dd*OppF2 + F4dd*OppF4 + tenDq*OpptenDq + zeta_3d*Oppldots + Bz*(2*OppSz+OppLz) | ||
| + | |||
| + | -- We normally do not include core-valence interactions between filed and partially filled | ||
| + | -- shells as it drops out anyhow. For the XAS we thus need to define a " | ||
| + | -- (more complete) Hamiltonian. | ||
| + | XASHamiltonian = Hamiltonian + zeta_2p * Oppcldots + F0pd * OppUpdF0 + F2pd * OppUpdF2 + G1pd * OppUpdG1 + G3pd * OppUpdG3 | ||
| + | |||
| + | -- we now also need the Hamiltonian for the configuration with a core hole in the 3s state | ||
| + | -- here it is: | ||
| + | Hamiltonian3s = Hamiltonian + F0sd * OppUsdF0 + G2sd * OppUsdG2 | ||
| + | |||
| + | |||
| + | -- We saw in the previous example that NiO has a ground-state doublet with occupation | ||
| + | -- t2g^6 eg^2 and S=1 (S^2=S(S+1)=2). The next state is 1.1 eV higher in energy and thus | ||
| + | -- unimportant for the ground-state upto temperatures of 10 000 Kelvin. We thus restrict | ||
| + | -- the calculation to the lowest 3 eigenstates. | ||
| + | Npsi=3 | ||
| + | -- We need a filling of 6 electrons in the 2p shell | ||
| + | -- 2 electrons in the 3s shell | ||
| + | -- and 8 electrons in the 3d shell | ||
| + | StartRestrictions = {NF, NB, {" | ||
| + | |||
| + | -- And calculate the lowest 3 eigenfunctions | ||
| + | psiList = Eigensystem(Hamiltonian, | ||
| + | |||
| + | -- In order to get some information on these eigenstates it is good to plot expectation values | ||
| + | -- We first define a list of all the operators we would like to calculate the expectation value of | ||
| + | oppList={Hamiltonian, | ||
| + | |||
| + | -- next we loop over all operators and all states and print the expectation value | ||
| + | print(" | ||
| + | for i = 1,#psiList do | ||
| + | for j = 1,#oppList do | ||
| + | expectationvalue = Chop(psiList[i]*oppList[j]*psiList[i]) | ||
| + | io.write(string.format(" | ||
| + | end | ||
| + | io.write(" | ||
| + | end | ||
| + | |||
| + | -- spectra XAS | ||
| + | XASSpectra = CreateSpectra(XASHamiltonian, | ||
| + | XASSpectra.Print({{" | ||
| + | |||
| + | -- spectra FY | ||
| + | FYSpectra = CreateFluorescenceYield(XASHamiltonian, | ||
| + | FYSpectra.Print({{" | ||
| + | |||
| + | -- spectra RIXS | ||
| + | RIXSSpectra = CreateResonantSpectra(XASHamiltonian, | ||
| + | RIXSSpectra.Print({{" | ||
| + | |||
| + | print(" | ||
| + | |||
| + | -- and make some 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 out ' | ||
| + | set terminal postscript portrait enhanced color " | ||
| + | |||
| + | unset colorbox | ||
| + | |||
| + | energyshift=857.6 | ||
| + | energyshiftM1=110.8 | ||
| + | |||
| + | set multiplot | ||
| + | set size 0.5,0.55 | ||
| + | set origin 0,0 | ||
| + | |||
| + | set ylabel " | ||
| + | set xlabel " | ||
| + | |||
| + | set yrange [852:860] | ||
| + | set xrange [energyshiftM1-0.5: | ||
| + | |||
| + | plot "< | ||
| + | |||
| + | set origin 0.5,0 | ||
| + | |||
| + | set yrange [869:877] | ||
| + | set xrange [energyshiftM1-0.5: | ||
| + | |||
| + | plot "< | ||
| + | |||
| + | unset multiplot | ||
| + | |||
| + | set out ' | ||
| + | set size 1.0, 1.0 | ||
| + | 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 [energyshiftM1-0.5: | ||
| + | |||
| + | ofset = 0.25 | ||
| + | scale=10 | ||
| + | |||
| + | 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(" | ||
| + | </ | ||
| + | ### | ||
| + | |||
| + | ### | ||
| + | Just like in the case of $L_{2, | ||
| + | |||
| + | {{: | ||
| + | |||
| + | The first shows on the left the XAS spectra in black and the integrated RIXS spectra (FY) in blue. The core-core RIXS is shown on the right. | ||
| + | |||
| + | {{: | ||
| + | |||
| + | The second plot shows the same RIXS spectra, but now as an intensity map. | ||
| + | ### | ||
| + | |||
| + | ### | ||
| + | The output of the script is: | ||
| + | <file Quanty_Output RIXS_L23M1.out> | ||
| + | < | ||
| + | -2.721 | ||
| + | -2.721 | ||
| + | -2.721 | ||
| + | Start of LanczosTriDiagonalizeKrylovRR | ||
| + | Start of LanczosTriDiagonalizeKrylovRR | ||
| + | Finished calculating the spectra now start plotting. | ||
| + | This might take more time than the calculation | ||
| + | </ | ||
| + | ### | ||
| + | |||
| + | ===== Table of contents ===== | ||
| + | {{indexmenu> | ||