\(\renewcommand\AA{\unicode{x212B}}\)
DiscusMultipleScatteringCorrection v1¶
Summary¶
Calculates a multiple scattering correction using a Monte Carlo method
See Also¶
MayersSampleCorrection, CarpenterSampleCorrection, VesuvioCalculateMS
This algorithm is also known as: Muscat
Properties¶
Name |
Direction |
Type |
Default |
Description |
---|---|---|---|---|
InputWorkspace |
Input |
Mandatory |
The name of the input workspace. The input workspace must have X units of Momentum (k) for elastic calculations and units of energy transfer (DeltaE) for inelastic calculations. This is used to supply the sample details, the detector positions and the x axis range to calculate corrections for |
|
StructureFactorWorkspace |
Input |
Mandatory |
The name of the workspace containing S’(q) or S’(q, w). For elastic calculations, the input workspace must contain a single spectrum and have X units of momentum transfer. A workspace group containing one workspace per component can also be supplied if a calculation is being run on a workspace with a sample environment specified |
|
OutputWorkspace |
Output |
WorkspaceGroup |
Mandatory |
Name for the WorkspaceGroup that will be created. Each workspace in the group contains a calculated weight for a particular number of scattering events. The number of scattering events varies from 1 up to the number supplied in the NumberOfScatterings parameter. The group will also include an additional workspace for a calculation with a single scattering event where the absorption post scattering has been set to zero |
ScatteringCrossSection |
Input |
A workspace containing the scattering cross section as a function of k, \(\sigma_s(k)\). Note - this parameter would normally be left empty which results in the tabulated cross section data being used instead which implies no wavelength dependence |
||
NumberOfSimulationPoints |
Input |
number |
Optional |
The number of points on the input workspace x axis for which a simulation is attempted |
NeutronPathsSingle |
Input |
number |
1000 |
The number of “neutron” paths to generate for single scattering |
NeutronPathsMultiple |
Input |
number |
1000 |
The number of “neutron” paths to generate for multiple scattering |
SeedValue |
Input |
number |
123456789 |
Seed the random number generator with this value |
NumberScatterings |
Input |
number |
2 |
Number of scatterings |
Interpolation |
Input |
string |
Linear |
Method of interpolation used to compute unsimulated values. Allowed values: [‘Linear’, ‘CSpline’] |
SparseInstrument |
Input |
boolean |
False |
Enable simulation on special instrument with a sparse grid of detectors interpolating the results to the real instrument. |
NumberOfDetectorRows |
Input |
number |
5 |
Number of detector rows in the detector grid of the sparse instrument. |
NumberOfDetectorColumns |
Input |
number |
10 |
Number of detector columns in the detector grid of the sparse instrument. |
ImportanceSampling |
Input |
boolean |
False |
Enable importance sampling on the Q value chosen on multiple scatters based on Q.S(Q) |
MaxScatterPtAttempts |
Input |
number |
5000 |
Maximum number of tries made to generate a scattering point within the sample. Objects with holes in them, e.g. a thin annulus can cause problems if this number is too low. If a scattering point cannot be generated by increasing this value then there is most likely a problem with the sample geometry. |
SimulateEnergiesIndependently |
Input |
boolean |
False |
For inelastic calculation, whether the results for adjacent energy transfer bins are simulated separately. Currently applies to Direct geometry only |
NormalizeStructureFactors |
Input |
boolean |
False |
Enable normalization of supplied structure factor(s). May be required when running a calculation involving more than one material where the normalization of the default S(Q)=1 structure factor doesn’t match the normalization of a supplied non-isotropic structure factor |
Description¶
This algorithm calculates a Multiple Scattering correction using a Monte Carlo integration method. The method uses a structure function for the sample to determine the probability of a particular momentum transfer (q) and energy transfer (\(\omega\)) value for each scattering event and it doesn’t therefore rely on an assumption that the scattering is isotropic.
The structure function that the algorithm takes as input is a linear combination of the coherent and incoherent structure factors:
\(S'(Q, \omega) = \frac{1}{\sigma_b}(\sigma_{coh} S(Q, \omega) + \sigma_{inc} S_s(Q, \omega))\)
If the sample is a perfectly coherent scatterer then \(S'(Q, \omega) = S(Q, \omega)\)
The algorithm is based on code which was originally written in Fortran as part of the Discus program [1]. The code was subsequently resurrected and improved by Spencer Howells under the Muscat name and was included in the QENS MODES package [2], [3] These original programs calculated multiple scattering corrections for inelastic instruments but an elastic diffraction version of the code was also created and results from that program are included in this paper by Mancinelli [4].
Theory¶
The theory is outlined here for an inelastic calculation. The calculation performed for an elastic instrument is a special case of this with \(\omega=0\).
The algorithm calculates a set of dimensionless weights \(J_n\) describing the probability of detection at an angle \(\theta\) after n scattering events given a total incident flux \(I_0\) and a transmitted flux of T:
\(T_n(\theta,k_{in}, \omega) = J_n(\theta,k_{in}, \omega) I_0(k_{in})\)
The quantity \(J_n\) is calculated by performing the following integration:
The variables \(l_i^{max}\) represent the maximum path length before the next scatter given a particular phi and theta value (Q). Each \(l_i\) is actually a function of all of the earlier values for the \(l_i\), \(\phi\), \(Q\) and \(\omega\) variables ie \(l_i = l_i(l_1, l_2, ..., l_{i-1}, \phi_1, \phi_2, ..., \phi_i, Q_1, Q_2, ..., Q_i, \omega_1, \omega_2, ..., \omega_i)\). The integration over the variables \(Q\) and \(\omega\) is done over the kinematically accessible region \(D(k_i)\). The integral is done over \(\omega\) first and the range is defined by the minimum value on the \(\omega\) axis in the \(S(Q, \omega)\) profile and a maximum which is equal to the total energy loss of the pre-scatter neutron prior to the ith scatter. The limits on the q integration are then calculated as follows (and are also a function of i). These formulae for the q limits reduce to 0 and 2k for elastic.
\(q_i^{min} = |k_{i+1} - k_i|\)
\(q_i^{max} = q_i^{min} + 2 min(k_i, k_{i+1})\)
The following substitutions are then performed in order to make it more convenient to evaluate as a Monte Carlo integral:
\(t_i = \frac{1-e^{-\mu_T l_i}}{1-e^{-\mu_T l_i^{max}}}\)
\(u_i = \frac{\phi_i}{2 \pi}\)
\(2 k_i^2 = \frac{\sigma_s}{\sigma_s(k_i)} I(k_i)\) where \(I(k) = \iint \limits_{D(k)} Q S(Q, \omega) dQ d\omega\)
Using the new variables the integral is:
This is evaluated as a Monte Carlo integration by selecting random values for the variables \(t_i\) and \(u_i\) between 0 and 1. The integral over \(Q\omega\) space is performed by integrating a slightly modified \(S(Q,\omega)\) function over a rectangular region. \(S_{kin}(Q,\omega)\) equals zero if \(Q\) and \(\omega\) are outside the kinematically accessible region. The rectangular region spans the full length of the \(\omega\) axis in the \(S(Q,\omega)\) profile and goes from zero to the maximum possible \(Q_i\) for a particular \(k_i\) in the q direction.
A simulated path is traced through the sample to enable the \(l_i^{\ max}\) values to be calculated. The path is traced by calculating the \(l_i\), \(\theta\) and \(\phi\) values as follows from the random variables. The code keeps a note of the start coordinates of the current path segment and updates it when moving to the next segment using these randomly selected lengths and directions:
\(l_i = -\frac{1}{\mu_T}ln(1-(1-e^{-\mu_T l_i^{\ max}})t_i)\)
\(\cos\theta_i = (k_i^2 + k_{i+1}^2 - Q_i^2)/2 k_i k_{i+1}\)
\(\phi_i = 2 \pi u_i\)
The final Monte Carlo integration that is actually performed by the code is as follows where N is the number of scenarios:
where the integration ranges over the rectangular \(Q \omega\) region are defined as follows:
\(\Delta\omega = \omega^{max}-\omega^{min}\)
\(\Delta Q_i = k_i + \frac{2m}{\hbar}\sqrt{\frac{\hbar^2 k_i^2}{2m} - \omega_{min}}\)
This is similar to the formulation described in the Mancinelli paper except there is no random variable to decide whether a particular scattering event is coherent or incoherent.
The integral \(I(k)\) is evaluated deterministically up front at a set of k values and interpolated as required.
The factor for the final track segment can also be normalised by setting NormalizeStructureFactors=true
which replaces \(\sigma_s\) with \(2k_n^2 \sigma_s(k_n)/I(k_n)\). This feature wasn’t in the original Discus implementation.
The results for different \(\omega\) values can be calculated by simulating tracks separately for each \(\omega\) value or the same tracks can be reused with the multiple weights for the final track segment being calculated to achieve the required range of overall energy transfers.
Discus used the latter approach which results in the results for different \(\omega\) being correlated. This choice is controlled using the SimulateEnergiesIndependently
parameter
Importance Sampling¶
The algorithm includes an option to use importance sampling to improve the results for elastic instrument when running with S(Q) profiles containing spikes. Without this option enabled, the contribution from rare, high values in the structure factor is only visible at a very high number of scenarios.
The importance sampling is achieved using a further change of variables as follows:
\(v_i = P(Q_i) = \frac{I(Q_i)}{I(2k)}\) where \(I(x) = \int_{0}^{x} Q S(Q) dQ\)
With this approach the Q value for each segment is chosen as follows based on a \(v_i\) value randomly selected between 0 and 1:
\(Q_i = P^{-1}(v_i)\)
\(\cos\theta_i\) is determined from \(Q_i\) as before. The change of variables gives the following integral for \(J_n\):
Finally, the equivalent Monte Carlo integration that the algorithm performs with importance sampling enabled is:
The importance sampling has also been implemented for inelastic instruments by flatting out the 2D \(S(Q, \omega)\) profile into a 1D array. A 1D coordinate is created which is the actual Q value added onto the maximum Q from the preceding \(\omega\) row: \(Q'(Q,\omega_i) = Q + Q_{max}(\omega_{i-1})\) With this approach there is no interpolation performed between different \(\omega\) values. It’s not clear whether the importance sampling is useful for inelastic calculations since the area where the multiple scattering correction tends to be largest relative to the signal is away from the peak in \(S(Q, \omega)\).
Support for sample environment¶
The calculation can include scattering from the sample environment (e.g. can) in the Monte Carlo simulation. The term “segment” has previously been used to refer to a straight neutron path between two scattering events. For the purpose of this description the term “link” will be used to refer to a subsection of a segment that lies within a single material.
The modified calculation is illustrated here with an example of a sample contained in a can where a track may contain three different links (can, then sample, then can). If the selected scatter point occurs somewhere in the third link, the quantity \(t_i\) is redefined as:
This can be generally expressed as follows where n is the number of sample environment components:
Based on this the length of the ith segment can be derived from a \(t_i\) that has been randomly selected between 0 and 1 as follows where again the expression is for the specific case of a track containing three different links:
…and more generally (although perhaps less helpfully in terms of explaining how the code works):
It can be seen that the formula (1) can be solved for \(l_i\) by calculating the quantity on the right hand side and then sequentially subtracting \(\mu_i l_i^{max}\) from it for increasing i while keeping the running total >=0. The value of \(i\) when you can’t subtract any more \(\mu_i l_i^{max}\) identifies the component containing the scatter. Dividing by \(\mu_i\) at this point gives you the length into that component that the track reaches.
The other modification to the calculation to support scattering in the sample environment is that a different structure factor \(S(Q,\omega)\) , \(I(k)\) and scattering cross section \(\sigma_s\) is required for each material. The component containing each scatter is derived from the \(l_i\) calculation and is used to look up the material.
Outputs¶
The algorithm outputs a workspace group containing the following workspaces:
Several workspaces called
Scatter_n
where n is the number of scattering events considered. Each workspace contains “per detector” weights as a function of momentum or energy transfer for a specific number of scattering events. The number of scattering events ranges between 1 and the number specified in the NumberOfScatterings parameterSeveral workspaces called
Scatter_n_Integrated
which are integrals of theScatter_n
workspaces across the x axis (Momentum for elastic and DeltaE for inelastic)A workspace called
Scatter_1_NoAbsorb
is also created for a scenario where neutrons are scattered once, absorption is assumed to be zero and re-scattering after the simulated scattering event is assumed to be zero. This is the quantity \(J_{1}^{*}\) described in the Discus manualA workspace called
Scatter_2_n_Summed
which is the sum of theScatter_n
workspaces for n > 1A workspace called
Scatter_1_n_Summed
which is the sum of theScatter_n
workspaces for n >= 1A workspace called
Ratio_Single_To_All
which is theScatter_1
workspace divided byScatter_1_n_Summed
The output can be applied to a workspace containing a real sample measurement in one of two ways:
subtraction method. The additional intensity contributed by multiple scattering to either a raw measurement or a vanadium corrected measurement can be calculated from the weights output from this algorithm. The additional intensity can then be subtracted to give an idealised “single scatter” intensity. For example, the additional intensity measured at a detector due to multiple scattering is given by \((\sum_{n=2}^{\infty} J_n) E(\lambda) I_0(\lambda) \Delta \Omega\) where \(E(\lambda)\) is the detector efficiency, \(I_0(\lambda)\) is the incident intensity and \(\Delta \Omega\) is the solid angle subtended by the detector. The factors \(E(\lambda) I_0(\lambda) \Delta \Omega\) can be obtained from a Vanadium run - although to take advantage of the “per detector” multiple scattering weights, the preparation of the Vanadium data will need to take place “per detector” instead of on focussed datasets
factor method. The correction can be applied by multiplying the real sample measurement by \(J_1/\sum_{n=1}^{\infty} J_n\). This approach avoids having to create a suitably normalised intensity from the weights and the method is also more tolerant of any normalisation inaccuracies in the S(Q) profile
The multiple scattering correction should be applied before applying an absorption correction.
The Discus manual describes a further method of applying an attenuation correction and a multiple scattering correction in one step using a variation of the factor method. To achieve this the real sample measurement should be multipled by \(J_1^{*}/(\sum_{n=1}^{\infty} J_n\)). Note that this differs from the approach taken in other Mantid absorption correction algorithms such as MonteCarloAbsorption because of the properties of \(J_{1}^{*}\). \(J_{1}^{*}\) corrects for attenuation due to absorption before and after the simulated scattering event (which is the same as MonteCarloAbsorption) but it only corrects for attenuation due to scattering after the simulated scattering event. For this reason it’s not clear this feature from Discus is useful but it has been left in for historical reasons.
The sample shape (and optionally the sample environment shape) can be specified by running the algorithms SetSample or LoadSampleShape on the input workspace prior to running this algorithm.
The algorithm can take a long time to run on instruments with a lot of spectra andor a lot of bins in each spectrum. The run time can be reduced by enabling the following interpolation features:
the multiple scattering correction can be calculated on a subset of the bins in the input workspace by specifying a non-default value for NumberOfSimulationPoints. The other points will be calculated by interpolation
the algorithm can be performed on a subset of the detectors by setting SparseInstrument=True
Both of these interpolation features are described further in the documentation for the MonteCarloAbsorption algorithm
Usage¶
Example - elastic calculation on single spike S(Q) and an isotropic S(Q) for comparison
The double scatter profile shows a similar shape to the analytic result calculated in [5]:
Example - inelastic calculation on direct geometry (matches calculation in DISCUS paper [1] figure 1)
This is the equivalent plot from the original Discus Fortran program:
References¶
Categories: AlgorithmIndex | CorrectionFunctions
Source¶
C++ header: DiscusMultipleScatteringCorrection.h
C++ source: DiscusMultipleScatteringCorrection.cpp