probe - Instrument probe

NeutronProbe

Neutron probe.

PolarizedNeutronProbe

Polarized neutron probe

PolarizedNeutronQProbe

alias of PolarizedQProbe

PolarizedQProbe

Probe

Defines the incident beam used to study the material.

ProbeSet

QProbe

A pure Q, R probe

Qmeasurement_union

Determine the unique Q, dQ across all datasets.

XrayProbe

X-Ray probe.

load4

Load in four column data Q, R, dR, dQ.

make_probe

Return a reflectometry measurement object of the given resolution.

measurement_union

Determine the unique (T, dT, L, dL) across all datasets.

spin_asymmetry

Compute spin asymmetry for R++, R--.

Experimental probe.

The experimental probe describes the incoming beam for the experiment. Scattering properties of the sample are dependent on the type and energy of the radiation.

See Data Representation for details.

class refl1d.probe.NeutronProbe(T=None, dT=0, L=None, dL=0, data=None, intensity=1, background=0, back_absorption=1, theta_offset=0, sample_broadening=0, back_reflectivity=False, name=None, filename=None, dQ=None, resolution='normal')[source]

Bases: Probe

Neutron probe.

By providing a scattering factor calculator for X-ray scattering, model components can be defined by mass density and chemical composition.

Aguide = 270.0
property Q
Q_c(substrate=None, surface=None)
property Ro
static alignment_uncertainty(w, I, d=0)

Compute alignment uncertainty.

Parameters:

wfloat | degrees

Rocking curve full width at half max.

Ifloat | counts

Rocking curve integrated intensity.

d = 0: float | degrees

Motor step size

Returns:

dthetafloat | degrees

uncertainty in alignment angle

apply_beam(calc_Q, calc_R, resolution=True, interpolation=0)

Apply factors such as beam intensity, background, backabsorption, resolution to the data.

resolution is True if the resolution function should be applied to the reflectivity.

interpolation is the number of Q points to show between the nominal Q points of the probe. Use this to draw a smooth theory line between the data points. The resolution dQ is interpolated between the resolution of the surrounding Q points.

If an amplitude signal is provided, \(r\) will be scaled by \(\surd I + i \surd B / |r|\), which when squared will equal \(I |r|^2 + B\). The resolution function will be applied directly to the amplitude. Unlike intensity and background, the resulting \(|G \ast r|^2 \ne G \ast |r|^2\) for convolution operator \(\ast\), but it should be close.

property calc_Q
critical_edge(substrate=None, surface=None, n=51, delta=0.25)

Oversample points near the critical edge.

The critical edge is defined by the difference in scattering potential for the substrate and surface materials, or the reverse if back_reflectivity is true.

n is the number of \(Q\) points to compute near the critical edge.

delta is the relative uncertainty in the material density, which defines the range of values which are calculated.

Note: critical_edge() will remove the extra Q calculation points introduced by oversample().

The \(n\) points \(Q_i\) are evenly distributed around the critical edge in \(Q_c \pm \delta Q_c\) by varying angle \(\theta\) for a fixed wavelength \(< \lambda >\), the average of all wavelengths in the probe.

Specifically:

\[\begin{split}Q_c^2 &= 16 \pi (\rho - \rho_\text{incident}) \\ Q_i &= Q_c - \delta_i Q_c (i - (n-1)/2) \qquad \text{for} \; i \in 0 \ldots n-1 \\ \lambda_i &= < \lambda > \\ \theta_i &= \sin^{-1}(Q_i \lambda_i / 4 \pi)\end{split}\]

If \(Q_c\) is imaginary, then \(-|Q_c|\) is used instead, so this routine can be used for reflectivity signals which scan from back reflectivity to front reflectivity. For completeness, the angle \(\theta = 0\) is added as well.

property dQ
fresnel(substrate=None, surface=None)

Returns a Fresnel reflectivity calculator given the surface and and substrate. The calculated reflectivity includes The Fresnel reflectivity for the probe reflecting from a block of material with the given substrate.

Returns F = R(probe.Q), where R is magnitude squared reflectivity.

label(prefix=None, gloss='', suffix='')
log10_to_linear()

Convert data from log to linear.

Older reflectometry reduction code stored reflectivity in log base 10 format. Call probe.log10_to_linear() after loading this data to convert it to linear for subsequent display and fitting.

oversample(n=20, seed=1)

Generate an over-sampling of Q to avoid aliasing effects.

Oversampling is needed for thick layers, in which the underlying reflectivity oscillates so rapidly in Q that a single measurement has contributions from multiple Kissig fringes.

Sampling will be done using a pseudo-random generator so that accidental structure in the function does not contribute to the aliasing. The generator will usually be initialized with a fixed seed so that the point selection will not change from run to run, but a seed of None will choose a different set of points each time oversample is called.

The value n is the number of points that should contribute to each Q value when computing the resolution. These will be distributed about the nominal measurement value, but varying in both angle and energy according to the resolution function. This will yield more points near the measurement and fewer farther away. The measurement point itself will not be used to avoid accidental bias from uniform Q steps. Depending on the problem, a value of n between 20 and 100 should lead to stable values for the convolved reflectivity.

Note: oversample() will remove the extra Q calculation points introduced by critical_edge().

parameters()
plot(view=None, **kwargs)

Plot theory against data.

Need substrate/surface for Fresnel-normalized reflectivity

plot_Q4(**kwargs)

Plot the Q**4 reflectivity associated with the probe.

Note that Q**4 reflectivity has the intensity and background applied so that hydrogenated samples display more cleanly. The formula to reproduce the graph is:

\[R' = R / ( (100*Q)^{-4} I + B) \Delta R' = \Delta R / ( (100*Q)^{-4} I + B )\]

where \(B\) is the background.

plot_fft(theory=None, suffix='', label=None, substrate=None, surface=None, **kwargs)

FFT analysis of reflectivity signal.

plot_fresnel(substrate=None, surface=None, **kwargs)

Plot the Fresnel-normalized reflectivity associated with the probe.

Note that the Fresnel reflectivity has the intensity and background applied before normalizing so that hydrogenated samples display more cleanly. The formula to reproduce the graph is:

\[R' = R / (F(Q) I + B) \Delta R' = \Delta R / (F(Q) I + B)\]

where \(I\) is the intensity and \(B\) is the background.

plot_linear(**kwargs)

Plot the data associated with probe.

plot_log(**kwargs)

Plot the data associated with probe.

plot_logfresnel(*args, **kw)

Plot the log Fresnel-normalized reflectivity associated with the probe.

plot_residuals(theory=None, suffix='', label=None, plot_shift=None, **kwargs)
plot_resolution(suffix='', label=None, **kwargs)
plot_shift = 0
polarized = False
radiation = 'neutron'
residuals_shift = 0
resolution_guard()

Make sure each measured \(Q\) point has at least 5 calculated \(Q\) points contributing to it in the range \([-3\Delta Q, 3\Delta Q]\).

