SurfaceTopography.Nonuniform package

Submodules

SurfaceTopography.Nonuniform.Autocorrelation module

Height-difference autocorrelation functions for nonuniform line scans

SurfaceTopography.Nonuniform.Autocorrelation.height_height_autocorrelation(self, distances=None)

Compute the one-dimensional height-height autocorrelation function (ACF).

This function treats the nonuniform line scan as a piece-wise function of straight lines between the data points. The ACF is computed exactly for this piece-wise linear interpolation of the data.

Parameters:
  • self (NonuniformLineScan) – Container storing the nonuniform line scan.

  • r (array_like) – Array containing distances for which to compute the ACF. If no array is given, the function will automatically construct an array with equally spaced distances. (Default: None)

Returns:

  • distances (array) – Distances. (Units: length)

  • A (array) – Autocorrelation function. (Units: length**2)

SurfaceTopography.Nonuniform.Autocorrelation.height_difference_autocorrelation(self, reliable=True, algorithm='fft', distances=None, nb_interpolate=5, short_cutoff=<function mean>, resampling_method='bin-average', collocation='log', nb_points=None, nb_points_per_decade=10)

Compute the one-dimensional height-difference autocorrelation function (ACF).

This function treats the nonuniform line scan as a piece-wise function of straight lines between the data points. The ACF is computed exactly for this piece-wise linear interpolation of the data.

Parameters:
  • self (NonuniformLineScan) – Container storing the nonuniform line scan.

  • reliable (bool, optional) – Only return data deemed reliable. (Default: True)

  • algorithm (str) – Algorithm to compute autocorrelation. * ‘fft’: Interpolates the nonuniform line scan on a grid and then uses the FFT to compute the autocorrelation. Scales O(N log N) * ‘brute-force’: Brute-force computation using between line segements. Scale O(N^2 M) where M is the number of distance point. (Default: ‘fft’)

  • distances (array_like) – Array containing distances for which to compute the ACF. If no array is given, the function will automatically construct an array with equally spaced distances. Can be used only if ‘brute-force’ algorithm is used. (Default: None)

  • nb_interpolate (int) – Number of grid points to put between closest points on surface. Only used for ‘fft’ algorithm. (Default: 5)

  • short_cutoff (function, optional) – Function that determines how the short cutoff for the returned PSD is computed. If set to None, the full PSD of the interpolated data is returned. If the user passes a function, that function is called with an array that contains the distances between the points of the uniform line scan. Pass np.max, np.mean or np.min to use the maximum, mean or minimum of that distance as the cutoff.

  • resampling_method (str, optional) – Method can be None for no resampling (return on the grid of the data) ‘bin-average’ for simple bin averaging and ‘gaussian-process’ for Gaussian process regression. (Default: None)

  • collocation ({'log', 'quadratic', 'linear', array_like}, optional) – Resampling grid. Specifying ‘log’ yields collocation points equally spaced on a log scale, ‘quadratic’ yields bins with similar number of data points and ‘linear’ yields linear bins. Alternatively, it is possible to explicitly specify the bin edges. If bin_edges are explicitly specified, then rmax and nbins is ignored. (Default: ‘log’)

  • nb_points (int, optional) – Number of bins for averaging. Bins are automatically determined if set to None. (Default: None)

  • nb_points_per_decade (int, optional) – Number of points per decade for log-spaced collocation points. (Default: None)

Returns:

  • distances (array) – Distances. (Units: length)

  • acf (array) – Autocorrelation function. (Units: length**2)

SurfaceTopography.Nonuniform.BearingArea module

Bearing area curve, also known as Abbott-Firestone curve or cumulative distribution function of the surface heights.

class SurfaceTopography.Nonuniform.BearingArea.NonuniformBearingArea(x, h)

Bases: object

Accelerated bearing area calculation for nonuniform line scans.

__init__(x, h)
bounds(heights)

Compute bounds on the bearing area for a specific height. These bounds can be computed O(log(N)) time, where N is the number of elements in the line scan.

Parameters:

heights (float or np.ndarray) – Heights for which to compute the bearing area.

Returns:

  • lower_bound (float or np.ndarray) – Lower bound on fractional area above a threshold height.

  • upper_bound (float or np.ndarray) – Upper bound on fractional area above a threshold height.

property min
property max
SurfaceTopography.Nonuniform.BearingArea.bearing_area(self, heights=None)

Compute the bearing area for a specific height.

