Profile and IFT objects

class bioxtasraw.SASM.IFTM(p, r, err, i_orig, q_orig, err_orig, i_fit, parameters, i_extrap=[], q_extrap=[])

Bases: object

Inverse Fourier transform measurement (IFTM) object. Contains the P(r), r and error vectors, as well as the original data, the fit of the P(r) to the data, and all associated metadata.

Variables:
  • r (numpy.array) – The r vector of the P(r) function.

  • p (numpy.array) – The values of the P(r) function.

  • err (numpy.array) – The errors of the P(r) function.

  • q_orig (numpy.array) – The q vector of the input data.

  • i_orig (numpy.array) – The intensity vector of the input data.

  • err_orig (numpy.array) – The error vector of the input data.

  • i_fit (numpy.array) – The intensity vector of the fit to the input data.

  • q_extrap (numpy.array) – The q vector of the input data extrapolated to q=0.

  • i_extrap (numpy.array) – The intensity vector of the fit to the input data extrapolated to q=0.

__init__(p, r, err, i_orig, q_orig, err_orig, i_fit, parameters, i_extrap=[], q_extrap=[])

Constructor

Parameters:
  • p (numpy.array) – The input P(r) values.

  • r (numpy.array) – The input r values for the P(r) function.

  • err (numpy.array) – The input error values for the P(r) function.

  • i_orig (numpy.array) – The intensity values of the data used to do the IFT.

  • q_orig (numpy.array) – The q values of the data used to do the IFT.

  • err_orig (numpy.array) – The error values of the data used to do the IFT.

  • i_fit (numpy.array) – The intensity values of the fit of the P(r) function to the data.

  • parameters (dict) – A dictionary of the metadata. Should ontain at least {‘filename’: filename_with_no_path}

  • i_extrap (numpy.array, optional) – The intensity values of the fit of the P(r) function to the data extrapolated to q=0. If not provided, an empty array is used.

  • q_extrap (numpy.array, optional) – The q values of the input data extrapolated to q=0. If not provided, an empty array is used.

copy()

Creates a copy of the IFT.

Returns:

ift – The copied IFTM

Return type:

bioxtasraw.SASM.IFTM

extractAll()

Extracts the raw and scaled q, intensity, and error, the scale and offset values, the selected q range, and the parameters in a single dictionary.

Returns:

all_data – A dictionary with keys r_raw, p_raw, err_raw, i_orig_raw, q_orig_raw, err_orig_raw, i_fit_raw, i_extrap_raw, q_extrap_raw, and parameters, which correspond to those values from the IFTM.

Return type:

dict

getAllParameters()

Returns all of the metadata parameters associated with the IFT as a dictionary.

Returns:

parameters – The metadata associated with the IFT.

Return type:

dict

getLine()

Returns the plotted line for the P(r) function. Only used in the RAW GUI.

Returns:

line – The plotted line.

Return type:

matplotlib.lines.Line2D

getOffset()

Returns the offset for the P(r) function.

Returns:

offset – The offset.

Return type:

float

getParameter(key)

Gets a particular metadata parameter based on the provided key.

Parameters:

key (str) – A string that is a key in the parameters metadata dictionary.

Returns:

The parameter associated with the specified key. If the key is not in the parameter dictionary, None is returned.

Return type:

parameter

getScale()

Returns the scale factor for the P(r) function.

Returns:

scale – The scale factor.

Return type:

float

setAllParameters(new_parameters)

Sets the parameters dictionary, which contains the IFT metadata, to the new input value.

Parameters:

new_parameters (dict) – A dictionary containing the new parameters.

setParameter(key, value)

Sets a particular metadata parameter based on the provided key and value.

Parameters:
  • key (str) – The name of the new bit of metadata.

  • value (object) – The value of the new bit of metadata. Could be anything that is an acceptable value for a dictionary.

class bioxtasraw.SASM.SASM(i, q, err, parameters, q_err=None)

Bases: object

Small Angle Scattering Measurement (SASM) Object. Essentially a scattering profile with q, i, and uncertainty, plus a lot of metadata.

Variables:
  • q (numpy.array) – The scaled q vector, without the trimming specified by setQrange().

  • i (numpy.array) – The scaled intensity vector, without the trimming specified by setQrange().

  • err (numpy.array) – The scaled error vector, without the trimming specified by setQrange().

  • q_err (numpy.array) – The scaled q error vector, without the trimming specified by setQrange(). Typically only used with SANS data.

__init__(i, q, err, parameters, q_err=None)

Constructor

Parameters:
  • i (numpy.array) – The intensity vector.

  • q (numpy.array) – The q vector.

  • err (numpy.array) – The error vector.

  • parameters (dict) – A dictionary of metadata for the object. This should contain at least {‘filename’: filename_with_no_path}. Other reserved keys are: ‘counters’ : [(countername, value),…] Info from counter files ‘fileHeader’ : [(label, value),…] Info from the header in the loaded file