Not Implemented

restore_data()

Restore the original data after resynth.

resynth_data()

Generate new data according to the model R’ ~ N(R, dR).

The resynthesis step is a precursor to refitting the data, as is required for certain types of monte carlo error analysis. The first time it is run it will save the original R into Ro. If you reset R in the probe you will also need to reset Ro so that it is used for subsequent resynth analysis.

save(filename, theory, substrate=None, surface=None)

Save the data and theory to a file.

scattering_factors(material, density)[source]

Returns the scattering factors associated with the material given the range of wavelengths/energies used in the probe.

show_resolution = True
simulate_data(theory, noise=2.0)

Set the data for the probe to R + eps with eps ~ normal(dR^2).

theory is (Q, R),

If the percent noise is provided, set dR to R*noise/100 before simulating. noise defaults to 2% if no dR is present.

Note that measured data estimates uncertainty from the number of counts. This means that points above the true value will have larger uncertainty than points below the true value. This bias is not captured in the simulated data.

subsample(dQ)

Select points at most every dQ.

Use this to speed up computation early in the fitting process.

This changes the data object, and is not reversible.

The current algorithm is not picking the “best” Q value, just the nearest, so if you have nearby Q points with different quality statistics (as happens in overlapped regions from spallation source measurements at different angles), then it may choose badly. Simple solutions based on the smallest relative error dR/R will be biased toward peaks, and smallest absolute error dR will be biased toward valleys.

to_dict()

Return a dictionary representation of the parameters

view = 'log'
write_data(filename, columns=('Q', 'R', 'dR'), header=None)

Save the data to a file.

header is a string with trailing n containing the file header. columns is a list of column names from Q, dQ, R, dR, L, dL, T, dT.

The default is to write Q, R, dR data.

class refl1d.probe.PolarizedNeutronProbe(xs=None, name=None, Aguide=270.0, H=0, oversampling=None)[source]

Bases: object

Polarized neutron probe

xs (4 x NeutronProbe) is a sequence pp, pm, mp and mm.

Aguide (degrees) is the angle of the applied field relative to the plane of the sample, with angle 270º in the plane of the sample.

H (tesla) is the magnitude of the applied field

apply_beam(Q, R, resolution=True, interpolation=0)[source]

Apply factors such as beam intensity, background, backabsorption, and footprint to the data.

property calc_Q
fresnel(*args, **kw)[source]

Returns a Fresnel reflectivity calculator given the surface and and substrate. The calculated reflectivity includes The Fresnel reflectivity for the probe reflecting from a block of material with the given substrate.

Returns F = R(probe.Q), where R is magnitude squared reflectivity.

property mm
property mp
oversample(n=6, seed=1)[source]
parameters()[source]
plot(view=None, **kwargs)[source]

Plot theory against data.

Need substrate/surface for Fresnel-normalized reflectivity

plot_Q4(**kwargs)[source]
plot_SA(theory=None, label=None, plot_shift=None, **kwargs)[source]
plot_fresnel(**kwargs)[source]
plot_linear(**kwargs)[source]
plot_log(**kwargs)[source]
plot_logfresnel(**kwargs)[source]
plot_residuals(**kwargs)[source]
plot_resolution(**kwargs)[source]
property pm
polarized = True
property pp
restore_data()[source]

Restore the original data after resynth.

resynth_data()[source]

Generate new data according to the model R’ ~ N(R, dR).

The resynthesis step is a precursor to refitting the data, as is required for certain types of monte carlo error analysis. The first time it is run it will save the original R into Ro. If you reset R in the probe you will also need to reset Ro so that it is used for subsequent resynth analysis.

save(filename, theory, substrate=None, surface=None)[source]

Save the data and theory to a file.

scattering_factors(material, density)[source]

Returns the scattering factors associated with the material given the range of wavelengths/energies used in the probe.

select_corresponding(theory)[source]

Select theory points corresponding to the measured data.

Since we have evaluated theory at every Q, it is safe to interpolate measured Q into theory, since it will land on a node, not in an interval.

shared_beam(intensity=1, background=0, back_absorption=1, theta_offset=0, sample_broadening=0)[source]

Share beam parameters across all four cross sections.

New parameters are created for intensity, background, theta_offset, sample_broadening and back_absorption and assigned to the all cross sections. These can be replaced with an explicit parameter in an individual cross section if that parameter is independent.

show_resolution = None
simulate_data(theory, noise=2.0)[source]

Set the data for the probe to R + eps with eps ~ normal(dR^2).

theory is (Q, R),

If the percent noise is provided, set dR to R*noise/100 before simulating. noise defaults to 2% if no dR is present.

Note that measured data estimates uncertainty from the number of counts. This means that points above the true value will have larger uncertainty than points below the true value. This bias is not captured in the simulated data.

substrate = None
surface = None
to_dict()[source]

Return a dictionary representation of the parameters

view = None
property xs
refl1d.probe.PolarizedNeutronQProbe

alias of PolarizedQProbe

class refl1d.probe.PolarizedQProbe(xs=None, name=None, Aguide=270.0, H=0)[source]

Bases: PolarizedNeutronProbe

apply_beam(Q, R, resolution=True, interpolation=0)

Apply factors such as beam intensity, background, backabsorption, and footprint to the data.

property calc_Q
fresnel(*args, **kw)

Returns a Fresnel reflectivity calculator given the surface and and substrate. The calculated reflectivity includes The Fresnel reflectivity for the probe reflecting from a block of material with the given substrate.

Returns F = R(probe.Q), where R is magnitude squared reflectivity.

property mm
property mp
oversample(n=6, seed=1)
parameters()
plot(view=None, **kwargs)

Plot theory against data.

Need substrate/surface for Fresnel-normalized reflectivity

plot_Q4(**kwargs)
plot_SA(theory=None, label=None, plot_shift=None, **kwargs)
plot_fresnel(**kwargs)
plot_linear(**kwargs)
plot_log(**kwargs)
plot_logfresnel(**kwargs)
plot_residuals(**kwargs)
plot_resolution(**kwargs)
property pm
polarized = True
property pp
restore_data()

Restore the original data after resynth.

resynth_data()

Generate new data according to the model R’ ~ N(R, dR).

The resynthesis step is a precursor to refitting the data, as is required for certain types of monte carlo error analysis. The first time it is run it will save the original R into Ro. If you reset R in the probe you will also need to reset Ro so that it is used for subsequent resynth analysis.

save(filename, theory, substrate=None, surface=None)

Save the data and theory to a file.

scattering_factors(material, density)

Returns the scattering factors associated with the material given the range of wavelengths/energies used in the probe.

select_corresponding(theory)

Select theory points corresponding to the measured data.

Since we have evaluated theory at every Q, it is safe to interpolate measured Q into theory, since it will land on a node, not in an interval.

shared_beam(intensity=1, background=0, back_absorption=1, theta_offset=0, sample_broadening=0)

