Differences
This shows you the differences between two versions of the page.
Next revision | Previous revision | ||
documentation:tutorials:nio_crystal_field:nixs_m45 [2016/10/08 21:23] – created Maurits W. Haverkort | documentation:tutorials:nio_crystal_field:nixs_m45 [2019/02/21 08:22] (current) – Maurits W. Haverkort | ||
---|---|---|---|
Line 1: | Line 1: | ||
+ | {{indexmenu_n> | ||
+ | ====== nIXS $M_{4,5}$ ($d$-$d$ excitations) ====== | ||
+ | ### | ||
+ | Inelastic x-ray scattering IXS (non-resonant) nIXS or x-ray Raman scattering allows one to measure non-dipolar allowed transitions. A powerful technique to look at even $d$-$d$ transitions with well defined selection rules \cite{Haverkort: | ||
+ | ### | ||
+ | |||
+ | ### | ||
+ | This tutorial compares calculated spectra to experiment. In order to make the plots you need to download the experimental data. You can download them in a zip file here {{ : | ||
+ | ### | ||
+ | |||
+ | ### | ||
+ | This example shows low energy $d$-$d$ transitions in NiO. The input is: | ||
+ | <code Quanty nIXS_M45.Quanty> | ||
+ | -- This example calculates the d-d excitations in NiO using non-resonant Inelastic X-ray | ||
+ | -- Scattering. This is one of the most beautiful spectroscopy techniques as the selection | ||
+ | -- rules are very " | ||
+ | |||
+ | -- We use the A^2 term of the interaction to make transitions between states with photons | ||
+ | -- of much higher energy. These photons now cary non negligible momentum and one can make | ||
+ | -- transitions beyond the dipole limit. | ||
+ | |||
+ | -- Here we look at k=2 and k=4 transitions between the Ni 3d orbitals | ||
+ | |||
+ | -- we set the output to a minimum | ||
+ | Verbosity(0) | ||
+ | |||
+ | -- define the basis of one particle spin-orbitals | ||
+ | -- we only need the d orbitals in this case | ||
+ | NF=10 | ||
+ | NB=0 | ||
+ | IndexDn_3d={0, | ||
+ | IndexUp_3d={1, | ||
+ | |||
+ | -- define operators on this basis | ||
+ | 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(" | ||
+ | |||
+ | -- define the crystal-field operator | ||
+ | |||
+ | Akm = PotentialExpandedOnClm(" | ||
+ | OpptenDq = NewOperator(" | ||
+ | |||
+ | -- define number operators counting the number of eg and t2g electrons | ||
+ | |||
+ | Akm = PotentialExpandedOnClm(" | ||
+ | OppNeg = NewOperator(" | ||
+ | Akm = PotentialExpandedOnClm(" | ||
+ | OppNt2g = NewOperator(" | ||
+ | |||
+ | -- set some parameters (see PRB 85, 165113 (2012) for more information) | ||
+ | beta = 0.8 | ||
+ | U | ||
+ | F2dd = 11.142 * beta | ||
+ | F4dd = 6.874 * beta | ||
+ | F0dd = U+(F2dd+F4dd)*2/ | ||
+ | tenDq = 1.100 | ||
+ | zeta_3d = 0.081 | ||
+ | Bz = 0.000001 | ||
+ | |||
+ | -- create a parameter dependent Hamiltonian | ||
+ | Hamiltonian = F0dd*OppF0 + F2dd*OppF2 + F4dd*OppF4 + tenDq*OpptenDq + zeta_3d*Oppldots + Bz*(2*OppSz + OppLz) | ||
+ | |||
+ | |||
+ | -- 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 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 | ||
+ | |||
+ | -- 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 | ||
+ | |||
+ | -- dd transitions from 3d (index 7 in Rnl) to 3d (index 7 in Rnl) | ||
+ | -- <R(r) | j_k(q r) | R(r)> | ||
+ | function RjRdd (q) | ||
+ | Rj0R = 0 | ||
+ | Rj2R = 0 | ||
+ | Rj4R = 0 | ||
+ | dr = Rnl[3][1]-Rnl[2][1] | ||
+ | r0 = Rnl[2][1]-2*dr | ||
+ | for ir = 2, #Rnl, 1 do | ||
+ | r = r0 + ir * dr | ||
+ | Rj0R = Rj0R + Rnl[ir][7] * SphericalBesselJ(0, | ||
+ | Rj2R = Rj2R + Rnl[ir][7] * SphericalBesselJ(2, | ||
+ | Rj4R = Rj4R + Rnl[ir][7] * SphericalBesselJ(4, | ||
+ | end | ||
+ | return Rj0R, Rj2R, Rj4R | ||
+ | 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_dd(q, theta, phi) | ||
+ | Rj0R, Rj2R, Rj4R = RjRdd(q) | ||
+ | k=0 | ||
+ | A0 = ExpandOnClm(k, | ||
+ | T0 = NewOperator(" | ||
+ | k=2 | ||
+ | A2 = ExpandOnClm(k, | ||
+ | T2 = NewOperator(" | ||
+ | k=4 | ||
+ | A4 = ExpandOnClm(k, | ||
+ | T4 = NewOperator(" | ||
+ | T = T0+T2+T4 | ||
+ | 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=4.5 | ||
+ | |||
+ | print(" | ||
+ | |||
+ | -- define some transition operators | ||
+ | qtheta=0 | ||
+ | qphi=0 | ||
+ | Tq001 = TnIXS_dd(q, | ||
+ | |||
+ | qtheta=Pi/2 | ||
+ | qphi=Pi/4 | ||
+ | Tq110 = TnIXS_dd(q, | ||
+ | |||
+ | qtheta=acos(sqrt(1/ | ||
+ | qphi=Pi/4 | ||
+ | Tq111 = TnIXS_dd(q, | ||
+ | |||
+ | qtheta=acos(sqrt(9/ | ||
+ | qphi=acos(sqrt(1/ | ||
+ | Tq123 = TnIXS_dd(q, | ||
+ | |||
+ | -- calculate the spectra | ||
+ | nIXSSpectra = CreateSpectra(Hamiltonian, | ||
+ | |||
+ | -- 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 " | ||
+ | |||
+ | set yrange [0:6.5] | ||
+ | |||
+ | 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: | ||
+ | |||
+ | {{: | ||
+ | ### | ||
+ | |||
+ | ### | ||
+ | We calculate the spectrum in 4 different directions of momentum transfer. The experimental spectra (by Verbeni (2009) and Huotari (2008) //et al.//) are measured on a powder sample and thus do not show the strong momentum direction dependence. In previous measurements (Larson //et al.// (2007) and later in great detail by Hiroaka //et al.// (2009) this angular dependence has been observed, which agrees well with the theoretical predictions. | ||
+ | ### | ||
+ | |||
+ | ### | ||
+ | For completeness the output of the script is: | ||
+ | <file Quanty_Output nIXS_M45.out> | ||
+ | < | ||
+ | -2.444 | ||
+ | -2.444 | ||
+ | -2.444 | ||
+ | for q= 4.5 per a0 ( 8.5037675605428 per A) The ratio of k=0, k=2 and k=4 transition strength is: | ||
+ | </ | ||
+ | ### | ||
+ | |||
+ | |||
+ | ===== Table of contents ===== | ||
+ | {{indexmenu> |