Differences
This shows you the differences between two versions of the page.
Next revision | Previous revision | ||
documentation:tutorials:nio_crystal_field:fy_l23m45 [2016/10/08 19:14] – created Maurits W. Haverkort | documentation:tutorials:nio_crystal_field:fy_l23m45 [2018/04/26 15:21] (current) – Maurits W. Haverkort | ||
---|---|---|---|
Line 1: | Line 1: | ||
+ | {{indexmenu_n> | ||
+ | ====== FY $L_{2, | ||
+ | ### | ||
+ | The absorption cross section is in principle measured using transmission. Transmission experiments in the soft-x-ray regime can be difficult as the absorption is quite high. Alternatively one can measure the reflectivity, | ||
+ | ### | ||
+ | |||
+ | ### | ||
+ | An alternative is to measure the fluorescence yield. Although not proportional to the absorption cross section \cite{Kurian: | ||
+ | ### | ||
+ | |||
+ | ### | ||
+ | In the following example we calculate the excitation of a $2p$ electron into the $3d$ shell of Ni in NiO. ($L_{2,3}$ edge). We integrate over the decay of a $3d$ electron into the $2p$ orbital (removing an electron from the $3d$-shell i.e. $M_{4,5}$) We thus look at the $L_{2, | ||
+ | ### | ||
+ | |||
+ | ### | ||
+ | The input is: | ||
+ | <code Quanty FY_L23M45.Quanty> | ||
+ | -- In this example we will calculate the fluorescence yield spectra | ||
+ | -- One makes an excitation from 2p to 3d and then looks at a specific decay channel | ||
+ | -- Or at the sum over all channels | ||
+ | -- the spectra are integrated over the energy of the decay channel which allows for | ||
+ | -- extreme efficient calculation of these spectra. | ||
+ | -- Note that most detectors will not be equally sensitive to all possible photon energies | ||
+ | -- and one thus would always measure some weighted sum over the different decay channels | ||
+ | |||
+ | -- this file calculates the Ni L23M45 spectra. | ||
+ | -- (L23, i.e. we excite from 2p to 3d) | ||
+ | -- (M45, we decay from the 3d shell, into the 2p shell) | ||
+ | |||
+ | -- we minimize the output by setting the verbosity to 0 | ||
+ | Verbosity(0) | ||
+ | |||
+ | -- In order to do crystal-field theory for NiO we need to define a Ni d-shell. | ||
+ | -- A d-shell has 10 elements and we label again the even spin-orbitals to be spin down | ||
+ | -- and the odd spin-orbitals to be spin up. In order to calculate 2p to 3d excitations we | ||
+ | -- also need a Ni 2p shell. We thus have a total of 10+6=16 fermions, 6 Ni-2p and 10 Ni-3d | ||
+ | -- spin-orbitals | ||
+ | NF=16 | ||
+ | NB=0 | ||
+ | IndexDn_2p={0, | ||
+ | IndexUp_2p={1, | ||
+ | IndexDn_3d={6, | ||
+ | IndexUp_3d={7, | ||
+ | |||
+ | -- 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(" | ||
+ | |||
+ | -- new for core level spectroscopy are operators that define the interaction acting on the | ||
+ | -- Ni-2p shell. There is actually only one of these interactions, | ||
+ | -- spin-orbit interaction | ||
+ | |||
+ | Oppcldots= NewOperator(" | ||
+ | |||
+ | -- we also need to define the Coulomb interaction between the Ni 2p- and Ni 3d-shell | ||
+ | -- Again the interaction (e^2/ | ||
+ | -- between two shells we need to consider two cases. For the direct interaction a 2p electron | ||
+ | -- scatters of a 3d electron into a 2p and 3d electron. The radial integrals involve | ||
+ | -- the square of a 2p radial wave function at coordinate 1 and the square of a 3d radial | ||
+ | -- wave function at coordinate 2. The transfer of angular momentum can either be 0 or 2. | ||
+ | -- These processes are called direct and the resulting Slater integrals are F[0] and F[2]. | ||
+ | -- The second proces involves a 2p electron scattering of a 3d electron into the 3d shell | ||
+ | -- and at the same time the 3d electron scattering into a 2p shell. These exchange processes | ||
+ | -- involve radial integrals over the product of a 2p and 3d radial wave function. The transfer | ||
+ | -- of angular momentum in this case can be 1 or 3 and the Slater integrals are called G1 and G3. | ||
+ | |||
+ | -- In Quanty you can enter these processes by labeling 4 indices for the orbitals, once | ||
+ | -- the 2p shell with spin up, 2p shell with spin down, 3d shell with spin up and 3d shell with | ||
+ | -- spin down. Followed by the direct Slater integrals (F0 and F2) and the exchange Slater | ||
+ | -- integrals (G1 and G3) | ||
+ | |||
+ | -- Here we define the operators separately and later sum them with appropriate prefactors | ||
+ | |||
+ | OppUpdF0 = NewOperator(" | ||
+ | OppUpdF2 = NewOperator(" | ||
+ | OppUpdG1 = NewOperator(" | ||
+ | OppUpdG3 = 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. | ||
+ | |||
+ | -- x polarized light is defined as x = Cos[phi]Sin[theta] = sqrt(1/2) ( C_1^{(-1)} - C_1^{(1)}) | ||
+ | Akm = {{1, | ||
+ | TXASx = 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(" | ||
+ | -- z polarized light is defined as z = Cos[theta] = C_1^{(0)} | ||
+ | Akm = {{1,0,1}} | ||
+ | TXASz = 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() | ||
+ | |||
+ | -- the 3d to 2p dipole transition is the conjugate transpose of the 2p to 3d dipole transition | ||
+ | TXASxdag = ConjugateTranspose(TXASx) | ||
+ | TXASydag = ConjugateTranspose(TXASy) | ||
+ | TXASzdag = ConjugateTranspose(TXASz) | ||
+ | TXASldag = ConjugateTranspose(TXASl) | ||
+ | TXASrdag = ConjugateTranspose(TXASr) | ||
+ | |||
+ | -- 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 | ||
+ | -- 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 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 | ||
+ | -- in order to make sure we have a filling of 8 | ||
+ | -- electrons we need to define some restrictions | ||
+ | -- We need to restrict the occupation of the Ni-2p shell to 6 for the ground state and for | ||
+ | -- the Ni 3d-shell to 8. | ||
+ | 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 | ||
+ | |||
+ | -- here we calculate the x-ray absorption spectra, not the main task of this file, but nice to do so we can compare | ||
+ | XASSpectra = CreateSpectra(XASHamiltonian, | ||
+ | XASSpectra.Print({{" | ||
+ | |||
+ | -- and we calculate the FY spectra | ||
+ | FYSpectra = CreateFluorescenceYield(XASHamiltonian, | ||
+ | FYSpectra.Print({{" | ||
+ | |||
+ | -- in order to plot both the XAS and FY spectra we can define a gnuplot script | ||
+ | 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 xlabel "E (eV)" font " | ||
+ | set ylabel " | ||
+ | |||
+ | set out ' | ||
+ | set terminal postscript portrait enhanced color " | ||
+ | set yrange [0:0.6] | ||
+ | |||
+ | set multiplot layout 3, 3 | ||
+ | |||
+ | plot " | ||
+ | " | ||
+ | " | ||
+ | " | ||
+ | plot " | ||
+ | " | ||
+ | " | ||
+ | " | ||
+ | plot " | ||
+ | " | ||
+ | " | ||
+ | " | ||
+ | |||
+ | plot " | ||
+ | " | ||
+ | " | ||
+ | " | ||
+ | plot " | ||
+ | " | ||
+ | " | ||
+ | " | ||
+ | plot " | ||
+ | " | ||
+ | " | ||
+ | " | ||
+ | |||
+ | plot " | ||
+ | " | ||
+ | " | ||
+ | " | ||
+ | plot " | ||
+ | " | ||
+ | " | ||
+ | " | ||
+ | plot " | ||
+ | " | ||
+ | " | ||
+ | " | ||
+ | |||
+ | unset multiplot | ||
+ | ]] | ||
+ | |||
+ | |||
+ | -- write the gnuplot script to a file | ||
+ | file = io.open(" | ||
+ | file: | ||
+ | file: | ||
+ | |||
+ | -- call gnuplot to execute the script | ||
+ | os.execute(" | ||
+ | -- and transform the ps to pdf | ||
+ | os.execute(" | ||
+ | </ | ||
+ | ### | ||
+ | |||
+ | ### | ||
+ | The script returns 9 plots with each 4 curves. The local ground-state of Ni in NiO is 3-fold degenerate in the paramagnetic phase ($S=1$) The different columns show the spectra for the states with different $S_z$. In the paramagnetic phase one should summ these 3 spectra, in a full magnetized sample one measurers either the left or the right column. The different rows the different incoming polarization. Top row z-polarized, | ||
+ | |||
+ | {{: | ||
+ | ### | ||
+ | |||
+ | ### | ||
+ | The script shows some information on the ground-state, | ||
+ | <file Quanty_Output FY_L23M45.out> | ||
+ | < | ||
+ | -2.721 | ||
+ | -2.721 | ||
+ | -2.721 | ||
+ | Start of LanczosTriDiagonalizeKrylovRR | ||
+ | Start of LanczosTriDiagonalizeKrylovRR | ||
+ | Start of LanczosTriDiagonalizeKrylovRR | ||
+ | Start of LanczosTriDiagonalizeKrylovRR | ||
+ | Start of LanczosTriDiagonalizeKrylovRR | ||
+ | Start of LanczosTriDiagonalizeKrylovRR | ||
+ | Start of LanczosTriDiagonalizeKrylovRR | ||
+ | Start of LanczosTriDiagonalizeKrylovRR | ||
+ | Start of LanczosTriDiagonalizeKrylovRR | ||
+ | </ | ||
+ | ### | ||
+ | |||
+ | ===== Table of contents ===== | ||
+ | {{indexmenu> |