Share beam parameters across all four cross sections.

New parameters are created for intensity, background, theta_offset, sample_broadening and back_absorption and assigned to the all cross sections. These can be replaced with an explicit parameter in an individual cross section if that parameter is independent.

show_resolution = None
simulate_data(theory, noise=2.0)

Set the data for the probe to R + eps with eps ~ normal(dR^2).

theory is (Q, R),

If the percent noise is provided, set dR to R*noise/100 before simulating. noise defaults to 2% if no dR is present.

Note that measured data estimates uncertainty from the number of counts. This means that points above the true value will have larger uncertainty than points below the true value. This bias is not captured in the simulated data.

substrate = None
surface = None
to_dict()

Return a dictionary representation of the parameters

view = None
property xs
class refl1d.probe.Probe(T=None, dT=0, L=None, dL=0, data=None, intensity=1, background=0, back_absorption=1, theta_offset=0, sample_broadening=0, back_reflectivity=False, name=None, filename=None, dQ=None, resolution='normal')[source]

Bases: object

Defines the incident beam used to study the material.

For calculation purposes, probe needs to return the values \(Q_\text{calc}\) at which the model is evaluated. This is normally going to be the measured points only, but for some systems, such as those with very thick layers, oversampling is needed to avoid aliasing effects.

A measurement point consists of incident angle, angular resolution, incident wavelength, FWHM wavelength resolution, reflectivity and uncertainty in reflectivity.

A probe is a set of points, defined by vectors for point attribute. For convenience, the attribute can be initialized with a scalar if it is constant throughout the measurement, but will be set to a vector in the probe. The attributes are initialized as follows:

Tfloat or [float] | degrees

Incident angle

dTfloat or [float] | degrees

FWHM angular divergence

Lfloat or [float] | Å

Incident wavelength

dLfloat or [float] | Å

FWHM wavelength dispersion

data([float], [float])

R, dR reflectivity measurement and uncertainty

dQ[float] or None | Å-1

1-$sigma$ Q resolution when it cannot be computed directly from angular divergence and wavelength dispersion.

resolution‘normal’ or ‘uniform’

Distribution function for Q resolution.

Measurement properties:

intensityfloat or Parameter

Beam intensity

backgroundfloat or Parameter

Constant background

back_absorptionfloat or Parameter

Absorption through the substrate relative to beam intensity. A value of 1.0 means complete transmission; a value of 0.0 means complete absorption.

theta_offsetfloat or Parameter

Offset of the sample from perfect alignment

sample_broadeningfloat or Parameter

Additional FWHM angular divergence from sample curvature. Scale 1-\(\sigma\) rms by \(2 \surd(2 \ln 2) \approx 2.35\) to convert to FWHM.

back_reflectivityTrue or False

True if the beam enters through the substrate

Measurement properties are fittable parameters. theta_offset in particular should be set using probe.theta_offset.dev(dT), with dT equal to the FWHM uncertainty in the peak position for the rocking curve, as measured in radians. Changes to theta_offset will then be penalized in the cost function for the fit as if it were another measurement. Use alignment_uncertainty() to compute dT from the shape of the rocking curve.

Sample broadening adjusts the existing Q resolution rather than recalculating it. This allows it the resolution to describe more complicated effects than a simple gaussian distribution of wavelength and angle will allow. The calculation uses the mean wavelength, angle and angular divergence. See resolution.dQ_broadening() for details.

intensity and back_absorption are generally not needed — scaling the reflected signal by an appropriate intensity measurement will correct for both of these during reduction. background may be needed, particularly for samples with significant hydrogen content due to its large isotropic incoherent scattering cross section.

View properties:

viewstring

One of ‘fresnel’, ‘logfresnel’, ‘log’, ‘linear’, ‘q4’, ‘residuals’

show_resolutionbool

True if resolution bars should be plotted with each point.

plot_shiftfloat

The number of pixels to shift each new dataset so datasets can be seen separately

residuals_shift :

The number of pixels to shift each new set of residuals so the residuals plots can be seen separately.

Normally view is set directly in the class rather than the instance since it is not specific to the view. Fresnel and Q4 views are corrected for background and intensity; log and linear views show the uncorrected data. The Fresnel reflectivity calculation has resolution applied.

Aguide = 270.0
property Q
Q_c(substrate=None, surface=None)[source]
property Ro
static alignment_uncertainty(w, I, d=0)[source]

Compute alignment uncertainty.

Parameters:

wfloat | degrees

Rocking curve full width at half max.

Ifloat | counts

Rocking curve integrated intensity.

d = 0: float | degrees

Motor step size

Returns:

dthetafloat | degrees

uncertainty in alignment angle

apply_beam(calc_Q, calc_R, resolution=True, interpolation=0)[source]

Apply factors such as beam intensity, background, backabsorption, resolution to the data.

resolution is True if the resolution function should be applied to the reflectivity.

interpolation is the number of Q points to show between the nominal Q points of the probe. Use this to draw a smooth theory line between the data points. The resolution dQ is interpolated between the resolution of the surrounding Q points.

If an amplitude signal is provided, \(r\) will be scaled by \(\surd I + i \surd B / |r|\), which when squared will equal \(I |r|^2 + B\). The resolution function will be applied directly to the amplitude. Unlike intensity and background, the resulting \(|G \ast r|^2 \ne G \ast |r|^2\) for convolution operator \(\ast\), but it should be close.

property calc_Q
critical_edge(substrate=None, surface=None, n=51, delta=0.25)[source]

Oversample points near the critical edge.

The critical edge is defined by the difference in scattering potential for the substrate and surface materials, or the reverse if back_reflectivity is true.

n is the number of \(Q\) points to compute near the critical edge.

delta is the relative uncertainty in the material density, which defines the range of values which are calculated.

Note: critical_edge() will remove the extra Q calculation points introduced by oversample().

The \(n\) points \(Q_i\) are evenly distributed around the critical edge in \(Q_c \pm \delta Q_c\) by varying angle \(\theta\) for a fixed wavelength \(< \lambda >\), the average of all wavelengths in the probe.

Specifically:

\[\begin{split}Q_c^2 &= 16 \pi (\rho - \rho_\text{incident}) \\ Q_i &= Q_c - \delta_i Q_c (i - (n-1)/2) \qquad \text{for} \; i \in 0 \ldots n-1 \\ \lambda_i &= < \lambda > \\ \theta_i &= \sin^{-1}(Q_i \lambda_i / 4 \pi)\end{split}\]

If \(Q_c\) is imaginary, then \(-|Q_c|\) is used instead, so this routine can be used for reflectivity signals which scan from back reflectivity to front reflectivity. For completeness, the angle \(\theta = 0\) is added as well.

property dQ
fresnel(substrate=None, surface=None)[source]

Returns a Fresnel reflectivity calculator given the surface and and substrate. The calculated reflectivity includes The Fresnel reflectivity for the probe reflecting from a block of material with the given substrate.