The bearing area as a function of height is also known as the Abbott-Firestone curve. If expressed as a fractional area, it is the cumulative distribution function of the surface heights. This function returns this fractional area.

Parameters:
  • self (NonuniformLineScan) – Line scan container object.

  • heights (float or np.ndarray, optional) – Heights for which to compute the bearing area. (Default: None)

Returns:

  • fractional_area (float or np.ndarray) – Fractional area above the threshold height, if height is given.

  • bearing_area (NonuniformBearingArea) – Instance of NonuniformBearingArea class that caches the bearing area calculation.

SurfaceTopography.Nonuniform.Converters module

Utility function for conversion between uniform and nonuniform representations.

class SurfaceTopography.Nonuniform.Converters.UniformlyInterpolatedLineScan(topography, nb_points=None, padding=0, nb_interpolate=None, pixel_size=None, info={})

Bases: DecoratedUniformTopography

Interpolate a topography onto a uniform grid.

__init__(topography, nb_points=None, padding=0, nb_interpolate=None, pixel_size=None, info={})

Convert a nonuniform line scan to a uniform line scan by explicit linear interpolation between the discrete coordinates of the nonuniform scan.

nb_points, nb_interpolate or pixel_size needs to be specified.

Parameters:
  • topography (SurfaceTopography or UniformLineScan) – SurfaceTopography to interpolate.

  • nb_points (int, optional) – Number of equidistant grid points. Number of points will be automatically computed from the mean grid spacing and nb_interpolate if set to None. (Default: None)

  • padding (int, optional) – Number of padding grid points, zeros appended to the data. (Default: 0)

  • nb_interpolate (int, optional) – Number of grid points to between closest points on surface. (Default: None)

  • pixel_size (float, optional) – Grid spacing. (Default: None)

property dim

Returns 1 for line scans and 2 for topography maps.

property physical_sizes

Return the physical sizes of the topography.

property is_periodic

Return whether the topography is periodically repeated at the boundaries.

property is_uniform
property nb_grid_pts

Return nb_grid_pts, i.e. number of pixels, of the topography.

property pixel_size
property area_per_pt
property has_undefined_data
positions()

Returns array containing the lateral positions.

heights()

Computes the rescaled profile.

SurfaceTopography.Nonuniform.Derivative module

Compute derivatives of nonuniform line scans.

SurfaceTopography.Nonuniform.Derivative.derivative(self, n, scale_factor=None, distance=None, interpolation='linear', progress_callback=None)

Compute derivative of nonuniform line-scan. Function assumes nonperiodic topographies.

First derivative: Central differences.

Second derivative: Expand \(h(x+\Delta x_+)\) and \((x-\Delta x_-)\) up to second order in the grid spacing \(\Delta x_+\) and \(\Delta x_+\). Then \(\Delta x_- f(x+\Delta x_+) + \Delta x_+ f(x+\Delta x_-)\) yields:

\[\frac{d^2h}{dx^2} \approx 2 \frac{\Delta x_-\left[f(x+\Delta x_+)-f(x)\right] + \Delta x_+\left[f(x+\Delta x_-)-f(x)\right]}{\Delta x_+\Delta x_-(\Delta x_++\Delta x_-)}\]
Parameters:
  • self (SurfaceTopography.NonuniformLineScan) – Object containing height information.

  • n (int) – Number of times the derivative is taken.

  • scale_factor (int or list of ints or list of tuples of ints, optional) – Scale factors are not supported by nonuniform line scans. (Default: None)

  • distance (float or list of floats, optional) – Explicit distance scale for computation of the derivative. Note that the distance specifies the overall length of the stencil of lowest truncation order, not the effective grid spacing used by this stencil. The scale factor is then given by distance / (n * px) where n is the order of the derivative and px the grid spacing. (Default: None)

  • interpolation (str, optional) – Only ‘linear’ interpolation is supported by nonuniform line scans. (Default: ‘linear’)

  • progress_callback (func, optional) – Function taking iteration and the total number of iterations as arguments for progress reporting. (Default: None)

Returns:

derivative – Array with derivative values. Length of array is reduced by \(n\) with respect to the input array for the \(n\)-th derivative.

Return type:

np.ndarray

SurfaceTopography.Nonuniform.Detrending module

Helper functions to compute trends of surfaces

SurfaceTopography.Nonuniform.Detrending.polyfit(self, deg)

Compute the detrending plane that, if subtracted, minimizes the rms height of the surface. The detrending plane is parameterized as a polynomial:

\[p(x) = \sum_{k=0}^n a_k x^k\]

The values of \(a_k\) are returned by this function.

The rms height of the surface is given by (see rms_height)

\[h_\text{rms}^2 = \frac{1}{3L} \sum_{i=0}^{N-2} \left( h_i^2 + h_{i+1}^2 + h_i h_{i+1} \right) \Delta x_i\]

where \(N\) is the total number of data points. Hence we need to solve the following minimization problem:

\[\min_{\{a_k\}} \left\{ \frac{1}{3L} \sum_{i=0}^{N-2} \left[ (h_i - p(x_i))^2 + (h_{i+1} - p(x_{i+1}))^2 + (h_i - p(x_i))(h_{i+1} - p(x_{i+1})) \right] \Delta x_i \right\}\]

This gives the system of linear equations (one for each \(k\))

\[\sum_{i=0}^{N-2} \left( \left[ 2h_i + h_{i+1} \right] x_i^k + \left[ 2h_{i+1} + h_i \right] x_{i+1}^k \right) \Delta x_i = \sum_{l=0}^n a_l \sum_{i=0}^{N-2} \left( 2x_i^{k+l} + 2x_{i+1}^{k+l} + x_i^k x_{i+1}^l + x_{i+1}^k x_i^l \right) \Delta x_i\]
Parameters:
  • self (NonuniformLineScan) – SurfaceTopography object containing height information.

  • deg (int) – Degree of polynomial \(n\).

Returns:

a – Array with coefficients \(a_k\).

Return type:

array

class SurfaceTopography.Nonuniform.Detrending.DetrendedNonuniformTopography(topography, detrend_mode='height', coeffs=None, info={})

Bases: DecoratedNonuniformTopography

Remove trends from a topography. This is achieved by fitting polynomials to the topography data to extract trend lines. The resulting topography is then detrended by substracting these trend lines.

__init__(topography, detrend_mode='height', coeffs=None, info={})
Parameters:
  • topography (SurfaceTopography) – SurfaceTopography to be detrended.

  • detrend_mode (str) – ‘mean’: center the topography to its mean, no trend correction. ‘median’: center the topography to its median, no trend correction. ‘rms-tilt’: adjust slope such that rms height is minimized. ‘mad-tilt’: adjust slope such that mad height is minimized. ‘slope’: adjust slope such that rms slope is minimized. ‘rms-curvature’: adjust slope and curvature such that rms height is minimized. (Default: ‘height’)

  • coeffs (array-like, optional) – Coefficients of the detrending plane. If not given, they are computed from the detrend_mode. (Default: None)

property coeffs
property detrend_mode
property is_periodic

SurfaceTopography stays periodic only after detrend mode “center”. Otherwise the detrended SurfaceTopography is non-periodic.

property x_range
positions()

Returns array containing the lateral positions.

heights()

Computes the combined profile.

stringify_plane(fmt=<function DetrendedNonuniformTopography.<lambda>>)
property curvatures

SurfaceTopography.Nonuniform.Interpolation module

SurfaceTopography.Nonuniform.Interpolation.interpolate_linear(self)

Returns a linear interpolation function based on the topography’s heights.

SurfaceTopography.Nonuniform.Interpolation.interpolate_cubic(self)

Returns a linear interpolation function based on the topography’s heights.

SurfaceTopography.Nonuniform.PowerSpectrum module

Power-spectral density for nonuniform topographies.

SurfaceTopography.Nonuniform.PowerSpectrum.sinc(x)

sinc function

SurfaceTopography.Nonuniform.PowerSpectrum.dsinc(x)

Derivative of the sinc function

SurfaceTopography.Nonuniform.PowerSpectrum.ft_rectangle(q)

Fourier transform of a rectangle, \(f(x) = 1\) for \(|x| < 1/2\),

\[\tilde{f}(q) = 2\sin(q/2)/q = sinc(q/2)\]
SurfaceTopography.Nonuniform.PowerSpectrum.ft_one_sided_triangle(q)

Fourier transform of the one-sided triangle, \(f(x) = x\) for \(|x| < 1/2\),

\[\tilde{f}(q) = i [\cos(q/2)/q - 2\sin(q/2)/q^2] = d sinc(q/2)/dq\]

Returned value does not contain the factor of \(i\), i.e. this function is real-valued.

