Dans_Diffraction.functions_scattering (version 0.9)
index
c:\users\grp66007\onedrive - diamond light source ltd\pythonprojects\dans_diffraction\dans_diffraction\functions_scattering.py

Module: functions_scattering.py
 
Functions:
intensity(structurefactor)
    Returns the squared structure factor
phase_factor(hkl, uvw)
    Return the complex phase factor:
phase_factor_qr(q, r)
    Return the complex phase factor:
scatteringbasis(q, azi_ref_q=(1, 0, 0), psi=0)
    Determine the scattering and polarisation vectors of a reflection based on energy, azimuth and polarisation.
scatteringcomponents(moment, q, azi_ref_q=(1, 0, 0), psi=0)
    Transform magnetic vector into components within the scattering plane
scatteringvectors(q, energy_kev, azi_ref_q=(1, 0, 0), psi=0, polarisation='s-p')
    Determine the scattering and polarisation vectors of a reflection based on energy, azimuth and polarisation.
sf_magnetic_neutron(q, r, occ, moment, magnetic_formfactor=None, debyewaller=False)
    Calculate the magnetic structure factor for the given HKL, using neutron magnetic form factor
sf_magnetic_neutron_polarised(q, r, occ, moment, incident_polarisation_vector=(1, 0, 0), magnetic_formfactor=None, debyewaller=False)
    Calculate the magnetic structure factor for the given HKL, using neutron magnetic form factor
sf_magnetic_xray(q, r, occ, moment, magnetic_formfactor=None, debyewaller=False)
    Calculate the non-resonant magnetic component of the structure factor
sf_magnetic_xray_beamline(q, r, occ, moment, energy_kev, magnetic_formfactor=None, debyewaller=False, azi_ref_q=(1, 0, 0), psi=0, polarisation='s-p')
    Calculate the non-resonant magnetic component of the structure factor
sf_magnetic_xray_polarised(q, r, occ, moment, incident_polarisation_vector=(1, 0, 0), magnetic_formfactor=None, debyewaller=False)
    Calculate the non-resonant magnetic component of the structure factor
sf_magnetic_xray_resonant(q, r, occ, moment, energy_kev, debyewaller=False, azi_ref_q=(1, 0, 0), psi=0, polarisation='s-p', f0=0, f1=1, f2=0)
    Calculate the non-resonant magnetic component of the structure factor
sf_magnetic_xray_resonant_alternate(q, r, occ, moment, energy_kev, debyewaller=False, polarisation='sp', azi_ref_q=(1, 0, 0), psi=0, f0=1, f1=1, f2=1)
    Calculate structure factors using resonant scattering factors in the dipolar approximation
structure_factor(scattering_factor, occupancy, debyewaller, phase)
    Return the complex structure factor:
xray_resonant_scattering_factor(q, moment, energy_kev, polarisation='sp', flm=(1, 1, 1), psi=0, azi_ref_q=(1, 0, 0))
    Calcualte fxres, the resonant x-ray scattering factor
 
 
By Dan Porter, PhD
Diamond
2018
 
Version 0.9
Last updated: 13/07/21
 
Version History:
11/11/18 0.1    Version History started.
13/07/21 0.9    Functions re-written and tested
 
@author: DGPorter

 
Modules
       
datetime
Dans_Diffraction.functions_crystallography
Dans_Diffraction.functions_general
numpy

 
Functions
       