Returns F = R(probe.Q), where R is magnitude squared reflectivity.

label(prefix=None, gloss='', suffix='')[source]
log10_to_linear()[source]

Convert data from log to linear.

Older reflectometry reduction code stored reflectivity in log base 10 format. Call probe.log10_to_linear() after loading this data to convert it to linear for subsequent display and fitting.

oversample(n=20, seed=1)[source]

Generate an over-sampling of Q to avoid aliasing effects.

Oversampling is needed for thick layers, in which the underlying reflectivity oscillates so rapidly in Q that a single measurement has contributions from multiple Kissig fringes.

Sampling will be done using a pseudo-random generator so that accidental structure in the function does not contribute to the aliasing. The generator will usually be initialized with a fixed seed so that the point selection will not change from run to run, but a seed of None will choose a different set of points each time oversample is called.

The value n is the number of points that should contribute to each Q value when computing the resolution. These will be distributed about the nominal measurement value, but varying in both angle and energy according to the resolution function. This will yield more points near the measurement and fewer farther away. The measurement point itself will not be used to avoid accidental bias from uniform Q steps. Depending on the problem, a value of n between 20 and 100 should lead to stable values for the convolved reflectivity.

Note: oversample() will remove the extra Q calculation points introduced by critical_edge().

parameters()[source]
plot(view=None, **kwargs)[source]

Plot theory against data.

Need substrate/surface for Fresnel-normalized reflectivity

plot_Q4(**kwargs)[source]

Plot the Q**4 reflectivity associated with the probe.

Note that Q**4 reflectivity has the intensity and background applied so that hydrogenated samples display more cleanly. The formula to reproduce the graph is:

\[R' = R / ( (100*Q)^{-4} I + B) \Delta R' = \Delta R / ( (100*Q)^{-4} I + B )\]

where \(B\) is the background.

plot_fft(theory=None, suffix='', label=None, substrate=None, surface=None, **kwargs)[source]

FFT analysis of reflectivity signal.

plot_fresnel(substrate=None, surface=None, **kwargs)[source]

Plot the Fresnel-normalized reflectivity associated with the probe.

Note that the Fresnel reflectivity has the intensity and background applied before normalizing so that hydrogenated samples display more cleanly. The formula to reproduce the graph is:

\[R' = R / (F(Q) I + B) \Delta R' = \Delta R / (F(Q) I + B)\]

where \(I\) is the intensity and \(B\) is the background.

plot_linear(**kwargs)[source]

Plot the data associated with probe.

plot_log(**kwargs)[source]

Plot the data associated with probe.

plot_logfresnel(*args, **kw)[source]

Plot the log Fresnel-normalized reflectivity associated with the probe.

plot_residuals(theory=None, suffix='', label=None, plot_shift=None, **kwargs)[source]
plot_resolution(suffix='', label=None, **kwargs)[source]
plot_shift = 0
polarized = False
residuals_shift = 0
resolution_guard()[source]

Make sure each measured \(Q\) point has at least 5 calculated \(Q\) points contributing to it in the range \([-3\Delta Q, 3\Delta Q]\).

Not Implemented

restore_data()[source]

Restore the original data after resynth.

resynth_data()[source]

Generate new data according to the model R’ ~ N(R, dR).

The resynthesis step is a precursor to refitting the data, as is required for certain types of monte carlo error analysis. The first time it is run it will save the original R into Ro. If you reset R in the probe you will also need to reset Ro so that it is used for subsequent resynth analysis.

save(filename, theory, substrate=None, surface=None)[source]

Save the data and theory to a file.

scattering_factors(material, density)[source]

Returns the scattering factors associated with the material given the range of wavelengths/energies used in the probe.

show_resolution = True
simulate_data(theory, noise=2.0)[source]

Set the data for the probe to R + eps with eps ~ normal(dR^2).

theory is (Q, R),

If the percent noise is provided, set dR to R*noise/100 before simulating. noise defaults to 2% if no dR is present.

Note that measured data estimates uncertainty from the number of counts. This means that points above the true value will have larger uncertainty than points below the true value. This bias is not captured in the simulated data.

subsample(dQ)[source]

Select points at most every dQ.

Use this to speed up computation early in the fitting process.

This changes the data object, and is not reversible.

The current algorithm is not picking the “best” Q value, just the nearest, so if you have nearby Q points with different quality statistics (as happens in overlapped regions from spallation source measurements at different angles), then it may choose badly. Simple solutions based on the smallest relative error dR/R will be biased toward peaks, and smallest absolute error dR will be biased toward valleys.

to_dict()[source]

Return a dictionary representation of the parameters

view = 'log'
write_data(filename, columns=('Q', 'R', 'dR'), header=None)[source]

Save the data to a file.

header is a string with trailing n containing the file header. columns is a list of column names from Q, dQ, R, dR, L, dL, T, dT.

The default is to write Q, R, dR data.

class refl1d.probe.ProbeSet(probes, name=None)[source]

Bases: Probe

Aguide = 270.0
property Q
Q_c(substrate=None, surface=None)
property Ro
static alignment_uncertainty(w, I, d=0)

Compute alignment uncertainty.

Parameters:

wfloat | degrees

Rocking curve full width at half max.

Ifloat | counts

Rocking curve integrated intensity.

d = 0: float | degrees

Motor step size

Returns:

dthetafloat | degrees

uncertainty in alignment angle

apply_beam(calc_Q, calc_R, interpolation=0, **kw)[source]

Apply factors such as beam intensity, background, backabsorption, resolution to the data.

resolution is True if the resolution function should be applied to the reflectivity.

interpolation is the number of Q points to show between the nominal Q points of the probe. Use this to draw a smooth theory line between the data points. The resolution dQ is interpolated between the resolution of the surrounding Q points.

If an amplitude signal is provided, \(r\) will be scaled by \(\surd I + i \surd B / |r|\), which when squared will equal \(I |r|^2 + B\). The resolution function will be applied directly to the amplitude. Unlike intensity and background, the resulting \(|G \ast r|^2 \ne G \ast |r|^2\) for convolution operator \(\ast\), but it should be close.

property calc_Q
critical_edge(substrate=None, surface=None, n=51, delta=0.25)

Oversample points near the critical edge.

The critical edge is defined by the difference in scattering potential for the substrate and surface materials, or the reverse if back_reflectivity is true.

n is the number of \(Q\) points to compute near the critical edge.

delta is the relative uncertainty in the material density, which defines the range of values which are calculated.

Note: critical_edge() will remove the extra Q calculation points introduced by oversample().

The \(n\) points \(Q_i\) are evenly distributed around the critical edge in \(Q_c \pm \delta Q_c\) by varying angle \(\theta\) for a fixed wavelength \(< \lambda >\), the average of all wavelengths in the probe.

Specifically:

\[\begin{split}Q_c^2 &= 16 \pi (\rho - \rho_\text{incident}) \\ Q_i &= Q_c - \delta_i Q_c (i - (n-1)/2) \qquad \text{for} \; i \in 0 \ldots n-1 \\ \lambda_i &= < \lambda > \\ \theta_i &= \sin^{-1}(Q_i \lambda_i / 4 \pi)\end{split}\]

If \(Q_c\) is imaginary, then \(-|Q_c|\) is used instead, so this routine can be used for reflectivity signals which scan from back reflectivity to front reflectivity. For completeness, the angle \(\theta = 0\) is added as well.

property dQ
fresnel(*args, **kw)[source]

Returns a Fresnel reflectivity calculator given the surface and and substrate. The calculated reflectivity includes The Fresnel reflectivity for the probe reflecting from a block of material with the given substrate.

Returns F = R(probe.Q), where R is magnitude squared reflectivity.

label(prefix=None, gloss='', suffix='')
log10_to_linear()

Convert data from log to linear.

Older reflectometry reduction code stored reflectivity in log base 10 format. Call probe.log10_to_linear() after loading this data to convert it to linear for subsequent display and fitting.

oversample(**kw)[source]

Generate an over-sampling of Q to avoid aliasing effects.

Oversampling is needed for thick layers, in which the underlying reflectivity oscillates so rapidly in Q that a single measurement has contributions from multiple Kissig fringes.

Sampling will be done using a pseudo-random generator so that accidental structure in the function does not contribute to the aliasing. The generator will usually be initialized with a fixed seed so that the point selection will not change from run to run, but a seed of None will choose a different set of points each time oversample is called.

The value n is the number of points that should contribute to each Q value when computing the resolution. These will be distributed about the nominal measurement value, but varying in both angle and energy according to the resolution function. This will yield more points near the measurement and fewer farther away. The measurement point itself will not be used to avoid accidental bias from uniform Q steps. Depending on the problem, a value of n between 20 and 100 should lead to stable values for the convolved reflectivity.

Note: oversample() will remove the extra Q calculation points introduced by critical_edge().

parameters()[source]
parts(theory)[source]
plot(theory=None, **kw)[source]

Plot theory against data.

Need substrate/surface for Fresnel-normalized reflectivity

plot_Q4(theory=None, **kw)[source]

Plot the Q**4 reflectivity associated with the probe.

Note that Q**4 reflectivity has the intensity and background applied so that hydrogenated samples display more cleanly. The formula to reproduce the graph is:

\[R' = R / ( (100*Q)^{-4} I + B) \Delta R' = \Delta R / ( (100*Q)^{-4} I + B )\]

where \(B\) is the background.

plot_fft(theory=None, suffix='', label=None, substrate=None, surface=None, **kwargs)

FFT analysis of reflectivity signal.

plot_fresnel(theory=None, **kw)[source]

Plot the Fresnel-normalized reflectivity associated with the probe.

Note that the Fresnel reflectivity has the intensity and background applied before normalizing so that hydrogenated samples display more cleanly. The formula to reproduce the graph is:

\[R' = R / (F(Q) I + B) \Delta R' = \Delta R / (F(Q) I + B)\]

where \(I\) is the intensity and \(B\) is the background.

plot_linear(theory=None, **kw)[source]

Plot the data associated with probe.

plot_log(theory=None, **kw)[source]

Plot the data associated with probe.

plot_logfresnel(theory=None, **kw)[source]

Plot the log Fresnel-normalized reflectivity associated with the probe.

plot_residuals(theory=None, **kw)[source]
plot_resolution(**kw)[source]
plot_shift = 0
polarized = False
residuals_shift = 0
resolution_guard()

Make sure each measured \(Q\) point has at least 5 calculated \(Q\) points contributing to it in the range \([-3\Delta Q, 3\Delta Q]\).

Not Implemented

restore_data()[source]

Restore the original data after resynth.

resynth_data()[source]

Generate new data according to the model R’ ~ N(R, dR).

The resynthesis step is a precursor to refitting the data, as is required for certain types of monte carlo error analysis. The first time it is run it will save the original R into Ro. If you reset R in the probe you will also need to reset Ro so that it is used for subsequent resynth analysis.

save(filename, theory, substrate=None, surface=None)[source]

Save the data and theory to a file.

scattering_factors(material, density)[source]

Returns the scattering factors associated with the material given the range of wavelengths/energies used in the probe.

shared_beam(intensity=1, background=0, back_absorption=1, theta_offset=0, sample_broadening=0)[source]

Share beam parameters across all segments.

New parameters are created for intensity, background, theta_offset, sample_broadening and back_absorption and assigned to the all segments. These can be replaced with an explicit parameter in an individual segment if that parameter is independent.

show_resolution = True
simulate_data(theory, noise=2.0)[source]

Set the data for the probe to R + eps with eps ~ normal(dR^2).

theory is (Q, R),

If the percent noise is provided, set dR to R*noise/100 before simulating. noise defaults to 2% if no dR is present.

Note that measured data estimates uncertainty from the number of counts. This means that points above the true value will have larger uncertainty than points below the true value. This bias is not captured in the simulated data.

stitch(same_Q=0.001, same_dQ=0.001)[source]

Stitch together multiple datasets into a single dataset.

Points within tol of each other and with the same resolution are combined by interpolating them to a common \(Q\) value then averaged using Gaussian error propagation.

Returns:

probe | Probe Combined data set.

Algorithm:

To interpolate a set of points to a common value, first find the common \(Q\) value:

\[\hat Q = \sum{Q_k} / n\]

Then for each dataset \(k\), find the interval \([i, i+1]\) containing the value \(Q\), and use it to compute interpolated value for \(R\):

\[\begin{split}w &= (\hat Q - Q_i)/(Q_{i+1} - Q_i) \\ \hat R &= w R_{i+1} + (1-w) R_{i+1} \\ \hat \sigma_{R} &= \sqrt{ w^2 \sigma^2_{R_i} + (1-w)^2 \sigma^2_{R_{i+1}} } / n\end{split}\]

Average the resulting \(R\) using Gaussian error propagation:

\[\begin{split}\hat R &= \sum{\hat R_k}/n \\ \hat \sigma_R &= \sqrt{\sum \hat \sigma_{R_k}^2}/n\end{split}\]
subsample(dQ)

Select points at most every dQ.

Use this to speed up computation early in the fitting process.

This changes the data object, and is not reversible.

The current algorithm is not picking the “best” Q value, just the nearest, so if you have nearby Q points with different quality statistics (as happens in overlapped regions from spallation source measurements at different angles), then it may choose badly. Simple solutions based on the smallest relative error dR/R will be biased toward peaks, and smallest absolute error dR will be biased toward valleys.

to_dict()[source]

Return a dictionary representation of the parameters

property unique_L
view = 'log'
write_data(filename, columns=('Q', 'R', 'dR'), header=None)

Save the data to a file.

header is a string with trailing n containing the file header. columns is a list of column names from Q, dQ, R, dR, L, dL, T, dT.