static closest(qlist, q)

A convenience function which returns the index of the nearest q point in qlist to the input q value.

Parameters:
  • qlist (np.array) – The q list to search in.

  • q (float) – The q value to search for in the qlist.

Returns:

index – The index of the q value in qlist closest to the input q.

Return type:

int

copy()

Creates a copy of the SASM, without scale information, using, for example, the scaled, trimmed intensity as the raw intensity for the new SASM.

The preferred method of copying a profile is to use copy.deepcopy on the profile.

Returns:

profile – The copied profile

Return type:

bioxtasraw.SASM.SASM

copy_no_metadata()

Creates a deep copy of the SAMS without the metadata, which will usually be faster.

extractAll()

Extracts the raw and scaled q, intensity, and error, the scale and offset values, the selected q range, and the parameters in a single dictionary.

Returns:

all_data – A dictionary with keys q_raw, i_raw, err_raw, q, i, err, scale_factor, offset_value, q_scale_factor, selected_qrange, and parameters, which correspond to those values from the SASM.

Return type:

dict

getAllParameters()

Returns all of the metadata parameters associated with the profile as a dictionary.

Returns:

parameters – The metadata associated with the profile.

Return type:

dict

getErr()

Gets the scaled, offset, trimmed error vector. Usually this is what you want to use to the get the error vector.

getI()

Gets the scaled, offset, trimmed intensity vector. Usually this is what you want to use to the get the intensity vector.

getIofQ(qref)

Gets the intensity at a specific q value (or the closest such value in the q vector).

Parameters:

qref (float) – The reference q to get the intensity at.

Returns:

intensity – The intensity at the q point nearest the provided qref.

Return type:

float

getIofQRange(q1, q2)

Gets the total integrated intensity in the q range from q1 to q2 (or the closest such values in the q vector).

Parameters:
  • q1 (float) – The starting q value in the q range

  • q2 (float) – The ending q value in the q range.

Returns:

total_intensity – The total intensity.

Return type:

float

getLine()

Returns the plotted line for the profile. Only used in the RAW GUI.

Returns:

line – The plotted line.

Return type:

matplotlib.lines.Line2D

getMeanI()

Gets the mean intensity of the intensity vector.

Returns:

mean_intensity – The mean intensity.

Return type:

float

getOffset()

Returns the offset for the profile.

Returns:

offset – The offset.

Return type:

float

getParameter(key)

Gets a particular metadata parameter based on the provided key.

Parameters:

key (str) – A string that is a key in the parameters metadata dictionary.

Returns:

The parameter associated with the specified key. If the key is not in the parameter dictionary, None is returned.

Return type:

parameter

getQ()

Gets the scaled, offset, trimmed q vector. Usually this is what you want to use to the get the q vector.

getQErr()

Gets the scaled, offset, trimmed q error vector. Usually this is what you want to use to the get the q error vector. Q errors are usually only available with SANS data. Returns None if no q error has been defined.

getQrange()

Returns the currently selected q range as described in setQrange().

Returns:

q_range – A tuple with 2 indices, the start and end of the selected q range, such that q[start:end] returns the desired q range.

Return type:

tuple

getRawErr()

Gets the raw error vector, without scaling or offset from scale() and offset() and without trimming based on setQrange().

Returns:

err_raw – The raw error vector.

Return type:

numpy.array

getRawI()

Gets the raw intensity vector, without scaling or offset from scale() and offset() and without trimming based on setQrange().

Returns:

i_raw – The raw intensity vector.

Return type:

numpy.array

getRawQ()

Gets the raw q vector, without scaling based on the scaleQ() and without trimming based on setQrange().

Returns:

q_raw – The raw q vector.

Return type:

numpy.array

getRawQErr()

Gets the raw q error vector, without scaling or offset from scale() and offset() and without trimming based on setQrange().

Returns:

q_err_raw – The raw error vector.

Return type:

numpy.array

getScale()

Returns the scale factor for the profile.

Returns:

scale – The scale factor.

Return type:

float

getTotalI()

Gets the total integrated intensity of the intensity vector.

Returns:

total_intensity – The total intensity.

Return type:

float

offset(offset_value)

Applies an absolute offset to the profile intensity. For example, if the offset is 1, then all the the intensities in the profile are increased by 1. Offset supersedes the existing offset, so if the current offset is 1, and an offset_value of 2 is provided, the resulting offset is 2 (not 3).

Parameters:

offset_value (float) – The offset to be applied to the profile intensity.

offsetRawIntensity(offset)

Offsets the raw intensity and error values. These are the intensity and error values without any offset applied. The raw offset factor is not tracked, so this cannot easily be undone, unlike the :func:scaleQ function.

Parameters:

