experiment - Reflectivity fitness function

Experiment

Theory calculator.

ExperimentBase

MixedExperiment

Support composite sample reflectivity measurements.

nice

Fix v to a value with a given number of digits of precision

plot_sample

Quick plot of a reflectivity sample and the corresponding reflectivity.

Experiment definition

An experiment combines the sample definition with a measurement probe to create a fittable reflectometry model.

class refl1d.experiment.Experiment(sample=None, probe=None, name=None, roughness_limit=0, dz=None, dA=None, step_interfaces=None, smoothness=None, interpolation=0)[source]

Bases: ExperimentBase

Theory calculator. Associates sample with data, Sample plus data. Associate sample with measurement.

The model calculator is specific to the particular measurement technique that was applied to the model.

Measurement properties:

probe is the measuring probe

Sample properties:

sample is the model sample step_interfaces use slabs to approximate gaussian interfaces roughness_limit limit the roughness based on layer thickness dz minimum step size for computed profile steps in Angstroms dA discretization condition for computed profiles

If step_interfaces is True, then approximate the interface using microslabs with step size dz. The microslabs extend throughout the whole profile, both the interfaces and the bulk; a value for dA should be specified to save computation time. If False, then use the Nevot-Croce analytic expression for the interface between slabs.

The roughness_limit value should be reasonably large (e.g., 2.5 or above) to make sure that the Nevot-Croce reflectivity calculation matches the calculation of the displayed profile. Use a value of 0 if you want no limits on the roughness, but be aware that the displayed profile may not reflect the actual scattering densities in the material.

The dz step size sets the size of the slabs for non-uniform profiles. Using the relation d = 2 pi / Q_max, we use a default step size of d/20 rounded to two digits, with 5 Å as the maximum default. For simultaneous fitting you may want to set dz explicitly using to round(pi/Q_max/10, 1) so that all models use the same step size.

The dA condition measures the uncertainty in scattering materials allowed when combining the steps of a non-uniform profile into slabs. Specifically, the area of the box containing the minimum and the maximum of the non-uniform profile within the slab will be smaller than dA. A dA of 10 gives coarse slabs. If dA is not provided then each profile step forms its own slab. The dA condition will also apply to the slab approximation to the interfaces.

interpolation indicates the number of points to plot in between existing points.

smoothness DEPRECATED This parameter is not used.

amplitude(resolution=False, interpolation=0)[source]

Calculate reflectivity amplitude at the probe points.

format_parameters()
interpolation = 0
is_reset()

Returns True if a model reset was triggered.

property ismagnetic

True if experiment contains magnetic materials

magnetic_slabs()[source]
magnetic_smooth_profile(dz=0.1)[source]

Return the nuclear and magnetic scattering potential for the sample.

magnetic_step_profile()[source]

Return the nuclear and magnetic scattering potential for the sample.

property name
nllf()

Return the -log(P(data|model)).

Using the assumption that data uncertainty is uncorrelated, with measurements normally distributed with mean R and variance dR**2, this is just sum( resid**2/2 + log(2*pi*dR**2)/2 ).

The current version drops the constant term, sum(log(2*pi*dR**2)/2).

numpoints()
parameters()[source]

Fittable parameters to sample and probe

penalty()[source]
plot(plot_shift=None, profile_shift=None, view=None)
plot_profile(plot_shift=None)[source]
plot_reflectivity(show_resolution=False, view=None, plot_shift=None)
probe: probe.Probe = None
profile_shift = 0
reflectivity(resolution=True, interpolation=0)[source]

Calculate predicted reflectivity.

If resolution is true include resolution effects.

residuals()
restore_data()

Restore original data after resynthesis.

resynth_data()

Resynthesize data with noise from the uncertainty estimates.

save(basename)
save_json(basename)

Save the experiment as a json file

save_profile(basename)
save_refl(basename)
save_staj(basename)[source]
simulate_data(noise=2.0)

Simulate a random data set for the model.

This sets R and dR according to the noise level given.

Parameters:

noise: float or array or None | %

dR/R uncertainty as a percentage. If noise is set to None, then use dR from the data if present, otherwise default to 2%.

slabs()[source]

Return the slab thickness, roughness, rho, irho for the rendered model.

Note

Roughness is for the top of the layer.

smooth_profile(dz=0.1)[source]

Return the scattering potential for the sample.

If dz is not given, use dz = 0.1 A.

step_profile()[source]

Return the step scattering potential for the sample, ignoring interfaces.

to_dict()[source]
update()

Called when any parameter in the model is changed.

This signals that the entire model needs to be recalculated.

update_composition()

When the model composition has changed, we need to lookup the scattering factors for the new model. This is only needed when an existing chemical formula is modified; new and deleted formulas will be handled automatically.

write_data(filename, **kw)

Save simulated data to a file

class refl1d.experiment.ExperimentBase[source]

Bases: object

format_parameters()[source]
interpolation = 0
is_reset()[source]