The default is to write Q, R, dR data.

class refl1d.probe.QProbe(Q, dQ, data=None, name=None, filename=None, intensity=1, background=0, back_absorption=1, back_reflectivity=False, resolution='normal')[source]

Bases: Probe

A pure Q, R probe

This probe with no possibility of tricks such as looking up the scattering length density based on wavelength, or adjusting for alignment errors.

Aguide = 270.0
property Q
Q_c(substrate=None, surface=None)
property Ro
static alignment_uncertainty(w, I, d=0)

Compute alignment uncertainty.

Parameters:

wfloat | degrees

Rocking curve full width at half max.

Ifloat | counts

Rocking curve integrated intensity.

d = 0: float | degrees

Motor step size

Returns:

dthetafloat | degrees

uncertainty in alignment angle

apply_beam(calc_Q, calc_R, resolution=True, interpolation=0)

Apply factors such as beam intensity, background, backabsorption, resolution to the data.

resolution is True if the resolution function should be applied to the reflectivity.

interpolation is the number of Q points to show between the nominal Q points of the probe. Use this to draw a smooth theory line between the data points. The resolution dQ is interpolated between the resolution of the surrounding Q points.

If an amplitude signal is provided, \(r\) will be scaled by \(\surd I + i \surd B / |r|\), which when squared will equal \(I |r|^2 + B\). The resolution function will be applied directly to the amplitude. Unlike intensity and background, the resulting \(|G \ast r|^2 \ne G \ast |r|^2\) for convolution operator \(\ast\), but it should be close.

property calc_Q
critical_edge(substrate=None, surface=None, n=51, delta=0.25)[source]

Oversample points near the critical edge.

The critical edge is defined by the difference in scattering potential for the substrate and surface materials, or the reverse if back_reflectivity is true.

n is the number of \(Q\) points to compute near the critical edge.

delta is the relative uncertainty in the material density, which defines the range of values which are calculated.

Note: critical_edge() will remove the extra Q calculation points introduced by oversample().

The \(n\) points \(Q_i\) are evenly distributed around the critical edge in \(Q_c \pm \delta Q_c\) by varying angle \(\theta\) for a fixed wavelength \(< \lambda >\), the average of all wavelengths in the probe.

Specifically:

\[\begin{split}Q_c^2 &= 16 \pi (\rho - \rho_\text{incident}) \\ Q_i &= Q_c - \delta_i Q_c (i - (n-1)/2) \qquad \text{for} \; i \in 0 \ldots n-1 \\ \lambda_i &= < \lambda > \\ \theta_i &= \sin^{-1}(Q_i \lambda_i / 4 \pi)\end{split}\]

If \(Q_c\) is imaginary, then \(-|Q_c|\) is used instead, so this routine can be used for reflectivity signals which scan from back reflectivity to front reflectivity. For completeness, the angle \(\theta = 0\) is added as well.

property dQ
fresnel(substrate=None, surface=None)

Returns a Fresnel reflectivity calculator given the surface and and substrate. The calculated reflectivity includes The Fresnel reflectivity for the probe reflecting from a block of material with the given substrate.

Returns F = R(probe.Q), where R is magnitude squared reflectivity.

label(prefix=None, gloss='', suffix='')
log10_to_linear()

Convert data from log to linear.

Older reflectometry reduction code stored reflectivity in log base 10 format. Call probe.log10_to_linear() after loading this data to convert it to linear for subsequent display and fitting.

oversample(n=20, seed=1)[source]

Generate an over-sampling of Q to avoid aliasing effects.

Oversampling is needed for thick layers, in which the underlying reflectivity oscillates so rapidly in Q that a single measurement has contributions from multiple Kissig fringes.

Sampling will be done using a pseudo-random generator so that accidental structure in the function does not contribute to the aliasing. The generator will usually be initialized with a fixed seed so that the point selection will not change from run to run, but a seed of None will choose a different set of points each time oversample is called.

The value n is the number of points that should contribute to each Q value when computing the resolution. These will be distributed about the nominal measurement value, but varying in both angle and energy according to the resolution function. This will yield more points near the measurement and fewer farther away. The measurement point itself will not be used to avoid accidental bias from uniform Q steps. Depending on the problem, a value of n between 20 and 100 should lead to stable values for the convolved reflectivity.

Note: oversample() will remove the extra Q calculation points introduced by critical_edge().

parameters()
plot(view=None, **kwargs)

Plot theory against data.

Need substrate/surface for Fresnel-normalized reflectivity

plot_Q4(**kwargs)

Plot the Q**4 reflectivity associated with the probe.

Note that Q**4 reflectivity has the intensity and background applied so that hydrogenated samples display more cleanly. The formula to reproduce the graph is:

\[R' = R / ( (100*Q)^{-4} I + B) \Delta R' = \Delta R / ( (100*Q)^{-4} I + B )\]

where \(B\) is the background.

plot_fft(theory=None, suffix='', label=None, substrate=None, surface=None, **kwargs)

FFT analysis of reflectivity signal.

plot_fresnel(substrate=None, surface=None, **kwargs)

Plot the Fresnel-normalized reflectivity associated with the probe.

Note that the Fresnel reflectivity has the intensity and background applied before normalizing so that hydrogenated samples display more cleanly. The formula to reproduce the graph is:

\[R' = R / (F(Q) I + B) \Delta R' = \Delta R / (F(Q) I + B)\]

where \(I\) is the intensity and \(B\) is the background.

plot_linear(**kwargs)

Plot the data associated with probe.

plot_log(**kwargs)

Plot the data associated with probe.

plot_logfresnel(*args, **kw)

Plot the log Fresnel-normalized reflectivity associated with the probe.

plot_residuals(theory=None, suffix='', label=None, plot_shift=None, **kwargs)
plot_resolution(suffix='', label=None, **kwargs)
plot_shift = 0
polarized = False
residuals_shift = 0
resolution_guard()

Make sure each measured \(Q\) point has at least 5 calculated \(Q\) points contributing to it in the range \([-3\Delta Q, 3\Delta Q]\).

Not Implemented

restore_data()

Restore the original data after resynth.

resynth_data()

Generate new data according to the model R’ ~ N(R, dR).

The resynthesis step is a precursor to refitting the data, as is required for certain types of monte carlo error analysis. The first time it is run it will save the original R into Ro. If you reset R in the probe you will also need to reset Ro so that it is used for subsequent resynth analysis.

save(filename, theory, substrate=None, surface=None)

Save the data and theory to a file.

scattering_factors(material, density)[source]

Returns the scattering factors associated with the material given the range of wavelengths/energies used in the probe.

show_resolution = True
simulate_data(theory, noise=2.0)

Set the data for the probe to R + eps with eps ~ normal(dR^2).

theory is (Q, R),

If the percent noise is provided, set dR to R*noise/100 before simulating. noise defaults to 2% if no dR is present.

Note that measured data estimates uncertainty from the number of counts. This means that points above the true value will have larger uncertainty than points below the true value. This bias is not captured in the simulated data.