offset (float) – The offset to be applied to the raw profile intensity and error values.

removeParameter(key)

Removes a particular metadata parameter based on the provided key.

Parameters:

key (str) – A string that is a key in the parameters metadata dictionary.

removeZingers(start_idx=0, window_length=10, stds=4.0)

Removes spikes (zingers) from radially averaged data based on smoothing of the intensity profile.

Parameters:
  • start_idx (int) – The initial index in the intensity to start the dezingering process.

  • window_length (int) – The size of the window used to search for an replace spikes, as the number of intensity points.

  • stds (The standard deviation threshold used to detect spikes.) –

reset()

Removes scale and offset values from the intensity, uncertainty, and q.

scale(scale_factor)

Applies an absolute scale to the profile intensity. The scale factor supersedes the existing scale factor. For example, suppose the scale factor is currently 2. If a scale of 4 is provided, the resulting scale factor is 4 (not 8). Scale factors are always positive.

Parameters:

scale_factor (float) – The scale factor to be applied to the the profile intensity and uncertainty.

scaleQ(q_scale_factor)

Scales the profile q values by a factor. The scale factor supersedes the existing scale factor. For example, suppose the scale factor is currently 2. If a scale of 4 is provided, the resulting scale factor if 4 (not 8). Useful for converting between 1/A and 1/nm, for example.

Parameters:

q_scale_factor (float) – The scale factor to be applied to the profile q values.

scaleRawIntensity(scale)

Scales the raw intensity and error values. These are the intensity and error values without any scale applied. The raw scale factor is not tracked, so this cannot easily be undone, unlike the :func:scaleQ function.

Parameters:

scale (float) – The scale factor to be applied to the raw profile intensity and error values.

scaleRawQ(scale_factor)

Scales the raw q values. These are the q values without any scale applied. The raw q scale factor is not tracked, so this cannot easily be undone, unlike the :func:scaleQ function.

Parameters:

scale_factor (float) – The scale factor to be applied to the raw profile q values.

scaleRelative(relscale)

Applies a relative scale to the profile intensity. If the scale factor is currently 1, then this is the same as scale(). Otherwise, this scales relative to the current scale factor. For example, suppose the scale factor is currently 2. If a relative scale of 2 is provided, the resulting scale factor if 4. Scale factors are always positive.

Parameters:

relscale (float) – The relative scale factor to be applied to the the profile intensity and uncertainty.

scaleRelativeQ(relscale)

Applies a relative scale to the profile q. If the scale factor is currently 1, then this is the same as scale(). Otherwise, this scales relative to the current scale factor. For example, suppose the scale factor is currently 2. If a relative scale of 2 is provided, the resulting scale factor if 4.

Parameters:

relscale (float) – The relative scale factor to be applied to the the profile intensity and uncertainty.

setAllParameters(new_parameters)

Sets the parameters dictionary, which contains the profile metadata, to the new input value.

Parameters:

new_parameters (dict) – A dictionary containing the new parameters.

setParameter(key, value)

Sets a particular metadata parameter based on the provided key and value.

Parameters:
  • key (str) – The name of the new bit of metadata.

  • value (object) – The value of the new bit of metadata. Could be anything that is an acceptable value for a dictionary.

setQrange(qrange)

Sets the q range used for the profile. Useful for trimming leading or trailing values of the q profile that are not useful data.

Parameters:

qrange (tuple or list) – A tuple or list with two items. The first item is the starting index of the q vector to be used, the second item is the ending index of the q vector to be used, such that q[start:end] returns the desired q range.

setRawErr(new_raw_err)

Sets the raw error vector. Will overwrite whatever error vector is already in the object! Typically only used during calibration.

Parameters:

new_raw_err (numpy.array) – The new error vector.

setRawI(new_raw_i)

Sets the raw q vector. Will overwrite whatever q vector is already in the object! Typically only used during calibration.

Parameters:

new_raw_i (numpy.array) – The new intensity vector.

setRawQ(new_raw_q)

Sets the raw intensity vector. Will overwrite whatever intensity vector is already in the object! Typically only used during calibration.

Parameters:

new_raw_q (numpy.array) – The new q vector.

setRawQErr(new_raw_q_err)

Sets the raw q error vector. Will overwrite whatever error vector is already in the object! Typically only used during calibration. Q errors are typically only found in SANS data, and currently are only carried and written out with the data set, they are not used in any processing.

Parameters:

new_raw_err (numpy.array) – The new error vector.

setScaleValues(scale_factor, offset_value, q_scale_factor)

A convenience method that lets you set the scale offset, and q scale values all at once.

Parameters:
  • scale_factor (float) – The scale factor to be applied to the the profile intensity and uncertainty.

  • offset_value (float) – The offset to be applied to the profile intensity.

  • q_scale_factor (float) – The scale factor to be applied to the profile q values.