Returns True if a model reset was triggered.

magnetic_slabs()[source]
magnetic_step_profile()[source]
property name
nllf()[source]

Return the -log(P(data|model)).

Using the assumption that data uncertainty is uncorrelated, with measurements normally distributed with mean R and variance dR**2, this is just sum( resid**2/2 + log(2*pi*dR**2)/2 ).

The current version drops the constant term, sum(log(2*pi*dR**2)/2).

numpoints()[source]
parameters()[source]
plot(plot_shift=None, profile_shift=None, view=None)[source]
plot_profile(plot_shift=0.0)[source]
plot_reflectivity(show_resolution=False, view=None, plot_shift=None)[source]
probe: probe.Probe = None
reflectivity(resolution=True, interpolation=0)[source]
residuals()[source]
restore_data()[source]

Restore original data after resynthesis.

resynth_data()[source]

Resynthesize data with noise from the uncertainty estimates.

save(basename)[source]
save_json(basename)[source]

Save the experiment as a json file

save_profile(basename)[source]
save_refl(basename)[source]
simulate_data(noise=2.0)[source]

Simulate a random data set for the model.

This sets R and dR according to the noise level given.

Parameters:

noise: float or array or None | %

dR/R uncertainty as a percentage. If noise is set to None, then use dR from the data if present, otherwise default to 2%.

slabs()[source]
smooth_profile(dz=0.1)[source]
step_profile()[source]
to_dict()[source]
update()[source]

Called when any parameter in the model is changed.

This signals that the entire model needs to be recalculated.

update_composition()[source]

When the model composition has changed, we need to lookup the scattering factors for the new model. This is only needed when an existing chemical formula is modified; new and deleted formulas will be handled automatically.

write_data(filename, **kw)[source]

Save simulated data to a file

class refl1d.experiment.MixedExperiment(samples=None, ratio=None, probe=None, name=None, coherent=False, interpolation=0, **kw)[source]

Bases: ExperimentBase

Support composite sample reflectivity measurements.

Sometimes the sample you are measuring is not uniform. For example, you may have one portion of you polymer brush sample where the brushes are close packed and able to stay upright, whereas a different section of the sample has the brushes lying flat. Constructing two sample models, one with brushes upright and one with brushes flat, and adding the reflectivity incoherently, you can then fit the ratio of upright to flat.

samples the layer stacks making up the models ratio a list of parameters, such as [3, 1] for a 3:1 ratio probe the measurement to be fitted or simulated

coherent is True if the length scale of the domains is less than the coherence length of the neutron, or false otherwise.

Statistics such as the cost functions for the individual profiles can be accessed from the underlying experiments using composite.parts[i] for the various samples.

amplitude(resolution=False)[source]
format_parameters()
interpolation = 0
is_reset()

Returns True if a model reset was triggered.

magnetic_slabs()
magnetic_step_profile()
property name
nllf()

Return the -log(P(data|model)).

Using the assumption that data uncertainty is uncorrelated, with measurements normally distributed with mean R and variance dR**2, this is just sum( resid**2/2 + log(2*pi*dR**2)/2 ).

The current version drops the constant term, sum(log(2*pi*dR**2)/2).

numpoints()
parameters()[source]
penalty()[source]
plot(plot_shift=None, profile_shift=None, view=None)
plot_profile(plot_shift=None)[source]
plot_reflectivity(show_resolution=False, view=None, plot_shift=None)
probe: probe.Probe = None
reflectivity(resolution=True, interpolation=0)[source]

Calculate predicted reflectivity.

This will be the weighted sum of the reflectivity from the individual systems. If coherent is set, then the coherent sum will be used, otherwise the incoherent sum will be used.

If resolution is true include resolution effects.

interpolation is the number of theory points to show between data points.

residuals()
restore_data()

Restore original data after resynthesis.

resynth_data()

Resynthesize data with noise from the uncertainty estimates.

save(basename)
save_json(basename)

Save the experiment as a json file

save_profile(basename)[source]
save_refl(basename)
save_staj(basename)[source]
simulate_data(noise=2.0)

Simulate a random data set for the model.

This sets R and dR according to the noise level given.

Parameters:

noise: float or array or None | %

dR/R uncertainty as a percentage. If noise is set to None, then use dR from the data if present, otherwise default to 2%.

slabs()
smooth_profile(dz=0.1)
step_profile()
to_dict()[source]
update()[source]

Called when any parameter in the model is changed.

This signals that the entire model needs to be recalculated.

update_composition()

When the model composition has changed, we need to lookup the scattering factors for the new model. This is only needed when an existing chemical formula is modified; new and deleted formulas will be handled automatically.

write_data(filename, **kw)

Save simulated data to a file

refl1d.experiment.nice(v, digits=2)[source]

Fix v to a value with a given number of digits of precision

refl1d.experiment.plot_sample(sample, instrument=None, roughness_limit=0)[source]

Quick plot of a reflectivity sample and the corresponding reflectivity.