subsample(dQ)

Select points at most every dQ.

Use this to speed up computation early in the fitting process.

This changes the data object, and is not reversible.

The current algorithm is not picking the “best” Q value, just the nearest, so if you have nearby Q points with different quality statistics (as happens in overlapped regions from spallation source measurements at different angles), then it may choose badly. Simple solutions based on the smallest relative error dR/R will be biased toward peaks, and smallest absolute error dR will be biased toward valleys.

to_dict()

Return a dictionary representation of the parameters

view = 'log'
write_data(filename, columns=('Q', 'R', 'dR'), header=None)

Save the data to a file.

header is a string with trailing n containing the file header. columns is a list of column names from Q, dQ, R, dR, L, dL, T, dT.

The default is to write Q, R, dR data.

refl1d.probe.Qmeasurement_union(xs)[source]

Determine the unique Q, dQ across all datasets.

class refl1d.probe.XrayProbe(T=None, dT=0, L=None, dL=0, data=None, intensity=1, background=0, back_absorption=1, theta_offset=0, sample_broadening=0, back_reflectivity=False, name=None, filename=None, dQ=None, resolution='normal')[source]

Bases: Probe

X-Ray probe.

By providing a scattering factor calculator for X-ray scattering, model components can be defined by mass density and chemical composition.

Aguide = 270.0
property Q
Q_c(substrate=None, surface=None)
property Ro
static alignment_uncertainty(w, I, d=0)

Compute alignment uncertainty.

Parameters:

wfloat | degrees

Rocking curve full width at half max.

Ifloat | counts

Rocking curve integrated intensity.

d = 0: float | degrees

Motor step size

Returns:

dthetafloat | degrees

uncertainty in alignment angle

apply_beam(calc_Q, calc_R, resolution=True, interpolation=0)

Apply factors such as beam intensity, background, backabsorption, resolution to the data.

resolution is True if the resolution function should be applied to the reflectivity.

interpolation is the number of Q points to show between the nominal Q points of the probe. Use this to draw a smooth theory line between the data points. The resolution dQ is interpolated between the resolution of the surrounding Q points.

If an amplitude signal is provided, \(r\) will be scaled by \(\surd I + i \surd B / |r|\), which when squared will equal \(I |r|^2 + B\). The resolution function will be applied directly to the amplitude. Unlike intensity and background, the resulting \(|G \ast r|^2 \ne G \ast |r|^2\) for convolution operator \(\ast\), but it should be close.

property calc_Q
critical_edge(substrate=None, surface=None, n=51, delta=0.25)

Oversample points near the critical edge.

The critical edge is defined by the difference in scattering potential for the substrate and surface materials, or the reverse if back_reflectivity is true.

n is the number of \(Q\) points to compute near the critical edge.

delta is the relative uncertainty in the material density, which defines the range of values which are calculated.

Note: critical_edge() will remove the extra Q calculation points introduced by oversample().

The \(n\) points \(Q_i\) are evenly distributed around the critical edge in \(Q_c \pm \delta Q_c\) by varying angle \(\theta\) for a fixed wavelength \(< \lambda >\), the average of all wavelengths in the probe.

Specifically:

\[\begin{split}Q_c^2 &= 16 \pi (\rho - \rho_\text{incident}) \\ Q_i &= Q_c - \delta_i Q_c (i - (n-1)/2) \qquad \text{for} \; i \in 0 \ldots n-1 \\ \lambda_i &= < \lambda > \\ \theta_i &= \sin^{-1}(Q_i \lambda_i / 4 \pi)\end{split}\]

If \(Q_c\) is imaginary, then \(-|Q_c|\) is used instead, so this routine can be used for reflectivity signals which scan from back reflectivity to front reflectivity. For completeness, the angle \(\theta = 0\) is added as well.

property dQ
fresnel(substrate=None, surface=None)

Returns a Fresnel reflectivity calculator given the surface and and substrate. The calculated reflectivity includes The Fresnel reflectivity for the probe reflecting from a block of material with the given substrate.

Returns F = R(probe.Q), where R is magnitude squared reflectivity.

label(prefix=None, gloss='', suffix='')
log10_to_linear()

Convert data from log to linear.

Older reflectometry reduction code stored reflectivity in log base 10 format. Call probe.log10_to_linear() after loading this data to convert it to linear for subsequent display and fitting.

oversample(n=20, seed=1)

Generate an over-sampling of Q to avoid aliasing effects.

Oversampling is needed for thick layers, in which the underlying reflectivity oscillates so rapidly in Q that a single measurement has contributions from multiple Kissig fringes.

Sampling will be done using a pseudo-random generator so that accidental structure in the function does not contribute to the aliasing. The generator will usually be initialized with a fixed seed so that the point selection will not change from run to run, but a seed of None will choose a different set of points each time oversample is called.

The value n is the number of points that should contribute to each Q value when computing the resolution. These will be distributed about the nominal measurement value, but varying in both angle and energy according to the resolution function. This will yield more points near the measurement and fewer farther away. The measurement point itself will not be used to avoid accidental bias from uniform Q steps. Depending on the problem, a value of n between 20 and 100 should lead to stable values for the convolved reflectivity.

Note: oversample() will remove the extra Q calculation points introduced by critical_edge().

parameters()
plot(view=None, **kwargs)

Plot theory against data.

Need substrate/surface for Fresnel-normalized reflectivity

plot_Q4(**kwargs)

Plot the Q**4 reflectivity associated with the probe.

Note that Q**4 reflectivity has the intensity and background applied so that hydrogenated samples display more cleanly. The formula to reproduce the graph is:

\[R' = R / ( (100*Q)^{-4} I + B) \Delta R' = \Delta R / ( (100*Q)^{-4} I + B )\]

where \(B\) is the background.

plot_fft(theory=None, suffix='', label=None, substrate=None, surface=None, **kwargs)

FFT analysis of reflectivity signal.

plot_fresnel(substrate=None, surface=None, **kwargs)

Plot the Fresnel-normalized reflectivity associated with the probe.

Note that the Fresnel reflectivity has the intensity and background applied before normalizing so that hydrogenated samples display more cleanly. The formula to reproduce the graph is:

\[R' = R / (F(Q) I + B) \Delta R' = \Delta R / (F(Q) I + B)\]

where \(I\) is the intensity and \(B\) is the background.

plot_linear(**kwargs)

Plot the data associated with probe.

plot_log(**kwargs)

Plot the data associated with probe.

plot_logfresnel(*args, **kw)

Plot the log Fresnel-normalized reflectivity associated with the probe.

plot_residuals(theory=None, suffix='', label=None, plot_shift=None, **kwargs)
plot_resolution(suffix='', label=None, **kwargs)
plot_shift = 0
polarized = False
radiation = 'xray'
residuals_shift = 0
resolution_guard()