SurfaceTopography.Nonuniform.PowerSpectrum.apply_window(x, y, window=None)
SurfaceTopography.Nonuniform.PowerSpectrum.power_spectrum(self, reliable=True, algorithm='fft', wavevectors=None, nb_interpolate=5, short_cutoff=<function mean>, window=None, resampling_method='bin-average', collocation='log', nb_points=None, nb_points_per_decade=10)

Compute power-spectral density (PSD) for a nonuniform topography. The topography is assumed to be given by a series of points connected by straight lines.

Parameters:
  • self (NonuniformLineScan) – Container storing the nonuniform line scan.

  • reliable (bool, optional) – Only return data deemed reliable. (Default: True)

  • algorithm (str) – Algorithm to compute autocorrelation. * ‘fft’: Interpolates the nonuniform line scan on a grid and then uses the FFT to compute the autocorrelation. Scales O(N log N) * ‘brute-force’: Brute-force computation using between line segements. Scale O(N^2 M) where M is the number of distance point. (Default: ‘fft’)

  • wavevectors (array, optional) – Wavevectors at which to compute the PSD. If omitted, wavevectors are equally spaced with a spacing that corresponds to \(2\pi/\lambda\) where \(\lambda\) is the shortest distance between two points in the x-array. (Default: None)

  • nb_interpolate (int, optional) – Number of grid points to put between closest points on surface. Only used for ‘fft’ algorithm. (Default: 5)

  • short_cutoff (function, optional) – Function that determines how the short cutoff for the returned PSD is computed. If set to None, the full PSD of the interpolated data is returned. If the user passes a function, that function is called with an array that contains the distances between the points of the uniform line scan. Pass np.max, np.mean or np.min to use the maximum, mean or minimum of that distance as the cutoff.

  • window (str, optional) – Name of the window function to apply before computing the PSD. Presently only supports Hann window (‘hann’) or no window (None or ‘None’). Default: no window for periodic Topographies, “hann” window for nonperiodic Topographies

  • resampling_method (str, optional) – Method can be None for no resampling (return on the grid of the data) ‘bin-average’ for simple bin averaging and ‘gaussian-process’ for Gaussian process regression. (Default: None)

  • collocation ({'log', 'quadratic', 'linear', array_like}, optional) – Resampling grid. Specifying ‘log’ yields collocation points equally spaced on a log scale, ‘quadratic’ yields bins with similar number of data points and ‘linear’ yields linear bins. Alternatively, it is possible to explicitly specify the bin edges. If bin_edges are explicitly specified, then rmax and nbins is ignored. (Default: ‘log’)

  • nb_points (int, optional) – Number of bins for averaging. Bins are automatically determined if set to None. (Default: None)

  • nb_points_per_decade (int, optional) – Number of points per decade for log-spaced collocation points. (Default: None)

Returns:

  • wavevectors (array) – Wavevector array at which the PSD has been computed. (Unit: length)

  • psd (array) – PSD values. (Unit: length**-3)

SurfaceTopography.Nonuniform.ScalarParameters module

Functions computing scalar roughness parameters

SurfaceTopography.Nonuniform.ScalarParameters.moment(topography, alpha)

Computes the n-th moment of the height fluctuation of the line scan of length \(L\):

\[\langle h^\alpha \rangle = \frac{1}{L} \int_0^L dx\, h^\alpha(x)\right\]

This function approximates the topography between data points as piece-wise linear. The piece-wise linear section between point \(i\) and point \(i+1\) contributes

\[\frac{1}{\alpha+1} \frac{x_{i+1} - x_i}{L}\frac{h_{i+1}^{\alpha+1} - h_i^{\alpha+1}}{h_{i+1} - h_i}\]

to the above integral.

Parameters:
  • topography (NonuniformLineScan) – SurfaceTopography object containing height information.

  • alpha (int) – Order of moment.

Returns:

moment – Root-mean square height.

Return type:

float or array

SurfaceTopography.Nonuniform.ScalarParameters.rms_height(self)

Computes root-mean square height fluctuation of the line scan:

\[h_\text{rms} = \left[\frac{1}{L} \int_0^L dx\, h^2(x)\right]^{1/2}\]

Function approximates topography between data points as piece-wise linear. The piece-wise linear section between point \(i\) and point \(i+1\) contributes

\[\int_{0}^{\Delta x_i} dx\, \left( h_i + \frac{\Delta h_i}{\Delta x_i} x \right)^2 = \frac{1}{3} \left( h_i^2 + h_{i+1}^2 + h_i h_{i+1} \right) \Delta x_i\]