autointensity(scattering_type, q, r, *args, **kwargs)
Choose a scattering type can calcuate the scattered intensity
:param scattering_type: named scattering function, see "get_scattering_function()"
:param q: array [n,3] reflection positions in A^-1
:param r: array [m,3] atomic positions in A
:param args: additional arguments to pass to choosen scattering function
:param kwargs: named arguments to pass to choosen scattering function
:return: float array [n]
autostructurefactor(scattering_type, q, r, *args, **kwargs)
Choose a scattering type can calcuate the structure factor
:param scattering_type:
:param q: array [n,3] reflection positions in A^-1
:param r: array [m,3] atomic positions in A
:param args: additional arguments to pass to choosen scattering function
:param kwargs: named arguments to pass to choosen scattering function
:return: complex array [n]
get_scattering_function(scattering_type)
Return function for given scattering type
Function will return structure factor function
:param scattering_type: str : scattering name as defined in SCATTERING_NAMES
:return: function
intensity(structurefactor)
Returns the squared structure factor
:param structurefactor: complex array [n] structure factor
:return: array [n]
options(occ=None, debyewaller=None, scattering_factor=None, moment=None, incident_polarisation_vector=(1, 0, 0), magnetic_formfactor=None, energy_kev=8, polarisation='sp', azi_ref_q=(1, 0, 0), psi=0, f0=0, f1=1, f2=0)
Create an input dict that will work with all structure factor (sf_) functions
:param occ: [m,1] array of atomic occupancies
:param debyewaller: [n,m] array of thermal factors for each atom and reflection
:param scattering_factor: array [n,m] or [n]: radiation dependent scattering factor/ form factor,/ scattering length
:param moment: [m,3] array of magnetic moment direction in orthogonal basis
:param incident_polarisation_vector: [1,3] direction of incident polarisation
:param magnetic_formfactor: [n,m] array of magnetic form factors for each atom and relection
:param energy_kev: float value of incident x-ray energy in keV
:param azi_ref_q: [1,3] azimuthal refence, in cartesian basis (Q)
:param psi: float value of the azimthal angle - the rotation out of the scattering plane.
:param polarisation: str definition of the polarisation can be: ['ss','sp','ps','pp'] with 's'=sigma, 'p'=pi
:param f0: float Flm value 0 (charge)
:param f1: float Flm value 1
:param f2: float Flm value 2
:return: dict
phase_factor(hkl, uvw)
Return the complex phase factor:
    phase_factor = exp(i.2.pi.HKL.UVW')
:param hkl: array [n,3] integer reflections
:param uvw: array [m,3] atomic positions in atomic basis units
:return: complex array [n,m]
phase_factor_qr(q, r)
Return the complex phase factor:
    phase_factor = exp(i.Q.R')
:param q: array [n,3] reflection positions in A^-1
:param r: array [m,3] atomic positions in A
:return: complex array [n,m]
scatteringbasis(q, azi_ref_q=(1, 0, 0), psi=0)
Determine the scattering and polarisation vectors of a reflection based on energy, azimuth and polarisation.
:param q: [1*3] reflection vector in a cartesian basis
:param azi_ref_q: [1,3] direction along which the azimuthal zero angle is determind
:param psi: float azimuthal angle about U3 in degrees
:return: U1, U2, U3
The basis is chosen such that Q defines the scattering plane, the sigma direction is normal to this plane,
the pi direction is always within this plane.
The azimuthal angle defines a rotation about the Q axis in a clockwise mannor, matching I16.
At an azimuth of 0degrees, U1 is perpendicular to Q, along the direction of azim_zero.
scatteringcomponents(moment, q, azi_ref_q=(1, 0, 0), psi=0)
Transform magnetic vector into components within the scattering plane
:param moment: [n*3] array of magnetic moments in a cartesian basis
:param q: [1*3] reflection vector in a cartesian basis
:param azi_ref_q: [1*3] azimuthal reference in a cartesian basis
:param psi: float azimuthal angle
:return: (z1, z2, z3) components of the magnetic moment along the reflection direction
scatteringvectors(q, energy_kev, azi_ref_q=(1, 0, 0), psi=0, polarisation='s-p')
Determine the scattering and polarisation vectors of a reflection based on energy, azimuth and polarisation.
:param q: [n,3] reflection vector in a cartesian basis
:param energy_kev: x-ray scattering energy in keV
:param azi_ref_q: [1,3] direction along which the azimuthal zero angle is determind
:param psi: float angle in degrees about the azimuth
:param polarisation: polarisation with respect to the scattering plane, options:
            'ss' : sigma-sigma polarisation
            'sp' : sigma-pi polarisation
            'ps' : pi-sigma polarisation
            'pp' : pi-pi polarisation
:return: kin, kout, ein, eout
Returned values are [n,3] arrays
    kin : [n,3] array of incident wavevectors
    kout: [n,3] array of scattered wavevectors
    ein : [n,3] array of incident polarisation
    eout: [n,3] array of scattered polarisation
 
The basis is chosen such that Q defines the scattering plane, sigma and pi directions are normal to this plane.
sf_atom(q, r, scattering_factor=None, occ=None, debyewaller=None, **kwargs)
:param q: array [n,3] reflection positions in A^-1
:param r: array [m,3] atomic positions in A
:param scattering_factor: array [n,m] or [n]: radiation dependent scattering factor/ form factor,/ scattering length
:param occ: array [m]: occupancy of each atom
:param debyewaller: array [n,m]: thermal vibration factor of each atom and reflection
:param kwargs: additional options[*unused]
:return: complex array [n]
sf_magnetic_neutron(q, r, moment, magnetic_formfactor=None, occ=None, debyewaller=None, **kwargs)
Calculate the magnetic structure factor for the given HKL, using neutron magnetic form factor
Assumes an unpolarised incident beam.
    Reference: G. L. Squires, Introduction to the Theory of Thermal Neutron Scattering (Cambridge University Press, 1997).
:param q: [n,3] array of reflections in cartesian coordinate system in units of inverse-A
:param r: [m,3] array of atomic positions in A
:param moment: [m,3] array of magnetic moment direction in orthogonal basis
:param magnetic_formfactor: [n,m] array of magnetic form factors for each atom and relection
:param occ: [m,1] array of atomic occupancies
:param debyewaller: [n,m] array of thermal factors for each atom and reflection, or False to omit
:param kwargs: additional options[*unused]
:return sf: [n] complex array of structure factors
sf_magnetic_neutron_polarised(q, r, moment, incident_polarisation_vector=(1, 0, 0), magnetic_formfactor=None, occ=None, debyewaller=None, **kwargs)
Calculate the magnetic structure factor for the given HKL, using neutron magnetic form factor
    Reference: G. L. Squires, Introduction to the Theory of Thermal Neutron Scattering (Cambridge University Press, 1997).
:param q: [n,3] array of hkl reflections
:param r: [m,3] array of atomic positions in r.l.u.
:param moment: [m,3] array of magnetic moment direction in orthogonal basis
:param incident_polarisation_vector: [1,3] direction of incident polarisation
:param magnetic_formfactor: [n,m] array of magnetic form factors for each atom and relection
:param occ: [m,1] array of atomic occupancies
:param debyewaller: [n,m] array of thermal factors for each atom and reflection, or False to omit
:param kwargs: additional options[*unused]
:return sf: [n] complex array of structure factors
sf_magnetic_xray(q, r, moment, magnetic_formfactor=None, occ=None, debyewaller=None, **kwargs)
Calculate the non-resonant magnetic component of the structure factor
:param q: [n,3] array of hkl reflections
:param r: [m,3] array of atomic positions in r.l.u.
:param moment: [m,3] array of magnetic moment direction in orthogonal basis
:param magnetic_formfactor: [n,m] array of magnetic form factors for each atom and relection
:param occ: [m,1] array of atomic occupancies
:param debyewaller: [n,m] array of thermal factors for each atom and reflection
:param kwargs: additional options[*unused]
:return sf: [n] complex array of structure factors
 
    f_non-res_mag = i.r0.(hw/mc^2).fD.[.5.L.A + S.B]
    B = e_o X e_i + (k_o X e_o) * k_o.e_i - (k_i X e_i) * k_i.e_o - (k_o X e_o) X (k_i X e_i)
- ignore orbital moment L
- fD = magnetic form factor
- S = spin moment
- k_i, k_o = wavevector in, out
- e_i, e_o = polarisation in, out
From Hill+McMorrow Acta Cryst. 1996 A52, 236-244 Equ. (2)
Book: "X-ray Scattering and Absorption by Magnetic Materials" by Loevesy and Collins. Ch 2. Eqn.2.21+1
No orbital component assumed
magnetic moments assumed to be in the same reference frame as the polarisation
sf_magnetic_xray_beamline(q, r, moment, energy_kev, magnetic_formfactor=None, occ=None, debyewaller=None, azi_ref_q=(1, 0, 0), psi=0, polarisation='s-p', **kwargs)
Calculate the non-resonant magnetic component of the structure factor
:param q: [n,3] array of hkl reflections
:param r: [m,3] array of atomic positions in r.l.u.
:param moment: [m,3] array of magnetic moment direction in orthogonal basis
:param energy_kev: float value of incident x-ray energy in keV
:param magnetic_formfactor: [n,m] array of magnetic form factors for each atom and relection
:param occ: [m,1] array of atomic occupancies
:param debyewaller: [n,m] array of thermal factors for each atom and reflection
:param azi_ref_q: [1,3] azimuthal refence, in cartesian basis (Q)
:param psi: [p] array of azimthal angles - the rotation out of the scattering plane.
:param polarisation: str definition of the polarisation can be: ['ss','sp','ps','pp'] with 's'=sigma, 'p'=pi
:param kwargs: additional options[*unused]
:return sf: [n, p] complex array of structure factors for different reflections and azimuths
 
    f_non-res_mag = i.r0.(hw/mc^2).fD.[.5.L.A + S.B]
    B = e_o X e_i + (k_o X e_o) * k_o.e_i - (k_i X e_i) * k_i.e_o - (k_o X e_o) X (k_i X e_i)
- ignore orbital moment L
- fD = magnetic form factor
- S = spin moment
- k_i, k_o = wavevector in, out
- e_i, e_o = polarisation in, out
 
From Hill+McMorrow Acta Cryst. 1996 A52, 236-244 Equ. (2)
Book: "X-ray Scattering and Absorption by Magnetic Materials" by Loevesy and Collins. Ch 2. Eqn.2.21+1
No orbital component assumed
magnetic moments assumed to be in the same reference frame as the polarisation
sf_magnetic_xray_polarised(q, r, moment, incident_polarisation_vector=(1, 0, 0), magnetic_formfactor=None, occ=None, debyewaller=None, **kwargs)
Calculate the non-resonant magnetic component of the structure factor
:param q: [n,3] array of hkl reflections
:param r: [m,3] array of atomic positions in r.l.u.
:param moment: [m,3] array of magnetic moment direction in orthogonal basis
:param incident_polarisation_vector: [1,3] direction of incident polarisation
:param magnetic_formfactor: [n,m] array of magnetic form factors for each atom and relection
:param occ: [m,1] array of atomic occupancies
:param debyewaller: [n,m] array of thermal factors for each atom and reflection
:param kwargs: additional options[*unused]
:return sf: [n] complex array of structure factors
 
    f_non-res_mag = i.r0.(hw/mc^2).fD.[.5.L.A + S.B]
    B = e_o X e_i + (k_o X e_o) * k_o.e_i - (k_i X e_i) * k_i.e_o - (k_o X e_o) X (k_i X e_i)
- ignore orbital moment L
- fD = magnetic form factor
- S = spin moment
- k_i, k_o = wavevector in, out
- e_i, e_o = polarisation in, out
From Hill+McMorrow Acta Cryst. 1996 A52, 236-244 Equ. (2)
Book: "X-ray Scattering and Absorption by Magnetic Materials" by Loevesy and Collins. Ch 2. Eqn.2.21+1
No orbital component assumed
magnetic moments assumed to be in the same reference frame as the polarisation
sf_magnetic_xray_resonant(q, r, moment, energy_kev, occ=None, debyewaller=None, azi_ref_q=(1, 0, 0), psi=0, polarisation='sp', f0=0, f1=1, f2=0, **kwargs)
Calculate the non-resonant magnetic component of the structure factor
:param q: [n,3] array of hkl reflections
:param r: [m,3] array of atomic positions in r.l.u.
:param moment: [m,3] array of magnetic moment direction in orthogonal basis
:param energy_kev: float value of incident x-ray energy in keV
:param occ: [m,1] array of atomic occupancies
:param debyewaller: [n,m] array of thermal factors for each atom and reflection
:param azi_ref_q: [1,3] azimuthal refence, in cartesian basis (Q)
:param psi: [p] array of azimthal angles - the rotation out of the scattering plane.
:param polarisation: str definition of the polarisation can be: ['ss','sp','ps','pp'] with 's'=sigma, 'p'=pi
:param f0: float Flm value 0 (charge)
:param f1: float Flm value 1
:param f2: float Flm value 2
:param kwargs: additional options[*unused]
:return sf: [n, p] complex array of structure factors for different reflections and azimuths
 
f_res_mag = [(e'.e)F0 - i(e'xe).Z*F1 + (e'.Z)*(e.Z)*F2]
 
From Hill+McMorrow Acta Cryst. 1996 A52, 236-244 Equ. (2)
Book: "X-ray Scattering and Absorption by Magnetic Materials" by Loevesy and Collins. Ch 2. Eqn.2.21+1
No orbital component assumed
magnetic moments assumed to be in the same reference frame as the polarisation
sf_magnetic_xray_resonant_alternate(q, r, moment, energy_kev, occ=None, debyewaller=None, polarisation='sp', azi_ref_q=(1, 0, 0), psi=0, f0=0, f1=1, f2=0, **kwargs)
Calculate structure factors using resonant scattering factors in the dipolar approximation
:param q: [n,3] array of hkl reflections
:param r: [m,3] array of atomic positions in r.l.u.
:param moment: [m,3] array of magnetic moment direction in orthogonal basis
:param energy_kev: float value of incident x-ray energy in keV
:param occ: [m,1] array of atomic occupancies
:param debyewaller: [n,m] array of thermal factors for each atom and reflection
:param azi_ref_q: [1,3] azimuthal refence, in cartesian basis (Q)
:param psi: [p] array of azimthal angles - the rotation out of the scattering plane.
:param polarisation: str definition of the polarisation can be: ['ss','sp','ps','pp'] with 's'=sigma, 'p'=pi
:param f0: float Flm value 0 (charge)
:param f1: float Flm value 1
:param f2: float Flm value 2
:param kwargs: additional options[*unused]
:return sf: [n, p] complex array of structure factors for different reflections and azimuths
 
  I = Scattering.xray_resonant(HKL,energy_kev,polarisation,F0,F1,F2)
Returns an array with the same length as HKL, giving the real intensity at each reflection.
    energy_kev = x-ray energy in keV
    polarisation = x-ray polarisation: 'ss',{'sp'},'ps','pp'
    f0/1/2 = Resonance factor Flm
    azim_zero = [h,k,l] vector parallel to the origin of the azimuth
    psi = azimuthal angle defining the scattering plane
 
Uses the E1E1 resonant x-ray scattering amplitude:
    fxr_n = (ef.ei)*f0 -i(ef X ei).z_n*f1 + (ef.z_n)(ei.z_n)f2
 
Where ei and ef are the initial and final polarisation states, respectively,
and z_n is a unit vector in the direction of the magnetic moment of the nth ion.
The polarisation states are determined to be one of the natural synchrotron
states, where sigma (s) is perpendicular to the scattering plane and pi (p) is
parallel to it.
        ( s-s  s-p )
        ( p-s  p-p )
 
From Hill+McMorrow Acta Cryst. 1996 A52, 236-244 Equ. (15)
sf_xray_dispersion(q, r, scattering_factor, occ=None, debyewaller=None, **kwargs)
Calculate the resonant x-ray structure factor
:param q: [n,3] array of hkl reflections
:param r: [m,3] array of atomic positions in r.l.u.
:param scattering_factor: array [n,m,e]: energy dependent complex atomic form factor
:param occ: [m,1] array of atomic occupancies
:param debyewaller: [n,m] array of thermal factors for each atom and reflection
:param kwargs: additional options[*unused]
:return sf: [n, e] complex array of structure factors
structure_factor(scattering_factor, occupancy, debyewaller, phase)
Return the complex structure factor:
    structure_factor = sum_i( f.occ.dw.phase )
:param scattering_factor: array [n,m] or [n]: radiation dependent scattering factor/ form factor,/ scattering length
:param occupancy: array [m]: occupancy of each atom
:param debyewaller: array [n,m]: thermal vibration factor of each atom and reflection
:param phase: complex array [n,m]: complex phase factor exp(-i.Q.R)
:return: complex array [n]
xray_resonant_scattering_factor(q, moment, energy_kev, polarisation='sp', flm=(1, 1, 1), psi=0, azi_ref_q=(1, 0, 0))
Calcualte fxres, the resonant x-ray scattering factor
  fxres = Scattering.xray_resonant_scattering_factor(HKL,energy_kev,polarisation,flm,azim_zero,psi)
energy_kev = x-ray energy in keV
    polarisation = x-ray polarisation: 'ss',{'sp'},'ps','pp'
    flm = (f0, f1, f2) Resonance factor Flm, f0/1/2 should be 0 or 1 each
    azim_zero = [h,k,l] vector parallel to the origin of the azimuth {[1,0,0]}
    psi = azimuthal angle defining the scattering plane {0}
 
:param q: [n*3] array of reflection coordinates in cartesian basis (Q)
:param moment: [mx3] array of magnetic moments in cartesian basis
:param energy_kev: float energy in keV
:param polarisation: polarisation condition: 'sp', 'ss', 'ps', 'pp'. s=sigma, p=pi
:param flm: (f0, f1, f2) Resonance factor Flm, f0/1/2 should be 0 or 1 each
:param azi_ref_q: azimuthal refence, in cartesian basis (Q)
:param psi: float, azimuthal angle
:return: fxres [n*1] array of resonant x-ray scattering factors
 
Uses the E1E1 resonant x-ray scattering amplitude:
    fxr_n = (ef.ei)*F0 -i(ef X ei).z_n*F1 + (ef.z_n)(ei.z_n)F2
 
Where ei and ef are the initial and final polarisation states, respectively,
and z_n is a unit vector in the direction of the magnetic moment of the nth ion.
The polarisation states are determined to be one of the natural synchrotron
states, where sigma (s) is perpendicular to the scattering plane and pi (p) is
parallel to it.
        ( s-s  s-p )
        ( p-s  p-p )
 
From Hill+McMorrow Acta Cryst. 1996 A52, 236-244 Equ. (15)

 
Data
        DEBUG_MODE = False
MAX_QR_ARRAY = 10000000.0
SCATTERING_TYPES = {'electron': ['electron', 'ele', 'e'], 'neutron': ['neutron', 'n', 'nuclear'], 'neutron magnetic': ['neutron magnetic', 'magnetic neutron', 'magnetic'], 'xray': ['xray', 'x', 'x-ray', 'thomson', 'charge'], 'xray dispersion': ['dispersion', 'xray dispersion'], 'xray fast': ['xray fast', 'xfast'], 'xray magnetic': ['xray magnetic', 'magnetic xray', 'spin xray', 'xray spin'], 'xray resonant': ['xray resonant', 'resonant', 'resonant xray', 'rxs']}
TIME_REPORT = True