Make sure each measured \(Q\) point has at least 5 calculated \(Q\) points contributing to it in the range \([-3\Delta Q, 3\Delta Q]\).

Not Implemented

restore_data()

Restore the original data after resynth.

resynth_data()

Generate new data according to the model R’ ~ N(R, dR).

The resynthesis step is a precursor to refitting the data, as is required for certain types of monte carlo error analysis. The first time it is run it will save the original R into Ro. If you reset R in the probe you will also need to reset Ro so that it is used for subsequent resynth analysis.

save(filename, theory, substrate=None, surface=None)

Save the data and theory to a file.

scattering_factors(material, density)[source]

Returns the scattering factors associated with the material given the range of wavelengths/energies used in the probe.

show_resolution = True
simulate_data(theory, noise=2.0)

Set the data for the probe to R + eps with eps ~ normal(dR^2).

theory is (Q, R),

If the percent noise is provided, set dR to R*noise/100 before simulating. noise defaults to 2% if no dR is present.

Note that measured data estimates uncertainty from the number of counts. This means that points above the true value will have larger uncertainty than points below the true value. This bias is not captured in the simulated data.

subsample(dQ)

Select points at most every dQ.

Use this to speed up computation early in the fitting process.

This changes the data object, and is not reversible.

The current algorithm is not picking the “best” Q value, just the nearest, so if you have nearby Q points with different quality statistics (as happens in overlapped regions from spallation source measurements at different angles), then it may choose badly. Simple solutions based on the smallest relative error dR/R will be biased toward peaks, and smallest absolute error dR will be biased toward valleys.

to_dict()

Return a dictionary representation of the parameters

view = 'log'
write_data(filename, columns=('Q', 'R', 'dR'), header=None)

Save the data to a file.

header is a string with trailing n containing the file header. columns is a list of column names from Q, dQ, R, dR, L, dL, T, dT.

The default is to write Q, R, dR data.

refl1d.probe.load4(filename, keysep=':', sep=None, comment='#', name=None, intensity=1, background=0, back_absorption=1, back_reflectivity=False, Aguide=270.0, H=0, theta_offset=None, sample_broadening=None, L=None, dL=None, T=None, dT=None, dR=None, FWHM=False, radiation=None, columns=None, data_range=(None, None), resolution='normal', oversampling=None)[source]

Load in four column data Q, R, dR, dQ.

The file is loaded with bumps.data.parse_multi. keysep defaults to ‘:’ so that header data looks like JSON key: value pairs. sep is None so that the data uses white-space separated columns. comment is the standard ‘#’ comment character, used for “# key: value” lines, for commenting out data lines using “#number number number number”, and for adding comments after individual data lines. The parser isn’t very sophisticated, so be nice.

intensity is the overall beam intensity, background is the overall background level, and back_absorption is the relative intensity of data measured at negative Q compared to positive Q data. These can be values or a bumps Parameter objects.

back_reflectivity is True if reflectivity was measured through the substrate. This allows you to arrange the model from substrate to surface regardless of whether you are measuring through the substrate or reflecting off the surface.

theta_offset indicates sample alignment. In order to use theta offset you need to be able to convert from Q to wavelength and angle by providing values for the wavelength or the angle, and the associated resolution.

L, dL in Angstroms can be used to recover angle and angular resolution for monochromatic sources where wavelength is fixed and angle is varying. These values can also be stored in the file header as:

# wavelength: 4.75  # Ang
# wavelength_resolution: 0.02  # Ang (1-sigma)

T, dT in degrees can be used to recover wavelength and wavelength dispersion for time of flight sources where angle is fixed and wavelength is varying, or you can store them in the header of the file:

# angle: 2  # degrees
# angular_resolution: 0.2  # degrees (1-sigma)

If both angle and wavelength are varying in the data, you can specify a separate value for each point, such the following:

# wavelength: [1, 1.2, 1.5, 2.0, ...]
# wavelength_resolution: [0.02, 0.02, 0.02, ...]

dR can be used to replace the uncertainty estimate for R in the file with \(\Delta R = R * \text{dR}\). This allows files with only two columns, Q and R to be loaded. Note that points with dR=0 are automatically set to the minimum dR>0 in the dataset.

Instead of constants, you can provide function, dT = lambda T: f(T), dL = lambda L: f(L) or dR = lambda Q, R, dR: f(Q, R, dR) for more complex relationships (with dR() returning 1-\(\sigma\) \(\Delta R\)).

sample_broadening in degrees FWHM adds to the angular_resolution. Scale 1-\(\sigma\) rms by \(2 \surd(2 \ln 2) \approx 2.34\) to convert to FWHM.

Aguide and H are parameters for polarized beam measurements indicating the magnitude and direction of the applied field.

Polarized data is represented using a multi-section data file, with blank lines separating each section. Each section must have a polarization keyword, with value “++”, “+-”, “-+” or “–“.

FWHM is True if dQ, dT, dL are given as FWHM rather than 1-\(\sigma\). dR is always 1-\(\sigma\). sample_broadening is always FWHM.

radiation is ‘xray’ or ‘neutron’, depending on whether X-ray or neutron scattering length density calculator should be used for determining the scattering length density of a material. Default is ‘neutron’

columns is a string giving the column order in the file. Default order is “Q R dR dQ”. Note: include dR and dQ even if the file only has two or three columns, but put the missing columns at the end.

data_range indicates which data rows to use. Arguments are the same as the list slice arguments, (start, stop, step). This follows the usual semantics of list slicing, L[start:stop:step], with 0-origin indices, stop is last plus one and step optional. Use negative numbers to count from the end. Default is (None, None) for the entire data set.

resolution is ‘normal’ (default) or ‘uniform’. Use uniform if you are merging Q points from a finely stepped energy sensitive measurement.

oversampling is None or a positive integer indicating how many points to add between data point to support sparse data with denser theory (for PolarizedNeutronProbe)

refl1d.probe.make_probe(**kw)[source]

Return a reflectometry measurement object of the given resolution.

refl1d.probe.measurement_union(xs)[source]

Determine the unique (T, dT, L, dL) across all datasets.

refl1d.probe.spin_asymmetry(Qp, Rp, dRp, Qm, Rm, dRm)[source]

Compute spin asymmetry for R++, R–.

Parameters:

Qp, Rp, dRpvector

Measured ++ cross section and uncertainty.

Qm, Rm, dRmvector

Measured – cross section and uncertainty.

If dRp, dRm are None then the returned uncertainty will also be None.

Returns:

Q, SA, dSAvector

Computed spin asymmetry and uncertainty.

Algorithm:

Spin asymmetry, \(S_A\), is:

\[S_A = \frac{R_{++} - R_{--}}{R_{++} + R_{--}}\]

Uncertainty \(\Delta S_A\) follows from propagation of error:

\[\Delta S_A^2 = \frac{4(R_{++}^2\Delta R_{--}^2+R_{--}^2\Delta R_{++})} {(R_{++} + R_{--})^4}\]