to the above integral, where \(\Delta x_i=x_{i+1}-x_i\) and \(\Delta h_i=h_{i+1}-h_i\).

Parameters:

self (NonuniformLineScan) – SurfaceTopography object containing height information.

Returns:

rms_height – Root-mean square height.

Return type:

float or array

SurfaceTopography.Nonuniform.ScalarParameters.rms_slope(self)

Computes root-mean square slope fluctuation of the line scan:

\[h^\prime_\text{rms} = \left[ \frac{1}{L} \int_0^L dx\, \left(\frac{\partial h}{\partial x}\right)^2 \right]^{1/2}\]

Function approximates topography between data points as piece-wise linear. The piece-wise linear section between point \(i\) and point \(i+1\) contributes

\[\int_{0}^{\Delta x_i} dx\, \left( \frac{\Delta h_i}{\Delta x_i} \right)^2 = \frac{\Delta h_i^2}{\Delta x_i}\]

where \(\Delta x_i=x_{i+1}-x_i\) and \(\Delta h_i=h_{i+1}-h_i\) to the above integral.

Parameters:

self (NonuniformLineScan) – SurfaceTopography object containing height information.

Returns:

rms_slope – Root-mean square slope.

Return type:

float

SurfaceTopography.Nonuniform.ScalarParameters.rms_curvature(self)

Computes root-mean square slope fluctuation of the line scan:

Parameters:

self (NonuniformLineScan) – SurfaceTopography object containing height information.

Returns:

rms_slope – Root-mean square slope.

Return type:

float

SurfaceTopography.Nonuniform.VariableBandwidth module

Variable bandwidth analysis for nonuniform topographies

SurfaceTopography.Nonuniform.VariableBandwidth.checkerboard_detrend_profile(self, subdivisions, tol=1e-06)

Perform tilt correction (and substract mean value) in each individual rectangle of a checkerboard decomposition of the surface. This is identical to subdividing the surface into individual, nonoverlapping rectangles and performing individual tilt corrections on them.

The main application of this function is to carry out a variable bandwidth analysis of the surface.

Parameters:
  • self (NonuniformLineScan) – Container storing the uniform topography map

  • subdivisions (int) – Number of subdivisions.

  • tol (float) – Tolerance for searching for existing data points at domain boundaries. (Default: 1e-6)

Returns:

subdivided_line_scans – List with new, subdivided and detrended line scans.

Return type:

list of NonuniformLineScan

SurfaceTopography.Nonuniform.VariableBandwidth.variable_bandwidth_from_profile(self, quantities='bh', reliable=True, resampling_method=None, nb_grid_pts_cutoff=4)

Perform a variable bandwidth analysis by computing the mean root-mean-square height within increasingly finer subdivisions of the line scan.

Parameters:
  • self (obj:NonuniformLineScan) – Container storing the uniform topography map

  • quantities (str, optional) –

    Specification of return tuple, where each string character stand for specific quantity. Possible quantities are

    • ’m’: Magnification (Unit: dimensionless)

    • ’b’: Bandwidth (Unit: length)

    • ’h’: RMS height (Unit: length)

    For example, ‘mbh’ return a tuple with the three entries magnification, bandwidth, rms height. (Default: ‘bh’)

  • reliable (bool, optional) – Only return data deemed reliable. (Default: True)

  • resampling_method (str, optional) – Can only be None for the variable bandwidth analysis. (Default: None)

  • nb_grid_pts_cutoff (int) – Minimum number of data points to allow for subdivision. The analysis will automatically analyze subdivision down to this nb_grid_pts.

Returns:

  • magnifications (array) – Array containing the magnifications.

  • bandwidths (array) – Array containing the bandwidths, here the physical_sizes of the subdivided topography.

  • rms_heights (array) – Array containing the rms height corresponding to the respective magnification.

SurfaceTopography.Nonuniform.common module

Bin for small common helper function and classes for nonuniform topographies.

SurfaceTopography.Nonuniform.common.bandwidth(self)

Computes lower and upper bound of bandwidth, i.e. of the wavelengths or length scales occurring on a topography. The lower end of the bandwidth is given by the mean of the spacing of the individual points on the line scan. The upper bound is given by the overall length of the line scan.

Returns:

  • lower_bound (float) – Lower bound of the bandwidth.

  • upper_bound (float) – Upper bound of the bandwidth.

Module contents

Module containing all functions operating on uniform topographies