API Reference

SQWDndBase class

sqw.@SQWDnDBase.IX_dataset_1d(w)

Convert 1D sqw object into IX_dataset_1d

>> wout = IX_dataset_1d (w)

sqw.@SQWDnDBase.IX_dataset_2d(w)

Convert 2D sqw object into IX_dataset_2d

>> wout = IX_dataset_2d (w)

sqw.@SQWDnDBase.IX_dataset_3d(w)

Convert 3D sqw object into IX_dataset_3d

>> wout = IX_dataset_3d (w)

class sqw.@SQWDnDBase.SQWDnDBase

SQWDnDBase Abstract SQW/DnD object base class Abstract class defining common API and atrributes of the SQW and DnD objects

sqw.@SQWDnDBase.acos(w1)

Implements acos(w1) for objects

>> w = acos(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.acosh(w1)

Implements acosh(w1) for objects

>> w = acosh(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.acot(w1)

Implements acot(w1) for objects

>> w = acot(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.acoth(w1)

Implements acoth(w1) for objects

>> w = acoth(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.acsc(w1)

Implements acsc(w1) for objects

>> w = acsc(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.acsch(w1)

Implements acsch(w1) for objects

>> w = acsch(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.asec(w1)

Implements asec(w1) for objects

>> w = asec(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.asech(w1)

Implements asech(w1) for objects

>> w = asech(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.asin(w1)

Implements asin(w1) for objects

>> w = asin(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.asinh(w1)

Implements asinh(w1) for objects

>> w = asinh(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.atan(w1)

Implements atan(w1) for objects

>> w = atan(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.atanh(w1)

Implements atanh(w1) for objects

>> w = atanh(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.binary_op_manager(w1, w2, binary_op)

Implements a binary operation for objects with a signal and a variance array.

>> w = binary_op_manager(w1, w2, binary_op)

All binary operations on Matlab double arrays are permitted (+, -, *, /, ) and are applied element by element to the signal and variance arrays.

w1, w2 Objects on which the binary operation is to be performed.

One of these can be a Matlab double (i.e. double precision) array, in which case the variance array is taken to be zero.

If w1, w2 are scalar objects with the same signal array sizes: - The operation is performed element-by-element.

If one of w1 or w2 is a double array (and the other is a scalar object): - If a scalar, apply to each element of the object signal. - If it is an array of the same size as the object signal

array, apply the operation element by element.

If one or both of w1 and w2 are arrays of objects: - If objects have same array sizes, the binary operation is

applied object element-by-object element.

  • If one of the objects is scalar (i.e. only one object),

then it is applied by the binary operation to each object in the other array.

If one of w1, w2 is an array of objects and the other is a double array: - If the double is a scalar, it is applied to every object

in the array.

  • If the double is an array with the same size as the object

array, then each element is applied as a scalar to the corresponding object in the object array.

  • If the double is an array with larger size than the object

array, then the array is resolved into a stack of arrays, where the stack has the same size as the object array, and the each array in the stack is applied to the corresponding object in the object array. [Note that for this operation to be valid, each object must have the same signal array size.]

binary_op Function handle to a binary operation. All binary operations

on Matlab double or single arrays are permitted (+, -, *, /, ).

w Output object or array of objects.

NOTES: This is a generic method - works for any class (including sigvar) so long as the methods below are defined on that class.

Requires that objects have the following methods to find the size of the public signal and variance arrays, create a sigvar object from those arrays, and set them from another sigvar object.

>> sz = sigvar_size(obj) % Returns size of public signal and variance

% arrays

>> w = sigvar(obj) % Create a sigvar object from the public

% signal and variance arrays

>> obj = sigvar_set(obj,w) % Set signal and variance in an object from

% those in a sigvar object

sqw.@SQWDnDBase.calculate_q_bins(win)

Calculate qh,qk,ql,en for the centres of the bins of an n-dimensional sqw or dnd dataset

>> [q,en]=calculate_q_bins(win)

win Input sqw/dnd object

q Components of momentum (in rlu) for each bin in the dataset for a single energy bin

Arrays are packaged as cell array of column vectors for convenience with fitting routines etc.

i.e. q{1}=qh, q{2}=qk, q{3}=ql

en Column vector of energy bin centres. If energy was an integration axis, then returns the

centre of the energy integration range

sqw.@SQWDnDBase.calculate_qw_bins(win, optstr)

Calculate qh,qk,ql,en for the centres of the bins of an n-dimensional sqw or dnd dataset

>> qw=calculate_qw_bins(win) >> qw=calculate_qw_bins(win,’boundaries’) >> qw=calculate_qw_bins(win,’edges’)

win Input sqw or dnd object

Optional arguments: ‘boundaries’ Return qh,qk,ql,en at verticies of bins, not centres ‘edges’ Return qh,qk,ql,en at verticies of the hyper cuboid that

encloses the plot axes

qw Components of momentum (in rlu) and energy for each bin in

the dataset Arrays are packaged as cell array of column vectors for convenience with fitting routines etc.

i.e. qw{1}=qh, qw{2}=qk, qw{3}=ql, qw{4}=en

Note that the centre of the integration range is used in

the calculation of qh,qk,ql,en even with the options ‘boundaries’ or ‘edges’

If one or both of the integration ranges is infinite, then

the value of the corresponding coordinate is taken as zero.

sqw.@SQWDnDBase.change_crystal(varargin)

Change the crystal lattice and orientation of an sqw object or array of objects

Most commonly:

>> wout = change_crystal (w, rlu_corr) % change lattice parameters and orientation

OR

>> wout = change_crystal (w, alatt) % change just length of lattice vectors >> wout = change_crystal (w, alatt, angdeg) % change all lattice parameters >> wout = change_crystal (w, alatt, angdeg, rotmat) % change lattice parameters and orientation >> wout = change_crystal (w, alatt, angdeg, u, v) % change lattice parameters and redefine u, v

w Input sqw object

rlu_corr Matrix to convert notional rlu in the current crystal lattice to

the rlu in the the new crystal lattice together with any re-orientation of the crystal. The matrix is defined by the matrix:

qhkl(i) = rlu_corr(i,j) * qhkl_0(j)

This matrix can be obtained from refining the lattice and

orientation with the function refine_crystal (type >> help refine_crystal for more details).

OR

alatt New lattice parameters [a,b,c] (Angstroms) angdeg New lattice angles [alf,bet,gam] (degrees) rotmat Rotation matrix that relates crystal Cartesian coordinate frame of the new

lattice as a rotation of the current crystal frame. Orthonormal coordinates in the two frames are related by

v_new(i)= rotmat(i,j)*v_current(j)

u, v Redefine the two vectors that were used to determine the scattering plane

These are the vectors at whatever disorientation angles dpsi, gl, gs (which cannot be changed).

wout Output sqw object with changed crystal lattice parameters and orientation

NOTE

The input data set(s) can be reset to their original orientation by inverting the input data e.g.

  • call with inv(rlu_corr)

  • call with the original alatt, angdeg, u and v

sqw.@SQWDnDBase.compact(win)

Squeezes the data range in an sqw or dnd object to eliminate empty bins

Syntax:

>> wout = compact(win)

win Input object

wout Output object, with length of axes reduced to yield the

smallest cuboid that contains the non-empty bins.

sqw.@SQWDnDBase.convert_bins_for_shoelace(win, wref)

Converts data in the d2d object win into appropriate format for the shoelace rebinning function (i.e. bin corners rather than bin boundaries). The co-ordinate system is defined by the d2d object wref. If wref is empty then the co-ordinate system remains unchanged.

RAE 23/9/09

sqw.@SQWDnDBase.cos(w1)

Implements cos(w1) for objects

>> w = cos(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.cosh(w1)

Implements cosh(w1) for objects

>> w = cosh(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.cot(w1)

Implements cot(w1) for objects

>> w = cot(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.coth(w1)

Implements coth(w1) for objects

>> w = coth(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.csc(w1)

Implements csc(w1) for objects

>> w = csc(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.csch(w1)

Implements csch(w1) for objects

>> w = csch(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.da(w, varargin)

Draw an area plot of a 2D sqw dataset or array of datasets

>> da(w) >> da(w,xlo,xhi) >> da(w,xlo,xhi,ylo,yhi) >> da(w,xlo,xhi,ylo,yhi,zlo,zhi)

Advanced use:

>> da(w,…,’name’,fig_name) % Draw with name = fig_name

>> da(w,…,’-noaspect’) % Do not change aspect ratio

% according to data axes unit lengths

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = da(w,…)

sqw.@SQWDnDBase.dd(w, varargin)

Draws a plot of markers, error bars and lines of a 1D sqw object or array of objects

>> dd(w) >> dd(w,xlo,xhi) >> dd(w,xlo,xhi,ylo,yhi)

Advanced use:

>> dd(w,…,’name’,fig_name) % draw with name = fig_name

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = dd(w,…)

sqw.@SQWDnDBase.de(w, varargin)

Draws a plot of error bars of a 1D sqw object or array of objects

>> de(w) >> de(w,xlo,xhi) >> de(w,xlo,xhi,ylo,yhi)

Advanced use:

>> de(w,…,’name’,fig_name) % draw with name = fig_name

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = de(w,…)

sqw.@SQWDnDBase.dh(w, varargin)

Draws a histogram plot of a 1D sqw object or array of objects

>> dh(w) >> dh(w,xlo,xhi) >> dh(w,xlo,xhi,ylo,yhi)

Advanced use:

>> dh(w,…,’name’,fig_name) % draw with name = fig_name

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = dh(w,…)

sqw.@SQWDnDBase.dimensions(w)

Find number of dimensions and extent along each dimension of the signal arrays. - If 0D sqw object, nd=0, sz=zeros(1,0) (nb: []==zeros(0,0)) - if 1D sqw object, nd=1, sz=n1 - If 2D sqw object, nd=2, sz=[n1,n2] - If 3D sqw object, nd=3, sz=[n1,n2,n3] even if n3=1 - If 4D sqw object, nd=4, sz=[n1,n2,n3,n4] even if n4=1

The convention is that size(sz)=1 x ndim

>> [nd,sz]=dimensions(w)

sqw.@SQWDnDBase.dimensions_match(w, nd_ref)

Check that the dimensions of an array of sqw objects are all the same

>> [ok,mess]=dimensions_match(w) >> [ok,mess]=dimensions_match(w,nref)

w sqw object or array of objects nd_ref [optional] If not given, check all sqw objects in the array

have the same dimensionality. If given, check that they match this dimensionality

ok True if all have the smae dimensionality (and match nref, if given) mess Empty if ok==true; error message if not nd Dimensionality

sqw.@SQWDnDBase.disp2sqw_eval(win, dispreln, pars, fwhh, varargin)

Calculate sqw for a model scattering function

>> wout = disp2sqw_eval(win,dispreln,pars,fwhh,varargin)

win Dataset, or array of datasets, that provides the axes and points

for the calculation

dispreln Handle to function that calculates the dispersion relation w(Q) and

spectral weight, s(Q) Must have form:

[w,s] = dispreln (qh,qk,ql,p)

where
qh,qk,ql Arrays containing the coordinates of a set of points

in reciprocal lattice units

p Vector of parameters needed by dispersion function

e.g. [A,js,gam] as intensity, exchange, lifetime

w Array of corresponding energies, or, if more than

one dispersion relation, a cell array of arrays.

s Array of spectral weights, or, if more than

one dispersion relation, a cell array of arrays.

More general form is:

[w,s] = dispreln (qh,qk,ql,p,c1,c2,..)

where
p Typically a vector of parameters that we might want

to fit in a least-squares algorithm

c1,c2,… Other constant parameters e.g. file name for look-up

table.

pars Arguments needed by the function. Most commonly, a vector of parameter

values e.g. [A,js,gam] as intensity, exchange, lifetime. If a more general set of parameters is required by the function, then package these into a cell array and pass that as pars. In the example above then pars = {p, c1, c2, …}

fwhh Parametrizes the resolution function. There are three

possible input values of fwhh:

double A single FWHM value determines the FWHM of the

Gaussian resolution function

function_handle A function that produces the FWHM value as a

function of energy transfer, it has to have the following simple header (where omega can be a row vector of energies:

dE = resfun(omega)

function_handle A function handle of a function with two input
parameters with the following header:

I = shapefun(Emat,omega)

where Emat is a matrix with dimensions of [nQ nE] and omega is a column vector with nQ elements. The shapefun produces a peakshape for every Q point centered at the given omega and normalized to one. The output I has the same dimensions as the input Emat.

Optional arguments: (varargin)

‘-al[l]’ Requests that the calculated sqw be returned over

the whole of the domain of the input dataset. If not given, then the function will be returned only at those points of the dataset that contain data. Applies only to input with no pixel information - it is ignored if full sqw object.

‘-av[erage]’ Requests that the calculated sqw be computed for the

average values of h,k,l of the pixels in a bin, not for each pixel individually. Reduces cost of expensive calculations.

Applies only to the case of sqw object with pixel information - it is

ignored if dnd type object.

wout Output dataset or array of datasets

sqw.@SQWDnDBase.dispersion(win, dispreln, pars)

Calculate dispersion relation for dataset or array of datasets.

>> wdisp = dispersion (win, dispreln, p) % dispersion only >> [wdisp,weight] = dispersion (win, dispreln, p) % dispersion and spectral weight

The output dataset (or array of data sets), wdisp, will retain only the Q axes, and the signal array(s) will contain the values of energy along the Q axes. If the dispersion relation returns the spectral weight, this will be placed in the error array (actually the square of the spectral weight is put in the error array). In the case when the dispersion has been calculated on a plane in momentum (i.e. wdisp is IX_datset_2d) then the plot function ps2 (for plot_surface2)

>> ps2(wdisp)

will plot a surface with the z axis as energy and coloured according to the spectral weight.

The dispersion relation is calculated at the bin centres (that is, the individual pixel information in a sqw input object is not used).

If the function that calculates dispersion relations produces more than one branch, then in the case of a single input dataset the output will be an array of datasets, one for each branch. If the input is an array of datasets, then only the first dispersion branch will be returned, so there is one output dataset per input dataset.

win Dataset that provides the axes and points for the calculation

If one of the plot axes is energy transfer, then the output dataset

will have dimensionality one less than the input dataset

dispreln Handle to function that calculates the dispersion relation w(Q)
Must have form:

[w,s] = dispreln (qh,qk,ql,p)

where
qh,qk,ql Arrays containing the coordinates of a set of points

in reciprocal lattice units

p Vector of parameters needed by dispersion function

e.g. [A,js,gam] as intensity, exchange, lifetime

w Array of corresponding energies, or, if more than

one dispersion relation, a cell array of arrays.

s Array of spectral weights, or, if more than

one dispersion relation, a cell array of arrays.

More general form is:

[w,s] = dispreln (qh,qk,ql,p,c1,c2,..)

where
p Typically a vector of parameters that we might want

to fit in a least-squares algorithm

c1,c2,… Other constant parameters e.g. file name for look-up

table.

p Arguments needed by the function. Most commonly, a vector of parameter

values e.g. [A,js,gam] as intensity, exchange, lifetime. If a more general set of parameters is required by the function, then package these into a cell array and pass that as pars. In the example above then pars = {p, c1, c2, …}

wdisp Output dataset or array of datasets of the same type as the input argument.

The output dataset (or array of data sets) will retain only the Q axes, the

the signal array(s) will contain the values of energy along the Q axes, and the error array will contain the square of the spectral weight.

If the function that calculates dispersion relations produces more than one

branch, then in the case of a single input dataset the output will be an array of datasets, one for each branch. If the input is an array of datasets, then only the first dispersion branch will be returned, so there is one output dataset per input dataset.

weight Mirror output: the signal is the spectral weight, and the error array

contains the square of the frequency.

e.g. If win is a 2D dataset with Q and E axes, then wdisp is a 1D dataset

with just the Q axis

sqw.@SQWDnDBase.dl(w, varargin)

Draws a line plot of a 1D sqw object or array of objects

>> dl(w) >> dl(w,xlo,xhi) >> dl(w,xlo,xhi,ylo,yhi)

Advanced use:

>> dl(w,…,’name’,fig_name) % draw with name = fig_name

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = dl(w,…)

sqw.@SQWDnDBase.dm(w, varargin)

Draws a marker plot of a 1D sqw object or array of objects

>> dm(w) >> dm(w,xlo,xhi) >> dm(w,xlo,xhi,ylo,yhi)

Advanced use:

>> dm(w,…,’name’,fig_name) % draw with name = fig_name

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = dm(w,…)

sqw.@SQWDnDBase.dp(w, varargin)

Draws a plot of markers and error bars for a 1D sqw object or array of objects

>> dp(w) >> dp(w,xlo,xhi) >> dp(w,xlo,xhi,ylo,yhi)

Advanced use:

>> dp(w,…,’name’,fig_name) % draw with name = fig_name

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = dp(w,…)

sqw.@SQWDnDBase.ds(w, varargin)

Draw a surface plot of a 2D sqw dataset or array of datasets

>> ds(w) >> ds(w,xlo,xhi) >> ds(w,xlo,xhi,ylo,yhi) >> ds(w,xlo,xhi,ylo,yhi,zlo,zhi)

Advanced use:

>> ds(w,…,’name’,fig_name) % Draw with name = fig_name

>> ds(w,…,’-noaspect’) % Do not change aspect ratio

% according to data axes unit lengths

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = ds(w,…)

sqw.@SQWDnDBase.ds2(w, varargin)

Draw a surface plot of a 2D sqw dataset or array of datasets

>> ds2(w) % Use error bars to set colour scale >> ds2(w,wc) % Signal in wc sets colour scale

% wc can be any object with a signal array with same % size as w, e.g. sqw object, IX_dataset_2d object, or % a numeric array. % - If w is an array of objects, then wc must contain % the same number of objects. % - If wc is a numeric array then w must be a scalar % object.

>> ds2(…,xlo,xhi) >> ds2(…,xlo,xhi,ylo,yhi) >> ds2(…,xlo,xhi,ylo,yhi,zlo,zhi)

Differs from ds in that the signal sets the z axis, and the colouring is set by the error bars, or by another object. This enables two related functions to be plotted (e.g. dispersion relation where the ‘signal’ array holds the energy and the error array holds the spectral weight).

Advanced use:

>> ds2(w,…,’name’,fig_name) % Draw with name = fig_name

>> ds2(w,…,’-noaspect’) % Do not change aspect ratio

% according to data axes unit lengths

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = ds2(…)

sqw.@SQWDnDBase.equal_to_tol(w1, w2, varargin)

Check if two sqw objects are equal to a given tolerance

>> ok = equal_to_tol (a, b) >> ok = equal_to_tol (a, b, tol) >> ok = equal_to_tol (…, keyword1, val1, keyword2, val2,…) >> [ok, mess] = equal_to_tol (…)

Class specific version of the generic equal_to_tol that by default
  1. assumes NaN are equivalent (see option ‘nan_equal’), and

  2. ignores the order of pixels within a bin as the order is irrelevant (change the default with option ‘reorder’)

In addition, it is possible to check the contents of just a random fraction of non-empty bins (see option ‘fraction’) in order to speed up the comparison of large objects.

w1,w2 Test objects (scalar objects, or arrays of objects with same sizes)

tol Tolerance criterion for numeric arrays (Default: [0,0] i.e. equality)
It has the form: [abs_tol, rel_tol] where

abs_tol absolute tolerance (>=0; if =0 equality required) rel_tol relative tolerance (>=0; if =0 equality required)

If either criterion is satified then equality within tolerance is accepted.

Examples:

[1e-4, 1e-6] absolute 1e-4 or relative 1e-6 required [1e-4, 0] absolute 1e-4 required [0, 1e-6] relative 1e-6 required [0, 0] equality required 0 equivalent to [0,0]

For backwards compatibility, a scalar tolerance can be given where the sign determines absolute or relative tolerance

+ve : absolute tolerance abserr = abs(a-b) -ve : relative tolerance relerr = abs(a-b)/max(abs(a),abs(b))

Examples:

1e-4 absolute tolerance, equivalent to [1e-4, 0] -1e-6 relative tolerance, equivalent to [0, 1e-6]

[To apply an absolute as well as a relative tolerance with a

scalar negative value, set the value of the legacy keyword

‘min_denominator’ (see below)]

Valid keywords are:

‘nan_equal’ Treat NaNs as equal (true or false; default=true)

‘ignore_str’ Ignore the length and content of strings or cell arrays

of strings (true or false; default=false)

‘reorder’ Ignore the order of pixels within each bin
(true or false; default=true)

Only applies if sqw-type object

‘fraction’ Compare pixels in only a fraction of the non-empty bins
(0<= fracton <= 1; default=1 i.e. test all bins)

Only applies if sqw-type object

The reorder and fraction options are available because the order of the

pixels within the pix array for a given bin is unimportant. Reordering takes time, however, so the option to test on a few bins is given.

sqw.@SQWDnDBase.exp(w1)

Implements exp(w1) for objects

>> w = exp(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.func_eval(win, func_handle, pars, varargin)

Evaluate a function at the plotting bin centres of sqw object or array of sqw object Syntax:

>> wout = func_eval (win, func_handle, pars) >> wout = func_eval (win, func_handle, pars, [‘all’]) >> wout = func_eval (win, func_handle, pars, ‘outfile’, ‘output.sqw’)

If function is called on sqw-type object (i.e. has pixels), the pixels’ signal is also modified and evaluated

win Dataset or array of datasets; the function will be evaluated

at the bin centres along the plot axes

func_handle Handle to the function to be evaluated at the bin centres
Must have form:

y = my_function (x1,x2,… ,xn,pars)

or, more generally:

y = my_function (x1,x2,… ,xn,pars,c1,c2,…)

  • x1,x2,.xn Arrays of x coordinates along each of the n dimensions

  • pars Parameters needed by the function

  • c1,c2,… Any further arguments needed by the function e.g.

    they could be the filenames of lookup tables for resolution effects)

e.g. y=gauss2d(x1,x2,[ht,x0,sig])

y=gauss4d(x1,x2,x3,x4,[ht,x1_0,x2_0,x3_0,x4_0,sig1,sig2,sig3,sig4])

pars Arguments needed by the function.
  • Most commonly just a numeric array of parameters

  • If a more general set of parameters is needed by the function, then wrap as a cell array {pars, c1, c2, …}

Parameters

file (outfile     If present, the output of func_eval will be written to the) – of the given name/path. If numel(win) > 1, outfile must be omitted or a cell array of file paths with equal number of elements as win.

Additional allowed options:
‘all’ Requests that the calculated function be returned over

the whole of the domain of the input dataset. If not given, then the function will be returned only at those points of the dataset that contain data.

Applies only to input with no pixel information - this option is ignored if

the input is a full sqw object.

wout Output objects or array of objects

e.g.

>> wout = func_eval (w, @gauss4d, [ht,x1_0,x2_0,x3_0,x4_0,sig1,sig2,sig3,sig4])

where the function gauss appears on the matlab path

function y = gauss4d (x1, x2, x3, x4, pars) y = (pars(1)/(sig*sqrt(2*pi))) * …

sqw.@SQWDnDBase.get_proj_and_pbin(w)

Reverse engineer the projection and binning of an sqw object

>> [proj, pbin] = get_proj_and_pbin (w)

w sqw object

proj Projection as a projaxes object pbin Cell array, a row length 4, of the binning description of the

sqw object

sqw.@SQWDnDBase.log(w1)

Implements log(w1) for objects

>> w = log(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.log10(w1)

Implements log10(w1) for objects

>> w = log10(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.mask(win, mask_array)

Remove the bins indicated by the mask array

>> wout = mask (win, mask_array)

win Input sqw object

mask_array Array of 1 or 0 (or true or false) that indicate
which points to retain (true to retain, false to ignore)

Numeric or logical array of same number of elements

as the data_.

Note: mask will be applied to the stored data array

according as the projection axes, not the display axes. Thus permuting the display axes does not alter the effect of masking the data_.

wout Output dataset.

Original author: T.G.Perring

sqw.@SQWDnDBase.mask_pixels(win, mask_array)

Remove the pixels indicated by the mask array

>> wout = mask_pixels (win, mask_array) % Mask array >> wout = mask_pixels (win, wmask) % Mask according to pixel array

% contents

win Input sqw object

mask_array Array of 1 or 0 (or true or false) that indicate
which pixels to retain (true to retain, false to ignore)

Numeric or logical array of same number of pixels

as the data.

Note: mask will be applied to the stored data array

according as the projection axes, not the display axes. Thus permuting the display axes does not alter the effect of masking the data.

OR

wmask sqw object in which the signal in individual pixels is
interpreted as a mask array:

=1 (or true) to retain =0 (or false) to remove

wmask must have the same dimensionality, number of bins

along each dimension, and number of pixels in each bins as the array to be masked.

wout Output dataset.

sqw.@SQWDnDBase.mask_points(win, varargin)

Determine the points to keep on the basis of ranges and mask array.

>> sel = mask_points (win, ‘keep’, xkeep, ‘remove’, xremove, ‘mask’, mask)

or any selection (in any order) of the keyword-argument pairs e.g.

>> sel = mask_points (win, ‘mask’, mask, ‘remove’, xremove)

win Input sqw object

xkeep Ranges of display axes to retain for fitting. A range is specified by an array

of numbers which define a hypercube. For example in case of two dimensions:

[xlo, xhi, ylo, yhi]

or in the case of n-dimensions:

[x1_lo, x1_hi, x2_lo, x2_hi,…, xn_lo, xn_hi]

e.g. 1D: [50,70]

2D: [1,2,130,160]

More than one range can be defined in rows,

[Range_1; Range_2; Range_3;…; Range_m]

where each of the ranges are given in the format above.

xremove Ranges of display axes to remove from fitting.

mask Mask array of same number of elements as data array: 1 to keep, 0 to remove

Note: mask will be applied to the stored data array

according as the projection axes, not the display axes. Thus permuting the display axes does not alter the effect of masking the data_. The mask array works consistently with the input required by the mask method.

sel Mask array of same shape as data_. true for bins to keep, false to discard.

Advanced use: in addition the following two arguments, if present, suppress failure or the display of informational messges. Instead, the messages are returned to be used as desired.

ok =true if worked, =false if error

mess messages: if ok=true then informational or warning, if ok=false then the error message

sqw.@SQWDnDBase.mask_random_fraction_pixels(win, npix)

reduce the number of pixels randomly in a sqw object

The function uses the mask_pixels() function to keep only a fixed fraction of pixels randomly chosen. Useful when doing numerically intensive simulations of sqw objects with sqw_eval or fit_sqw, to speed things up

wout = mask_random_fraction_pixels(win,npix_frac)

win Input sqw object

npix_frac Fraction of pixels in win.data_.pix array to keep.

These are chosen at random. If win is an array then npix can either be a scalar, in which case all outputs will have the same number of retained pixels, or an array of the same size as win, in which case each mask is applied separately.

wout Output dataset.

sqw.@SQWDnDBase.mask_random_pixels(win, npix)

reduce the number of pixels randomly in a sqw object

The function uses the mask_pixels() function to keep only a fixed amount of pixels randomly chosen. Useful when doing numerically intensive simulations of sqw objects with sqw_eval or fit_sqw, to speed things up

wout = mask_random_pixels(win,npix)

win Input sqw object

npix Number of pixels in win.data_.pix array to keep.

These are chosen at random. If win is an array then npix can either be a scalar, in which case all outputs will have the same number of retained pixels, or an array of the same size as win, in which case each mask is applied separately.

wout Output dataset.

sqw.@SQWDnDBase.minus(w1, w2)

Implements w1 - w2 for objects

>> w = w1 - w2

w1, w2 Objects on which the binary operation is to be performed.

One of these can be a Matlab double (i.e. double precision) array, in which case the variance array is taken to be zero.

If w1, w2 are scalar objects with the same signal array sizes: - The operation is performed element-by-element.

If one of w1 or w2 is a double array (and the other is a scalar object): - If a scalar, apply to each element of the object signal. - If it is an array of the same size as the object signal

array, apply the operation element by element.

If one or both of w1 and w2 are arrays of objects: - If objects have same array sizes, the binary operation is

applied object element-by-object element.

  • If one of the objects is scalar (i.e. only one object),

then it is applied by the binary operation to each object in the other array.

If one of w1, w2 is an array of objects and the other is a double array: - If the double is a scalar, it is applied to every object

in the array.

  • If the double is an array with the same size as the object

array, then each element is applied as a scalar to the corresponding object in the object array.

  • If the double is an array with larger size than the object

array, then the array is resolved into a stack of arrays, where the stack has the same size as the object array, and the each array in the stack is applied to the corresponding object in the object array. [Note that for this operation to be valid, each object must have the same signal array size.]

w Output object or array of objects.

sqw.@SQWDnDBase.mldivide(w1, w2)

Implements w1 w2 for objects

>> w = w1 w2

w1, w2 Objects on which the binary operation is to be performed.

One of these can be a Matlab double (i.e. double precision) array, in which case the variance array is taken to be zero.

If w1, w2 are scalar objects with the same signal array sizes: - The operation is performed element-by-element.

If one of w1 or w2 is a double array (and the other is a scalar object): - If a scalar, apply to each element of the object signal. - If it is an array of the same size as the object signal

array, apply the operation element by element.

If one or both of w1 and w2 are arrays of objects: - If objects have same array sizes, the binary operation is

applied object element-by-object element.

  • If one of the objects is scalar (i.e. only one object),

then it is applied by the binary operation to each object in the other array.

If one of w1, w2 is an array of objects and the other is a double array: - If the double is a scalar, it is applied to every object

in the array.

  • If the double is an array with the same size as the object

array, then each element is applied as a scalar to the corresponding object in the object array.

  • If the double is an array with larger size than the object

array, then the array is resolved into a stack of arrays, where the stack has the same size as the object array, and the each array in the stack is applied to the corresponding object in the object array. [Note that for this operation to be valid, each object must have the same signal array size.]

w Output object or array of objects.

sqw.@SQWDnDBase.mpower(w1, w2)

Implements w1 ^ w2 for objects

>> w = w1 ^ w2

w1, w2 Objects on which the binary operation is to be performed.

One of these can be a Matlab double (i.e. double precision) array, in which case the variance array is taken to be zero.

If w1, w2 are scalar objects with the same signal array sizes: - The operation is performed element-by-element.

If one of w1 or w2 is a double array (and the other is a scalar object): - If a scalar, apply to each element of the object signal. - If it is an array of the same size as the object signal

array, apply the operation element by element.

If one or both of w1 and w2 are arrays of objects: - If objects have same array sizes, the binary operation is

applied object element-by-object element.

  • If one of the objects is scalar (i.e. only one object),

then it is applied by the binary operation to each object in the other array.

If one of w1, w2 is an array of objects and the other is a double array: - If the double is a scalar, it is applied to every object

in the array.

  • If the double is an array with the same size as the object

array, then each element is applied as a scalar to the corresponding object in the object array.

  • If the double is an array with larger size than the object

array, then the array is resolved into a stack of arrays, where the stack has the same size as the object array, and the each array in the stack is applied to the corresponding object in the object array. [Note that for this operation to be valid, each object must have the same signal array size.]

w Output object or array of objects.

sqw.@SQWDnDBase.mrdivide(w1, w2)

Implements w1 / w2 for objects

>> w = w1 / w2

w1, w2 Objects on which the binary operation is to be performed.

One of these can be a Matlab double (i.e. double precision) array, in which case the variance array is taken to be zero.

If w1, w2 are scalar objects with the same signal array sizes: - The operation is performed element-by-element.

If one of w1 or w2 is a double array (and the other is a scalar object): - If a scalar, apply to each element of the object signal. - If it is an array of the same size as the object signal

array, apply the operation element by element.

If one or both of w1 and w2 are arrays of objects: - If objects have same array sizes, the binary operation is

applied object element-by-object element.

  • If one of the objects is scalar (i.e. only one object),

then it is applied by the binary operation to each object in the other array.

If one of w1, w2 is an array of objects and the other is a double array: - If the double is a scalar, it is applied to every object

in the array.

  • If the double is an array with the same size as the object

array, then each element is applied as a scalar to the corresponding object in the object array.

  • If the double is an array with larger size than the object

array, then the array is resolved into a stack of arrays, where the stack has the same size as the object array, and the each array in the stack is applied to the corresponding object in the object array. [Note that for this operation to be valid, each object must have the same signal array size.]

w Output object or array of objects.

sqw.@SQWDnDBase.mtimes(w1, w2)

Implements w1 * w2 for objects

>> w = w1 * w2

w1, w2 Objects on which the binary operation is to be performed.

One of these can be a Matlab double (i.e. double precision) array, in which case the variance array is taken to be zero.

If w1, w2 are scalar objects with the same signal array sizes: - The operation is performed element-by-element.

If one of w1 or w2 is a double array (and the other is a scalar object): - If a scalar, apply to each element of the object signal. - If it is an array of the same size as the object signal

array, apply the operation element by element.

If one or both of w1 and w2 are arrays of objects: - If objects have same array sizes, the binary operation is

applied object element-by-object element.

  • If one of the objects is scalar (i.e. only one object),

then it is applied by the binary operation to each object in the other array.

If one of w1, w2 is an array of objects and the other is a double array: - If the double is a scalar, it is applied to every object

in the array.

  • If the double is an array with the same size as the object

array, then each element is applied as a scalar to the corresponding object in the object array.

  • If the double is an array with larger size than the object

array, then the array is resolved into a stack of arrays, where the stack has the same size as the object array, and the each array in the stack is applied to the corresponding object in the object array. [Note that for this operation to be valid, each object must have the same signal array size.]

w Output object or array of objects.

sqw.@SQWDnDBase.multifit_func(varargin)

Simultaneously fit function(s) to one or more sqw objects

>> myobj = multifit_func (w1, w2, …) % w1, w2 objects or arrays of objects

This creates a fitting object of class mfclass_Horace with the provided data, which can then be manipulated to add further data, set the fitting functions, initial parameter values etc. and fit or simulate the data. For details about how to do this <a href=”matlab:help(‘mfclass_Horace’);”>Click here</a>

For example:

>> myobj = multifit_func (w1, w2, …); % set the data

:

>> myobj = myobj.set_fun (@function_name, pars); % set forgraound function(s) >> myobj = myobj.set_bfun (@function_name, pars); % set background function(s)

:

>> myobj = myobj.set_free (pfree); % set which parameters are floating >> myobj = myobj.set_bfree (bpfree); % set which parameters are floating >> [wfit,fitpars] = myobj.fit; % perform fit

This method fits function(s) of the plot axes for both the foreground and the background function(s). The format of the fit functions depends on the number of plot axes for each sqw object. For examples see: <a href=”matlab:edit(‘example_1d_function’);”>example_1d_function</a> <a href=”matlab:edit(‘example_2d_function’);”>example_2d_function</a> <a href=”matlab:edit(‘example_3d_function’);”>example_3d_function</a>

See also multifit_sqw multifit_sqw_sqw

[Help for legacy use (2017 and earlier):

If you are still using the legacy version then it is strongly recommended that you change to the new operation. Help for the legacy operation can be <a href=”matlab:help(‘sqw/multifit_legacy_func’);”>found here</a>]

sqw.@SQWDnDBase.multifit_sqw(varargin)

Simultaneously fit function(s) of S(Q,w)to one or more sqw objects

>> myobj = multifit_sqw (w1, w2, …) % w1, w2 objects or arrays of objects

This creates a fitting object of class mfclass_Horace_sqw with the provided data, which can then be manipulated to add further data, set the fitting functions, initial parameter values etc. and fit or simulate the data. For details about how to do this <a href=”matlab:help(‘mfclass_Horace_sqw’);”>Click here</a>

For example:

>> myobj = multifit_sqw (w1, w2, …); % set the data

:

>> myobj = myobj.set_fun (@function_name, pars); % set forgraound function(s) >> myobj = myobj.set_bfun (@function_name, pars); % set background function(s)

:

>> myobj = myobj.set_free (pfree); % set which parameters are floating >> myobj = myobj.set_bfree (bpfree); % set which parameters are floating >> [wfit,fitpars] = myobj.fit; % perform fit

This method fits model(s) for S(Q,w) as the foreground function(s), and function(s) of the plot axes for the background function(s)

For the format of foreground fit functions: <a href=”matlab:edit(‘example_sqw_spin_waves’);”>Damped spin waves</a> <a href=”matlab:edit(‘example_sqw_flat_mode’);”>Dispersionless excitations</a>

The format of the background fit functions depends on the number of plot axes for each sqw object. For examples see: <a href=”matlab:edit(‘example_1d_function’);”>example_1d_function</a> <a href=”matlab:edit(‘example_2d_function’);”>example_2d_function</a> <a href=”matlab:edit(‘example_3d_function’);”>example_3d_function</a>

See also multifit multifit_sqw_sqw

[Help for legacy use (2017 and earlier):

If you are still using the legacy version then it is strongly recommended that you change to the new operation. Help for the legacy operation can be <a href=”matlab:help(‘sqw/multifit_legacy_sqw’);”>found here</a>]

sqw.@SQWDnDBase.multifit_sqw_sqw(varargin)

Simultaneously fit function(s) of S(Q,w)to one or more sqw objects

>> myobj = multifit_sqw_sqw (w1, w2, …) % w1, w2 objects or arrays of objects

This creates a fitting object of class mfclass_Horace_sqw_sqw with the provided data, which can then be manipulated to add further data, set the fitting functions, initial parameter values etc. and fit or simulate the data. For details about how to do this <a href=”matlab:help(‘mfclass_Horace_sqw_sqw’);”>Click here</a>

For example:

>> myobj = multifit_sqw_sqw (w1, w2, …); % set the data

:

>> myobj = myobj.set_fun (@function_name, pars); % set forgraound function(s) >> myobj = myobj.set_bfun (@function_name, pars); % set background function(s)

:

>> myobj = myobj.set_free (pfree); % set which parameters are floating >> myobj = myobj.set_bfree (bpfree); % set which parameters are floating >> [wfit,fitpars] = myobj.fit; % perform fit

This method fits function(s) of S(Q,w) as both the foreground and the background function(s). For the format of the fit functions: <a href=”matlab:edit(‘example_sqw_spin_waves’);”>Damped spin waves</a> <a href=”matlab:edit(‘example_sqw_flat_mode’);”>Dispersionless excitations</a>

See also multifit multifit_sqw

[Help for legacy use (2017 and earlier):

If you are still using the legacy version then it is strongly recommended that you change to the new operation. Help for the legacy operation can be <a href=”matlab:help(‘sqw/multifit_legacy_sqw_sqw’);”>found here</a>]

sqw.@SQWDnDBase.pa(w, varargin)

Overplot an area plot of a 2D sqw dataset or array of datasets

>> pa(w)

Advanced use:
>> pa(w,’name’,fig_name) % overplot on figure with name = fig_name

% or figure with given figure number or handle

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = pa(w,…)

sqw.@SQWDnDBase.paoc(w)

Overplot an area plot of a 2D sqw dataset or array of datasets on the current figure

>> paoc(w)

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = paoc(w)

sqw.@SQWDnDBase.pd(w, varargin)

Overplot markers, error bars and lines for a 1D sqw object or array of objects on an existing plot

>> pd(w)

Advanced use:
>> pd(w,’name’,fig_name) % overplot on figure with name = fig_name

% or figure with given figure number or handle

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = pd(w,…)

sqw.@SQWDnDBase.pdoc(w)

Overplot markers, error bars and lines for a 1D sqw object or array of objects on the current plot

>> pdoc(w)

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = pdoc(w)

sqw.@SQWDnDBase.pe(w, varargin)

Overplot error bars for a 1D sqw object or array of objects on an existing plot

>> pe(w)

Advanced use:
>> pe(w,’name’,fig_name) % overplot on figure with name = fig_name

% or figure with given figure number or handle

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = pe(w,…)

sqw.@SQWDnDBase.peoc(w)

Overplot error bars for a 1D sqw object or array of objects on the current plot

>> peoc(w)

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = peoc(w)

sqw.@SQWDnDBase.ph(w, varargin)

Overplot histogram for a 1D sqw object or array of objects on an existing plot

>> ph(w)

Advanced use:
>> ph(w,’name’,fig_name) % overplot on figure with name = fig_name

% or figure with given figure number or handle

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = ph(w,…)

sqw.@SQWDnDBase.phoc(w)

Overplot histogram for a 1D sqw object or array of objects on the current plot

>> phoc(w)

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = phoc(w)

sqw.@SQWDnDBase.pl(w, varargin)

Overplot line for a 1D sqw object or array of objects on an existing plot

>> pl(w)

Advanced use:
>> pl(w,’name’,fig_name) % overplot on figure with name = fig_name

% or figure with given figure number or handle

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = pl(w,…)

sqw.@SQWDnDBase.ploc(w)

Overplot line for a 1D sqw object or array of objects on the current plot

>> ploc(w)

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = ploc(w)

sqw.@SQWDnDBase.plot(w, varargin)

Plot 1D, 2D or 3D sqw object or array of objects

>> plot(w) >> plot(w,opt1,opt2,…) % plot with optional arguments

Equivalent to:

>> dp(w) % 1D dataset >> dp(w,…)

>> da(w) % 2D dataset >> da(w,…)

>> sliceomatic(w) % 3D dataset >> sliceomatic(w,…)

For details of optional parameters type >> help sqw/dp, >> help sqw/da, or >> help sqw/sliceomatic as appropriate

sqw.@SQWDnDBase.plotover(w, varargin)

Overplot 1D, 2D or 3D sqw object or array of objects

>> plotover(w) >> plotover(w,opt1,opt2,…) % plot with optional arguments

Equivalent to:

>> pp(w) % 1D dataset >> pp(w,…)

>> pa(w) % 2D dataset >> pa(w,…)

For details of optional parameters type >> help sqw/pp, >> help sqw/pa, as appropriate

sqw.@SQWDnDBase.plus(w1, w2)

Implements w1 + w2 for objects

>> w = w1 + w2

w1, w2 Objects on which the binary operation is to be performed.

One of these can be a Matlab double (i.e. double precision) array, in which case the variance array is taken to be zero.

If w1, w2 are scalar objects with the same signal array sizes: - The operation is performed element-by-element.

If one of w1 or w2 is a double array (and the other is a scalar object): - If a scalar, apply to each element of the object signal. - If it is an array of the same size as the object signal

array, apply the operation element by element.

If one or both of w1 and w2 are arrays of objects: - If objects have same array sizes, the binary operation is

applied object element-by-object element.

  • If one of the objects is scalar (i.e. only one object),

then it is applied by the binary operation to each object in the other array.

If one of w1, w2 is an array of objects and the other is a double array: - If the double is a scalar, it is applied to every object

in the array.

  • If the double is an array with the same size as the object

array, then each element is applied as a scalar to the corresponding object in the object array.

  • If the double is an array with larger size than the object

array, then the array is resolved into a stack of arrays, where the stack has the same size as the object array, and the each array in the stack is applied to the corresponding object in the object array. [Note that for this operation to be valid, each object must have the same signal array size.]

w Output object or array of objects.

sqw.@SQWDnDBase.pm(w, varargin)

Overplot markers for a 1D sqw object or array of objects on an existing plot

>> pm(w)

Advanced use:
>> pm(w,’name’,fig_name) % overplot on figure with name = fig_name

% or figure with given figure number or handle

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = pm(w,…)

sqw.@SQWDnDBase.pmoc(w)

Overplot markers for a 1D sqw object or array of objects on the current plot

>> pmoc(w)

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = pmoc(w)

sqw.@SQWDnDBase.pp(w, varargin)

Overplot markers and error bars for a 1D sqw object or array of objects on an existing plot

>> pp(w)

Advanced use:
>> pp(w,’name’,fig_name) % overplot on figure with name = fig_name

% or figure with given figure number or handle

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = pp(w,…)

sqw.@SQWDnDBase.ppoc(w)

Overplot markers and error bars for a 1D sqw object or array of objects on the current plot

>> ppoc(w)

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = ppoc(w)

sqw.@SQWDnDBase.ps(w, varargin)

Overplot a surface plot of a 2D sqw dataset or array of datasets

>> ps(w)

Advanced use:
>> ps(w,’name’,fig_name) % overplot on figure with name = fig_name

% or figure with given figure number or handle

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = ps(w,…)

sqw.@SQWDnDBase.ps2(w, varargin)

Overplot a surface plot of a 2D sqw dataset or array of datasets

>> ps2(w) % Use error bars to set colour scale >> ps2(w,wc) % Signal in wc sets colour scale

% wc can be any object with a signal array with same % size as w, e.g. sqw object, IX_dataset_2d object, or % a numeric array. % - If w is an array of objects, then wc must contain % the same number of objects. % - If wc is a numeric array then w must be a scalar % object.

Differs from ds in that the signal sets the z axis, and the colouring is set by the error bars, or by another object. This enables two related functions to be plotted (e.g. dispersion relation where the ‘signal’ array holds the energy and the error array holds the spectral weight).

Advanced use:
>> ps(w,’name’,fig_name) % overplot on figure with name = fig_name

% or figure with given figure number or handle

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = ps2(w,…)

sqw.@SQWDnDBase.ps2oc(w, varargin)

Overplot a surface plot of a 2D sqw dataset or array of datasets on the current figure The colour scale comes from a second source

>> ps2oc(w) % Use error bars to set colour scale >> ps2oc(w,wc) % Signal in wc sets colour scale

% wc can be any object with a signal array with same % size as w, e.g. sqw object, IX_dataset_2d object, or % a numeric array. % - If w is an array of objects, then wc must contain % the same number of objects. % - If wc is a numeric array then w must be a scalar % object.

Differs from ds in that the signal sets the z axis, and the colouring is set by the error bars, or by another object. This enables two related functions to be plotted (e.g. dispersion relation where the ‘signal’ array holds the energy and the error array holds the spectral weight).

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = ps2oc(w)

sqw.@SQWDnDBase.psoc(w)

Overplot a surface plot of a 2D sqw dataset or array of datasets on the current figure

>> psoc(w)

Return figure, axes and plot handles:

>> [fig_handle, axes_handle, plot_handle] = psoc(w)

sqw.@SQWDnDBase.recompute_bin_data(w)

Given sqw_type object, recompute w.data.s and w.data.e from the contents of pix array

>> wout=recompute_bin_data(w)

sqw.@SQWDnDBase.save(w, varargin)

Save a sqw object or array of sqw objects to file

>> save (w) % prompt for file >> save (w, file) % give file >> save (w, file,loader) % save file using specific data loader

(-update option, is provided, will be ignored)

>> save (w, file,’-update’) % if the target file exist, update it to

latest format if this is possible. If update is possible, pixels in file will not be overwritten.

Input:

w sqw object file [optional] File for output. if none given, then prompted for a file

Note that if w is an array of sqw objects then file must be a cell array of filenames of the same size.

Output:

sqw.@SQWDnDBase.sec(w1)

Implements sec(w1) for objects

>> w = sec(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.sech(w1)

Implements sech(w1) for objects

>> w = sech(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.sigvar_get(w)

Get signal and variance from object, and a logical array of which values to ignore

>> [s,var,mask_null] = sigvar_get (w)

sqw.@SQWDnDBase.sin(w1)

Implements sin(w1) for objects

>> w = sin(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.sinh(w1)

Implements sinh(w1) for objects

>> w = sinh(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.sliceomatic(w, varargin)

Plots 3D sqw object using sliceomatic

>> sliceomatic (w) >> sliceomatic (w, ‘isonormals’, true) % to enable isonormals

>> sliceomatic (w,…,’-noaspect’) % Do not change aspect ratio

% according to data axes unit lengths

To get handles to the graphics figure:

>> [figureHandle_, axesHandle_, plotHandle_] = sliceomatic(w,…)

NOTES:

  • Ensure that the slice color plotting is in ‘texture’ mode -

    On the ‘AllSlices’ menu click ‘Color Texture’. No indication will be made on this menu to show that it has been selected, but you can see the result if you right-click on an arrow indicating a slice on the graphics window.

  • To set the default for future Sliceomatic sessions -

    On the ‘Object_Defaults’ menu select ‘Slice Color Texture’

sqw.@SQWDnDBase.sliceomatic_overview(w, varargin)

Plots 3D sqw object using sliceomatic with view straight down one of the axes

>> sliceomatic_overview (w) % down third (vertical) axis >> sliceomatic_overview (w, axis) % down axis of choice (axis=1,2 or 3)

>> sliceomatic_overview (w,… ‘isonormals’, true) % to enable isonormals

>> sliceomatic_overview (w,…,’-noaspect’) % Do not change aspect ratio

% according to data axes unit lengths

To get handles to the graphics figure:

>> [figureHandle_, axesHandle_, plotHandle_] = sliceomatic(w,…)

Do a sliceomatic plot, but set the axes so that we look straight down the chosen axis, so that when the slider is moved we get a series of what appear to be 2d slices.

NOTES:

  • Ensure that the slice color plotting is in ‘texture’ mode -

    On the ‘AllSlices’ menu click ‘Color Texture’. No indication will be made on this menu to show that it has been selected, but you can see the result if you right-click on an arrow indicating a slice on the graphics window.

  • To set the default for future Sliceomatic sessions -

    On the ‘Object_Defaults’ menu select ‘Slice Color Texture’

sqw.@SQWDnDBase.smooth(win, varargin)

Smooth method - gataway to dnd object smoothing only.

sqw.@SQWDnDBase.smooth_units(win, varargin)

Smooth method - gataway to dnd object smoothing only.

sqw.@SQWDnDBase.sqrt(w1)

Implements sqrt(w1) for objects

>> w = sqrt(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.sqw_eval(win, sqwfunc, pars, varargin)

Calculate sqw for a model scattering function

>> wout = sqw_eval(win, sqwfunc, p) >> wout = sqw_eval(___, ‘-all’) >> wout = sqw_eval(___, ‘all’, true) >> wout = sqw_eval(___, ‘-average’) >> wout = sqw_eval(___, ‘average’, true) >> sqw_eval(___, ‘outfile’, outfile) >> wout = sqw_eval(___, ‘outfile’, outfile) >> sqw_eval(__, ‘outfile’, outfile, ‘filebacked’, true) >> wout = sqw_eval(__, ‘filebacked’, true)

win Dataset (or array of datasets) that provides the axes and points

for the calculation

sqwfunc Handle to function that calculates S(Q, w)
Most commonly used form is:

weight = sqwfunc (qh, qk, ql, en, p)

where

qh,qk,ql,en Arrays containing the coordinates of a set of points p Vector of parameters needed by dispersion function

e.g. [A, js, gam] as intensity, exchange, lifetime

weight Array containing calculated spectral weight

More general form is:

weight = sqwfunc (qh, qk, ql, en, p, c1, c2, ..)

where
p Typically a vector of parameters that we might want

to fit in a least-squares algorithm

c1, c2, … Other constant parameters e.g. file name for look-up

table

pars Arguments needed by the function. Most commonly, a vector of parameter

values e.g. [A, js, gam] as intensity, exchange, lifetime. If a more general set of parameters is required by the function, then package these into a cell array and pass that as pars. In the example above then pars = {p, c1, c2, …}

outfile If present, the outputs will be written to the file of the given

name/path. If numel(win) > 1, outfile must either be omitted or be a cell array of file paths with equal number of elements as win.

all If true, requests that the calculated sqw be returned over

the whole of the domain of the input dataset. If false, then the function will be returned only at those points of the dataset that contain data_.

Applies only to input with no pixel information - it is ignored if

full sqw object. [default = false]

average If true, requests that the calculated sqw be computed for the

average values of h, k, l of the pixels in a bin, not for each pixel individually. Reduces cost of expensive calculations. Applies only to the case of sqw object with pixel information - it is ignored if dnd type object. [default = false]

filebacked If true, the result of the function will be saved to file and

the output will be a file path. If no outfile is specified, a unique path within tempdir() will be generated. Default is false.

Note: all optional string input parameters can be truncated up to minimal

difference between them e.g. routine would accept ‘al’ and ‘av’, ‘ave’, ‘aver’ etc….

wout If filebacked is false, an sqw object or array of sqw objects.

If filebacked is true, a file path or cell array of file paths. Output argument must be specified if outfile not given.

sqw.@SQWDnDBase.sqw_eval_nopix_(obj, sqwfunc, all_bins, pars)

SQW_EVAL_NOPIX_

Helper function for sqw eval executed on a pixel-less object (i.e. DnD or SQW with no pixels Called by sqw_eval_ defined in sqw/DnDBase

obj Dataset (or array of datasets) that provides the axes and points

for the calculation

sqwfunc Handle to function that calculates S(Q,w) all_bins Boolean flag wither to apply function to all bins or only those contaiing data pars Arguments needed by the function.

sqw.@SQWDnDBase.tan(w1)

Implements tan(w1) for objects

>> w = tan(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.tanh(w1)

Implements tanh(w1) for objects

>> w = tanh(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.uminus(w1)

Implements uminus(w1) for objects

>> w = uminus(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.uplus(w1)

Implements uplus(w1) for objects

>> w = uplus(w1)

w1 Input object or array of objects on which to apply the

unary operator.

w Output object or array of objects.

sqw.@SQWDnDBase.xye(w, null_value)

Get the bin centres, intensity and error bar for a 1D, 2D, 3D or 4D dataset

>> S = xye(w) >> S = xye(w, null_value)

w sqw object or array of objects (1D, 2D, 3D or 4D) null_value Numeric value to substitute for the intensity in bins

with no data. Default: NaN

S Structure with the following fields:

x If a 1D sqw object, a column vector with the bin centres

Otherwise, a cell array (row vector) of arrays, each array containing the bin centres for the plot axes. Each array of bin centres is a column vector.

y n1 x n2 x… array of intensities, where n1 is the number of

bins along the first plot axis, n2 the number of bins along the second plot axis etc

e n1 x n2 x… array of error bars, where n1 is the number of

bins along the first plot axis, n2 the number of bins along the second plot axis etc

Experiment class

class sqw.@Experiment.Experiment(varargin)

EXPERIMENT Container object for all data describing the Experiment

Experiment(varargin)

Create a new Experiment object.

obj = Experiment() obj = Experiment(detector_array[s], instrument[s], sample[s])

Required:

detector_array Detector array (IX_detector_array objects) instrument Instrument (Concrete class inheriting IX_inst) sample Sample data (IX_sample object)

Each argument can be a single object or array of objects.

classVersion(~)

define version of the class to store in mat-files and nxsqw data format. Each new version would presumably read the older version, so version substitution is based on this number

static combine_experiments(exp_cellarray, allow_equal_headers, drop_subzone_headers)
take cellarray of experiments (e.g., generated from each runfile build

during gen_sqw generation) and combine then together into single Experiment info class

This is the HACK, providing only basic functionality. Previous header-s on the basis of sqw_header and part, present in write_nsqw_to_sqw implementation offers much more.

TODO: Do proper optimization on the way. See

sqw_header.header_combine(header,allow_equal_headers,drop_subzone_headers)

TODO: use allow_equal_headers,drop_subzone_headers variables

appropriately

TODO: repeat at least the logic within sqw_header helper class

and write_nsqw_to_sqw combine/check headers operation

convert_to_old_headers(header_num)

convert Experiment into the structure suitable to be stored in old binary sqw files (up to version 3.xxx)

this structure is also used in number of places of the old code where, e.g., structure sorting is implemented but this usage is deprecated and will be removed in a future.

get_aver_experiment()

some, presumably average, run-data. Naive implementation, all data are the same

get_unique_instruments()

compatibility fields with old binary file formats TODO: needs proper implementation

get_unique_samples()

compatibility fields with old binary file formats TODO: needs proper implementation

header_average()

very crude implementation for the header, average over all runs.

indepFields(~)

get independent fields, which fully define the state of the serializable object.

instruments = None

Mirrors of private properties

is_same_ebins()

return true if all energy bins of all runs are the same

static loadobj(S)

boilerplate loadobj method, calling generic method of save-able class

n_runs = None

return the number of runs, this class contains

set_efix_emode(efix, emode)

change efix and (optionally) emode in all experiment descriptions if emode is absent or described by any character string, the emode is kept unchanged

PixelData class

class sqw.PixelData.@PixelData.PixelData(arg, mem_alloc)

PixelData Provides an interface for access to pixel data

This class provides getters and setters for each data column in an SQW pixel array. You can access the data using the attributes listed below, using the get_data() method (to retrieve column data) or using the get_pixels() method (retrieve row data).

Construct this class with an 9 x N array, a file path to an SQW object or an instance of sqw_binfile_common.

>> pix_data = PixelData(data); >> pix_data = PixelData(‘/path/to/sqw.sqw’); >> pix_data = PixelData(‘/path/to/sqw.sqw’, mem_alloc); >> pix_data = PixelData(faccess_obj); >> pix_data = PixelData(faccess_obj, mem_alloc);

Constructing via a file or sqw_binfile_common will create a file-backed data object. No pixel data will be loaded from the file on construction. Data will be loaded when a getter is called e.g. pix_data.signal. Data will be loaded in pages such that the data held in memory will not exceed the size (in bytes) specified by private attribute page_memory_size_ - this can be set on construction (see mem_alloc above).

The file-backed operations work by loading “pages” of data into memory as required. If editing pixels, to avoid losing changes, if a page has been edited and the next page is then loaded, the “dirty” page will be written to a tmp file. This class’s getters will then retrieve data from the tmp file if that data is requested from the “dirty” page. Note that “dirty” pages are written to tmp files as floats, but stored in memory as double. This means data is truncated when moving pages, hence pixel data should not be relied upon being accurate to double precision.

Usage:

>> pix_data = PixelData(data) >> signal = pix_data.signal;

or equivalently:

>> pix_data = PixelData(); >> pix_data.data = data; >> signal = pix_data.get_data(‘signal’);

To retrieve multiple fields of data, e.g. run_idx and energy_idx, for pixels 1 to 10:

>> pix_data = PixelData(data); >> signal = pix_data.get_data({‘run_idx’, ‘energy_idx’}, 1:10);

To retrieve data for pixels 1, 4 and 10 (returning another PixelData object):

>> pix_data = PixelData(data); >> pixel_subset = pix_data.get_pixels([1, 4, 10])

To sum the signal of a file-backed object where the page size is less than amount of data in the file:

>> pix = PixelData(‘my_data.sqw’) >> signal_sum = 0; >> while pix.has_more() >> signal_sum = signal_sum + pix.signal; >> pix.advance(); >> end

Properties:
u1, u2, u3 - The 1st, 2nd and 3rd dimensions of the Crystal

Cartesian coordinates in projection axes, units are per Angstrom (1 x n arrays)

dE - The energy transfer value for each pixel in meV (1 x n array) coordinates - The coords in projection axes of the pixel data [u1, u2, u3, dE] (4 x n array) q_coordinates - The spacial coords in projection axes of the pixel data [u1, u2, u3] (3 x n array) run_idx - The run index the pixel originated from (1 x n array) detector_idx - The detector group number in the detector listing for the pixels (1 x n array) energy_idx - The energy bin numbers (1 x n array) signal - The signal array (1 x n array). variance - The variance on the signal array (variance i.e. error bar squared) (1 x n array)

num_pixels - The number of pixels in the data block. pix_range - [2x4] array of the range of pixels coordinates in Crystal Cartesian coordinate system.

data - The raw pixel data - usage of this attribute is discouraged, the structure

of the return value is not guaranteed.

page_size - The number of pixels in the currently loaded page.

DATA_POINT_SIZE = '8'

num bytes in a double

DEFAULT_PAGE_SIZE = 'realmax'

this gives no paging by default

PixelData(arg, mem_alloc)

Construct a PixelData object from the given data. Default construction initialises the underlying data as an empty (9 x 0) array.

>> obj = PixelData(ones(9, 200))

>> obj = PixelData(200) % initialise 200 pixels with underlying data set to zero

>> obj = PixelData(file_path) % initialise pixel data from an sqw file

>> obj = PixelData(faccess_reader) % initialise pixel data from an sqw file reader

>> obj = PixelData(faccess_reader, mem_alloc) % set maximum memory allocation

arg A 9 x n matrix, where each row corresponds to a pixel and
the columns correspond to the following:

col 1: u1 col 2: u2 col 3: u3 col 4: dE col 5: run_idx col 6: detector_idx col 7: energy_idx col 8: signal col 9: variance

arg An integer specifying the desired number of pixels. The underlying

data will be filled with zeros.

arg A path to an SQW file.

arg An instance of an sqw_binfile_common file reader.

mem_alloc The maximum amount of memory allocated to hold pixel

data in bytes. If pixels cannot all be held in memory at one time, they will be loaded from the file (specified by ‘arg’) when they are required. This argument does nothing if the class is constructed with in-memory data. (Optional)

advance(varargin)

Load the next page of pixel data from the file backing the object

This function will throw a PIXELDATA:advance error if attempting to advance past the final page of data in the file.

This function does nothing if the pixel data is not file-backed.

>> obj.advance() >> obj.advance(‘nosave’, true)

nosave Keyword argument. Set to true to discard changes to cache.

(default: false)

current_page_number The new page and total number of pages advance will walk through to complete the algorithm

base_page_size = None

The number of pixels that can fit in one page of data

static cat(varargin)

Concatenate the given PixelData objects’ pixels. This function performs a straight-forward data concatenation.

>> joined_pix = PixelData.cat(pix_data1, pix_data2);

varargin A cell array of PixelData objects

obj A PixelData object containing all the pixels in the inputted

PixelData objects

coordinates = None

The coordinates of the pixels in the projection axes, i.e.: u1,

copy()

Make an independent copy of this object This method simply constructs a new PixelData instance by calling the constructor with the input object as an argument. Because of this, any properties that need to be explicitly copied must be copied within this classes “copy-constructor”.

dE = None

The array of energy deltas of the pixels (1 x n array) [meV]

data = None

The full raw pixel data block. Usage of this attribute exposes

delete()

Class destructor to delete any temporary files

detector_idx = None

The detector group number in the detector listing for the pixels (1 x n array)

energy_idx = None

The energy bin numbers (1 x n array)

file_path = None

The file that the pixel data has been read from, empty if no file

has_more()

Returns true if there are subsequent pixels stored in the file that are not held in the current page

>> has_more = pix.has_more();

is_filebacked()

Return true if the pixel data is backed by a file or files. Returns false if all pixel data is held in memory

isempty()

Return true if the PixelData object holds no pixel data

static loadobj(S)

Load a PixelData object from a .mat file

>> obj = PixelData.loadobj(S)
S A data, produeced by saveobj operation and stored

in .mat file

obj An instance of PixelData object or array of objects

move_to_first_page()

Reset the object to point to the first page of pixel data in the file and clear the current cache

This function does nothing if pixels are not file-backed.

num_pixels = None

The number of pixels in the data block

page_size = None

The number of pixels in the current page

pix_range = None

The range of pixels coordinates in Crystal Cartesian

q_coordinates = None

The spatial dimensions of the Crystal Cartesian

run_idx = None

The run index the pixel originated from (1 x n array)

set_range(pix_range)

Function allows to set the pixels range (min/max values of pixels coordinates)

Use with caution!!! No checks that the set range is the correct range for pixels, holded by the class are performed, while subsequent algorithms may rely on pix range to be correct. A out-of memory write can occur during rebinning if the range is smaller, then the actual range.

Necessary to set up the pixel range when filebased pixels are modified by algorithm and correct range calculations are expensive

signal = None

The signal array (1 x n array)

struct()

convert object into saveable and serializable structure

u1 = None

The 1st dimension of the Crystal Cartesian orientation (1 x n array) [A^-1]

u2 = None

The 2nd dimension of the Crystal Cartesian orientation (1 x n array) [A^-1]

u3 = None

The 3rd dimension of the Crystal Cartesian orientation (1 x n array) [A^-1]

variance = None

The variance on the signal array

sqw.PixelData.@PixelData.append(obj, pix)

Join the pixels in the given PixelData object to the end of this object.

The pixels to append must all be in memory and you cannot append pixels if the inputted PixelData object has more pixels than are allowed in a single page.

pix A PixelData object containing the pixels to append

sqw.PixelData.@PixelData.compute_bin_data(obj, npix)

Compute the mean signal and variance given the number of contributing pixels for each bin Returns empty arrays if obj contains no pixels.

>> [mean_signal, mean_variance] = compute_bin_data(obj, npix)

npix The number of contributing pixels to each bin of the plot axes

mean_signal The average signal for each plot axis bin.

size(mean_signal) = size(npix)

mean_variance The average variance for each plot axis bin.

size(mean_variance) = size(npix)

sqw.PixelData.@PixelData.do_binary_op(obj, operand, binary_op, varargin)
DO_BINARY_OP perform a binary operation between this object and the given

operand

>> pix_diff = obj.do_binary_op(other_pix, @minus_single, ‘flip’, true)

>> pix_sum = obj.do_binary_op(signal_array, @plus_single, ‘npix’, npix)

>> pix_sum = obj.do_binary_op(signal_array, @plus_single, ‘npix’, npix, ‘flip’, true)

operand The second operand to use in the binary operation.
The operand must have one of the following types:
  • scalar double

  • double array, the size of the array must be equal to obj.num_pixels

  • object with fields ‘s’ and ‘e’ (e.g. dnd or sigvar)

  • another PixelData object with obj.num_pixels equal

binary_op Function handle pointing to the desired binary operation. The

function should take 2 objects with ‘.s’ and ‘.e’ attributes, e.g. a sigvar object

flip Flip the order of the operands, e.g. perform “this - operand” if

flip is false, perform “operand - this” if flip is true.

npix An array giving number of pixels in each bin. This argument should

have equal size to operand (assuming operand is numeric) and sum(npix, [], ‘all’) must be equal to obj.num_pixels

sqw.PixelData.@PixelData.do_unary_op(obj, unary_op)

Perform a unary operation on this object’s signal and variance arrays

unary_op Function handle pointing to the operation to perform. This

operation should take a sigvar object as an argument.

sqw.PixelData.@PixelData.equal_to_tol(obj, other_pix, varargin)

EQUAL_TO_TOL Check if two PixelData objects are equal to a given tolerance

pix The first pixel data object to compare.

other_pix The second pixel data object to compare.

tol Tolerance criterion for numeric arrays

(default = [0, 0] i.e. equality) It has the form: [abs_tol, rel_tol] where

abs_tol absolute tolerance (>=0; if =0 equality required) rel_tol relative tolerance (>=0; if =0 equality required)

If either criterion is satisfied then equality within tolerance is accepted.

Examples:

[1e-4, 1e-6] absolute 1e-4 or relative 1e-6 required [1e-4, 0] absolute 1e-4 required [0, 1e-6] relative 1e-6 required [0, 0] equality required 0 equivalent to [0,0]

A scalar tolerance can be given where the sign determines if the tolerance is absolute or relative:

+ve : absolute tolerance abserr = abs(a-b) -ve : relative tolerance relerr = abs(a-b)/max(abs(a),abs(b))

Examples:

1e-4 absolute tolerance, equivalent to [1e-4, 0] -1e-6 relative tolerance, equivalent to [0, 1e-6]

nan_equal Treat NaNs as equal (true or false; default=true).

name_a Explicit name of variable a for use in messages

Usually not required, as the name of a variable will be discovered. However, if the input argument is an array element e.g. my_variable{3} then the name is not discoverable in Matlab, and default ‘input_1’ will be used unless a different value is given with the keyword ‘name_a’. (default = ‘input_1’).

name_b Explicit name of variable b for use in messages.

The same comments apply as for ‘name_a’ except the default is ‘input_2’. (default = ‘input_2’).

sqw.PixelData.@PixelData.get_data(obj, pix_fields, varargin)

GET_DATA Retrieve data for a field, or fields, for the given pixel indices in the full pixel block. If no pixel indices are given, the full range of pixels is returned.

This method provides a convenient way of retrieving multiple fields of data from the pixel block. When retrieving multiple fields, the columns of data will be ordered corresponding to the order the fields appear in the inputted cell array.

>> sig_and_err = pix.get_data({‘signal’, ‘variance’})

retrieves the signal and variance over the whole range of pixels

>> run_det_id_range = pix.get_data({‘run_idx’, ‘detector_idx’}, 15640:19244);

retrieves the run and detector IDs for pixels 15640 to 19244

pix_fields The name of a field, or a cell array of field names abs_pix_indices The pixel indices to retrieve, if not given, get full range.

The syntax for these indices attempts to replicate indexing into a regular Matlab array. You can use logical indices as well as normal indices, and you can index into the array “out-of-order”. However, you cannot use end, but it is possible to achieve the same effect using the num_pixels property.

sqw.PixelData.@PixelData.get_pix_in_ranges(obj, abs_indices_starts, abs_indices_ends, recalculate_pix_range)

GET_PIX_IN_RANGES read pixels in the specified ranges Ranges are inclusive.

>> pix = get_pix_in_ranges([1, 12, 25], [6, 12, 27])

pix_starts Absolute indices of the starts of pixel ranges [Nx1 or 1xN array]. pix_ends Absolute indices of the ends of pixel ranges [Nx1 or 1xN array].

pix_out A PixelData object containing the pixels in the given ranges.

sqw.PixelData.@PixelData.get_pixels(obj, abs_pix_indices, varargin)

Retrieve the pixels at the given indices in the full pixel block, return a new PixelData object.

>> pix_out = pix.get_pixels(15640:19244) % retrieve pixels at indices 15640 to 19244

>> pix_out = pix.get_pixels([1, 0, 1]) % retrieve pixels at indices 1 and 3

The function attempts to mimic the behaviour you would see when indexing into a Matlab array. The difference being the returned object is a PixelData object and not an array.

This function may be useful if you want to extract data for a particular image bin.

abs_pix_indices A vector of positive integers or a vector of logicals.

The syntax for these indices attempts to replicate indexing into a regular Matlab array. You can use logical indices as well as normal indices, and you can index into the array “out-of-order”. However, you cannot use end, but it is possible to achieve the same effect using the num_pixels property.

Optional: ‘-ignore_range’ – if provided, new pix_object will not contain correct

pixel ranges

sqw.PixelData.@PixelData.mask(obj, mask_array, varargin)

MASK remove the pixels specified by the input logical array

You must specify exactly one return argument when calling this function.

mask_array A logical array specifying which pixels should be kept/removed

from the PixelData object. Must be of length equal to the number of pixels in ‘obj’ or equal in size to the ‘npix’ argument. A true/1 in the array indicates that the pixel at that index should be retained, a false/0 indicates the pixel should be removed.

npix (Optional)

Array of integers that specify how many times each value in mask_array should be replicated. This is useful for when masking all pixels contributing to a bin. Size must be equal to that of ‘mask_array’. E.g.:

mask_array = [ 0, 1, 1, 0, 1] npix = [ 3, 2, 2, 1, 2] full_mask = [0, 0, 0, 1, 1, 1, 1, 0, 1, 1]

The npix array must account for all pixels in the PixelData object i.e. sum(npix, ‘all’) == obj.num_pixels. It must also be the same dimensions as ‘mask_array’ i.e. all(size(mask_array) == size(npix)).

pix_out A PixelData object containing only non-masked pixels.

sqw.PixelData.@PixelData.move_to_page(obj, page_number, varargin)
Set the object to point at the given page number

This function does nothing if the object is not file-backed or is already on the given page

Inputs: page_number – page number to move to

Returns: page_number – the page this routine moved to total_num_pages – total number of pages, present in the file

sqw.PixelData.@PixelData.noisify(obj, varargin)

This is a method of the PixelData class. It is called from the noisify(sqw_object,[options]) function to allow the noisify functionality to be applied on a per-page basis to the sqw object’s pixel data. See noisify(sqw_object,…) for details of the input and the Herbert noisify function for details of how the input is used. This noisify adds random noise to the object’s signal array, and a fixed error bar to the variance array, paging as required. The options in varargin specify how that is done; see the above functions for details. For options where the signal absolute maximum is used in the noise scaling, a paged pre-scan of the signal provides the maximum over all pages.

obj The PixelData instance. varargs Options for random number distribution and noise magnitude

scaling.

pix_out If specified, returns a “noisified” copy of the input data

Otherwise it is a reference to the input object, which is “noisified” in place.

sqw.PixelData.@PixelData.recalc_pix_range(obj)

Recalculate pixels range in the situations, where the range for some reason appeared to be missing (i.e. loading pixels from old style files) or changed through private interface (for efficiency) and the internal integrity of the object has been violated.

returns obj for compatibility with recalc_pix_range method of combine_pixel_info class, which may be used instead of PixelData for the same purpose. recalc_pix_range is a normal Matlab value object (not a handle object), returning its changes in LHS

sqw.PixelData.@PixelData.set_data(obj, pix_fields, data, varargin)

SET_PIXELS Update the data on the given pixel data fields

The number of columns in ‘data’ must be equal to the number of fields in ‘pix_fields’. The number of rows in ‘data’ must be equal to the number of elements in ‘abs_pix_indices’.

Set the first 100 pixels’ signal and variance to zero

>> set_data({‘signal’, ‘variance’}, zeros(2, 100), 1:100);

pix_fields The name of a field, or a cell array of field names. data The data with which to set the given fields. abs_pix_indices The indices to set data on. If not specified all indices are

updated and ‘size(data, 2)’ must equal to obj.num_pixels.