xrayutilities.math package¶
Submodules¶
xrayutilities.math.algebra module¶
module providing analytic algebraic functions not implemented in scipy or any other dependency of xrayutilities. In particular the analytic solution of a quartic equation which is needed for the solution of the dynamic scattering equations.
xrayutilities.math.fit module¶
module with a function wrapper to scipy.optimize.leastsq for fitting of a 2D function to a peak or a 1D Gauss fit with the odr package
- xrayutilities.math.fit.fit_peak2d(x, y, data, start, drange, fit_function, maxfev=2000)[source]¶
fit a two dimensional function to a two dimensional data set e.g. a reciprocal space map.
- Parameters:
- xarray-like
first data coordinate (does not need to be regularly spaced)
- yarray-like
second data coordinate (does not need to be regularly spaced)
- dataarray-like
data set used for fitting (e.g. intensity at the data coordinates)
- startlist
set of starting parameters for the fit used as first parameter of function fit_function
- drangelist
limits for the data ranges used in the fitting algorithm, e.g. it is clever to use only a small region around the peak which should be fitted, i.e. [xmin, xmax, ymin, ymax]
- fit_functioncallable
function which should be fitted. Call signature must be
fit_function(x, y, *params) -> ndarray()
- Returns:
- fitparamlist
fitted parameters
- covarray-like
covariance matrix
- xrayutilities.math.fit.gauss_fit(xdata, ydata, iparams=None, maxit=300)[source]¶
Gauss fit function using odr-pack wrapper in scipy similar to https://github.com/tiagopereira/python_tips/wiki/Scipy%3A-curve-fitting
- Parameters:
- xdataarray-like
x-coordinates of the data to be fitted
- ydataarray-like
y-coordinates of the data which should be fit
- iparams: list, optional
initial paramters for the fit, determined automatically if not given
- maxitint, optional
maximal iteration number of the fit
- Returns:
- paramslist
the parameters as defined in function
Gauss1d(x, *param)
- sd_paramslist
For every parameter the corresponding errors are returned.
- itlimbool
flag to tell if the iteration limit was reached, should be False
- xrayutilities.math.fit.linregress(x, y)[source]¶
fast linregress to avoid usage of scipy.stats which is slow! NaN values in y are ignored by this function.
- Parameters:
- x, yarray-like
data coordinates and values
- Returns:
- ptuple
parameters of the linear fit (slope, offset)
- rsq: float
R^2 value
Examples
>>> (k, d), R2 = linregress([1, 2, 3], [3.3, 4.1, 5.05])
- xrayutilities.math.fit.multPeakFit(x, data, peakpos, peakwidth, dranges=None, peaktype='Gaussian', returnerror=False)[source]¶
function to fit multiple Gaussian/Lorentzian peaks with linear background to a set of data
- Parameters:
- xarray-like
x-coordinate of the data
- dataarray-like
data array with same length as x
- peakposlist
initial parameters for the peak positions
- peakwidthlist
initial values for the peak width
- drangeslist of tuples
list of tuples with (min, max) value of the data ranges to use. does not need to have the same number of entries as peakpos
- peaktype{‘Gaussian’, ‘Lorentzian’}
type of peaks to be used
- returnerrorbool
decides if the fit errors of pos, sigma, and amp are returned (default: False)
- Returns:
- poslist
peak positions derived by the fit
- sigmalist
peak width derived by the fit
- amplist
amplitudes of the peaks derived by the fit
- backgroundarray-like
background values at positions x
- if returnerror == True:
- sd_poslist
standard error of peak positions as returned by scipy.odr.Output
- sd_sigmalist
standard error of the peak width
- sd_amplist
standard error of the peak amplitude
- xrayutilities.math.fit.multPeakPlot(x, fpos, fwidth, famp, background, dranges=None, peaktype='Gaussian', fig='xu_plot', ax=None, fact=1.0)[source]¶
function to plot multiple Gaussian/Lorentz peaks with background values given by an array
- Parameters:
- xarray-like
x-coordinate of the data
- fposlist
positions of the peaks
- fwidthlist
width of the peaks
- famplist
amplitudes of the peaks
- backgroundarray-like
background values, same shape as x
- drangeslist of tuples
list of (min, max) values of the data ranges to use. does not need to have the same number of entries as fpos
- peaktype{‘Gaussian’, ‘Lorentzian’}
type of peaks to be used
- figint, str, or None
matplotlib figure number or name
- axmatplotlib.Axes
matplotlib axes as alternative to the figure name
- factfloat
factor to use as multiplicator in the plot
- xrayutilities.math.fit.peak_fit(xdata, ydata, iparams=None, peaktype='Gauss', maxit=300, background='constant', plot=False, func_out=False, debug=False)[source]¶
fit function using odr-pack wrapper in scipy similar to https://github.com/tiagopereira/python_tips/wiki/Scipy%3A-curve-fitting for Gauss, Lorentz or Pseudovoigt-functions
- Parameters:
- xdataarray_like
x-coordinates of the data to be fitted
- ydataarray_like
y-coordinates of the data which should be fit
- iparamslist, optional
initial paramters, determined automatically if not specified
- peaktype{‘Gauss’, ‘Lorentz’, ‘PseudoVoigt’,
‘PseudoVoigtAsym’, ‘PseudoVoigtAsym2’}, optional
type of peak to fit
- maxitint, optional
maximal iteration number of the fit
- background{‘constant’, ‘linear’}, optional
type of background function
- plotbool or str, optional
flag to ask for a plot to visually judge the fit. If plot is a string it will be used as figure name, which makes reusing the figures easier.
- func_outbool, optional
returns the fitted function, which takes the independent variables as only argument (f(x))
- Returns:
- paramslist
the parameters as defined in function Gauss1d/Lorentz1d/PseudoVoigt1d/ PseudoVoigt1dasym. In the case of linear background one more parameter is included!
- sd_paramslist
For every parameter the corresponding errors are returned.
- itlimbool
flag to tell if the iteration limit was reached, should be False
- fitfuncfunction, optional
the function used in the fit can be returned (see func_out).
xrayutilities.math.functions module¶
module with several common function needed in xray data analysis
- xrayutilities.math.functions.Debye1(x)[source]¶
function to calculate the first Debye function [1] as needed for the calculation of the thermal Debye-Waller-factor by numerical integration
\[D_1(x) = (1/x) \int_0^x t/(\exp(t)-1) dt\]- Parameters:
- xfloat
argument of the Debye function
- Returns:
- float
D1(x) float value of the Debye function
References
- xrayutilities.math.functions.Gauss1d(x, *p)[source]¶
function to calculate a general one dimensional Gaussian
- Parameters:
- xarray-like
coordinate(s) where the function should be evaluated
- plist
list of parameters of the Gaussian [XCEN, SIGMA, AMP, BACKGROUND] for information: SIGMA = FWHM / (2*sqrt(2*log(2)))
- Returns:
- array-like
the value of the Gaussian described by the parameters p at position x
Examples
Calling with a list of parameters needs a call looking as shown below (note the ‘*’) or explicit listing of the parameters
>>> Gauss1d(x, *p)
>>> import numpy >>> Gauss1d(numpy.linspace(0, 10, 10), 5, 1, 1e3, 0) array([3.72665317e-03, 5.19975743e-01, 2.11096565e+01, 2.49352209e+02, 8.56996891e+02, 8.56996891e+02, 2.49352209e+02, 2.11096565e+01, 5.19975743e-01, 3.72665317e-03])
- xrayutilities.math.functions.Gauss1dArea(*p)[source]¶
function to calculate the area of a Gauss function with neglected background
- Parameters:
- plist
list of parameters of the Gauss-function [XCEN, SIGMA, AMP, BACKGROUND]
- Returns:
- float
the area of the Gaussian described by the parameters p
- xrayutilities.math.functions.Gauss1d_der_p(x, *p)[source]¶
function to calculate the derivative of a Gaussian with respect the parameters p
for parameter description see Gauss1d
- xrayutilities.math.functions.Gauss1d_der_x(x, *p)[source]¶
function to calculate the derivative of a Gaussian with respect to x
for parameter description see Gauss1d
- xrayutilities.math.functions.Gauss2d(x, y, *p)[source]¶
function to calculate a general two dimensional Gaussian
- Parameters:
- x, yarray-like
coordinate(s) where the function should be evaluated
- plist
list of parameters of the Gauss-function [XCEN, YCEN, SIGMAX, SIGMAY, AMP, BACKGROUND, ANGLE]; SIGMA = FWHM / (2*sqrt(2*log(2))); ANGLE = rotation of the X, Y direction of the Gaussian in radians
- Returns:
- array-like
the value of the Gaussian described by the parameters p at position (x, y)
- xrayutilities.math.functions.Gauss2dArea(*p)[source]¶
function to calculate the area of a 2D Gauss function with neglected background
- Parameters:
- plist
list of parameters of the Gauss-function [XCEN, YCEN, SIGMAX, SIGMAY, AMP, ANGLE, BACKGROUND]
- Returns:
- float
the area of the Gaussian described by the parameters p
- xrayutilities.math.functions.Gauss3d(x, y, z, *p)[source]¶
function to calculate a general three dimensional Gaussian
- Parameters:
- x, y, zarray-like
coordinate(s) where the function should be evaluated
- plist
list of parameters of the Gauss-function [XCEN, YCEN, ZCEN, SIGMAX, SIGMAY, SIGMAZ, AMP, BACKGROUND];
SIGMA = FWHM / (2*sqrt(2*log(2)))
- Returns:
- array-like
the value of the Gaussian described by the parameters p at positions (x, y, z)
- xrayutilities.math.functions.Lorentz1d(x, *p)[source]¶
function to calculate a general one dimensional Lorentzian
- Parameters:
- xarray-like
coordinate(s) where the function should be evaluated
- plist
list of parameters of the Lorentz-function [XCEN, FWHM, AMP, BACKGROUND]
- Returns:
- array-like
the value of the Lorentian described by the parameters p at position (x, y)
- xrayutilities.math.functions.Lorentz1dArea(*p)[source]¶
function to calculate the area of a Lorentz function with neglected background
- Parameters:
- plist
list of parameters of the Lorentz-function [XCEN, FWHM, AMP, BACKGROUND]
- Returns:
- float
the area of the Lorentzian described by the parameters p
- xrayutilities.math.functions.Lorentz1d_der_p(x, *p)[source]¶
function to calculate the derivative of a Gaussian with respect the parameters p
for parameter description see Lorentz1d
- xrayutilities.math.functions.Lorentz1d_der_x(x, *p)[source]¶
function to calculate the derivative of a Gaussian with respect to x
for parameter description see Lorentz1d
- xrayutilities.math.functions.Lorentz2d(x, y, *p)[source]¶
function to calculate a general two dimensional Lorentzian
- Parameters:
- x, yarray-like
coordinate(s) where the function should be evaluated
- plist
list of parameters of the Lorentz-function [XCEN, YCEN, FWHMX, FWHMY, AMP, BACKGROUND, ANGLE]; ANGLE = rotation of the X, Y direction of the Lorentzian in radians
- Returns:
- array-like
the value of the Lorentian described by the parameters p at position (x, y)
- xrayutilities.math.functions.NormGauss1d(x, *p)[source]¶
function to calculate a normalized one dimensional Gaussian
- Parameters:
- xarray-like
coordinate(s) where the function should be evaluated
- plist
list of parameters of the Gaussian [XCEN, SIGMA]; for information: SIGMA = FWHM / (2*sqrt(2*log(2)))
- Returns:
- array-like
the value of the normalized Gaussian described by the parameters p at position x
- xrayutilities.math.functions.NormLorentz1d(x, *p)[source]¶
function to calculate a normalized one dimensional Lorentzian
- Parameters:
- xarray-like
coordinate(s) where the function should be evaluated
- plist
list of parameters of the Lorentzian [XCEN, FWHM]
- Returns:
- array-like
the value of the normalized Lorentzian described by the parameters p at position x
- xrayutilities.math.functions.PseudoVoigt1d(x, *p)[source]¶
function to calculate a pseudo Voigt function as linear combination of a Gauss and Lorentz peak
- Parameters:
- xarray-like
coordinate(s) where the function should be evaluated
- plist
list of parameters of the pseudo Voigt-function [XCEN, FWHM, AMP, BACKGROUND, ETA]; ETA: 0 …1 0 means pure Gauss and 1 means pure Lorentz
- Returns:
- array-like
the value of the PseudoVoigt described by the parameters p at position x
- xrayutilities.math.functions.PseudoVoigt1dArea(*p)[source]¶
function to calculate the area of a pseudo Voigt function with neglected background
- Parameters:
- plist
list of parameters of the Lorentz-function [XCEN, FWHM, AMP, BACKGROUND, ETA]; ETA: 0 …1 0 means pure Gauss and 1 means pure Lorentz
- Returns:
- float
the area of the PseudoVoigt described by the parameters p
- xrayutilities.math.functions.PseudoVoigt1d_der_p(x, *p)[source]¶
function to calculate the derivative of a PseudoVoigt with respect the parameters p
for parameter description see PseudoVoigt1d
- xrayutilities.math.functions.PseudoVoigt1d_der_x(x, *p)[source]¶
function to calculate the derivative of a PseudoVoigt with respect to x
for parameter description see PseudoVoigt1d
- xrayutilities.math.functions.PseudoVoigt1dasym(x, *p)[source]¶
function to calculate an asymmetric pseudo Voigt function as linear combination of asymmetric Gauss and Lorentz peak
- Parameters:
- xarray-like
coordinate(s) where the function should be evaluated
- plist
list of parameters of the pseudo Voigt-function [XCEN, FWHMLEFT, FWHMRIGHT, AMP, BACKGROUND, ETA]; ETA: 0 …1 0 means pure Gauss and 1 means pure Lorentz
- Returns:
- array-like
the value of the PseudoVoigt described by the parameters p at position x
- xrayutilities.math.functions.PseudoVoigt1dasym2(x, *p)[source]¶
function to calculate an asymmetric pseudo Voigt function as linear combination of asymmetric Gauss and Lorentz peak
- Parameters:
- xnaddray
coordinate(s) where the function should be evaluated
- plist
list of parameters of the pseudo Voigt-function [XCEN, FWHMLEFT, FWHMRIGHT, AMP, BACKGROUND, ETALEFT, ETARIGHT]; ETA: 0 …1 0 means pure Gauss and 1 means pure Lorentz
- Returns:
- array-like
the value of the PseudoVoigt described by the parameters p at position x
- xrayutilities.math.functions.PseudoVoigt2d(x, y, *p)[source]¶
function to calculate a pseudo Voigt function as linear combination of a Gauss and Lorentz peak in two dimensions
- Parameters:
- x, yarray-like
coordinate(s) where the function should be evaluated
- plist
list of parameters of the pseudo Voigt-function [XCEN, YCEN, FWHMX, FWHMY, AMP, BACKGROUND, ANGLE, ETA]; ETA: 0 …1 0 means pure Gauss and 1 means pure Lorentz
- Returns:
- array-like
the value of the PseudoVoigt described by the parameters p at position (x, y)
- xrayutilities.math.functions.TwoGauss2d(x, y, *p)[source]¶
function to calculate two general two dimensional Gaussians
- Parameters:
- x, yarray-like
coordinate(s) where the function should be evaluated
- plist
list of parameters of the Gauss-function [XCEN1, YCEN1, SIGMAX1, SIGMAY1, AMP1, ANGLE1, XCEN2, YCEN2, SIGMAX2, SIGMAY2, AMP2, ANGLE2, BACKGROUND]; SIGMA = FWHM / (2*sqrt(2*log(2))) ANGLE = rotation of the X, Y direction of the Gaussian in radians
- Returns:
- array-like
the value of the Gaussian described by the parameters p at position (x, y)
- xrayutilities.math.functions.heaviside(x)[source]¶
Heaviside step function for numpy arrays
- Parameters:
- x: scalar or array-like
argument of the step function
- Returns:
- int or array-like
Heaviside step function evaluated for all values of x with datatype integer
- xrayutilities.math.functions.kill_spike(data, threshold=2.0, offset=None)[source]¶
function to smooth single data points which differ from the average of the neighboring data points by more than the threshold factor or more than the offset value. Such spikes will be replaced by the mean value of the next neighbors.
Warning
Use this function carefully not to manipulate your data!
- Parameters:
- dataarray-like
1d numpy array with experimental data
- thresholdfloat or None
threshold factor to identify outlier data points. If None it will be ignored.
- offsetNone or float
offset value to identify outlier data points. If None it will be ignored.
- Returns:
- array-like
1d data-array with spikes removed
- xrayutilities.math.functions.multPeak1d(x, *args)[source]¶
function to calculate the sum of multiple peaks in 1D. the peaks can be of different type and a background function (polynom) can also be included.
- Parameters:
- xarray-like
coordinate where the function should be evaluated
- argslist
list of peak/function types and parameters for every function type two arguments need to be given first the type of function as string with possible values ‘g’: Gaussian, ‘l’: Lorentzian, ‘v’: PseudoVoigt, ‘a’: asym. PseudoVoigt, ‘p’: polynom the second type of arguments is the tuple/list of parameters of the respective function. See documentation of math.Gauss1d, math.Lorentz1d, math.PseudoVoigt1d, math.PseudoVoigt1dasym, and numpy.polyval for details of the different function types.
- Returns:
- array-like
value of the sum of functions at position x
- xrayutilities.math.functions.multPeak2d(x, y, *args)[source]¶
function to calculate the sum of multiple peaks in 2D. the peaks can be of different type and a background function (polynom) can also be included.
- Parameters:
- x, yarray-like
coordinates where the function should be evaluated
- argslist
list of peak/function types and parameters for every function type two arguments need to be given first the type of function as string with possible values ‘g’: Gaussian, ‘l’: Lorentzian, ‘v’: PseudoVoigt, ‘c’: constant the second type of arguments is the tuple/list of parameters of the respective function. See documentation of math.Gauss2d, math.Lorentz2d, math.PseudoVoigt2d for details of the different function types. The constant accepts a single float which will be added to the data
- Returns:
- array-like
value of the sum of functions at position (x, y)
xrayutilities.math.misc module¶
- xrayutilities.math.misc.center_of_mass(pos, data, background='none', full_output=False)[source]¶
function to determine the center of mass of an array
- Parameters:
- posarray-like
position of the data points
- dataarray-like
data values
- background{‘none’, ‘constant’, ‘linear’}
type of background, either ‘none’, ‘constant’ or ‘linear’
- full_outputbool
return background cleaned data and background-parameters
- Returns:
- float
center of mass position
- xrayutilities.math.misc.fwhm_exp(pos, data)[source]¶
function to determine the full width at half maximum value of experimental data. Please check the obtained value visually (noise influences the result)
- Parameters:
- posarray-like
position of the data points
- dataarray-like
data values
- Returns:
- float
fwhm value
xrayutilities.math.transforms module¶
- xrayutilities.math.transforms.ArbRotation(axis, alpha, deg=True)[source]¶
Returns a transform that represents a rotation around an arbitrary axis by the angle alpha. positive rotation is anti-clockwise when looking from positive end of axis vector
- Parameters:
- axislist or array-like
rotation axis
- alphafloat
rotation angle in degree (deg=True) or in rad (deg=False)
- degbool
determines the input format of ang (default: True)
- Returns:
- Transform
- class xrayutilities.math.transforms.AxisToZ(newzaxis)[source]¶
Bases:
CoordinateTransform
Creates a coordinate transformation to move a certain axis to the z-axis. The rotation is done along the great circle. The x-axis of the new coordinate frame is created to be normal to the new and original z-axis. The new y-axis is create in order to obtain a right handed coordinate system.
- class xrayutilities.math.transforms.AxisToZ_keepXY(newzaxis)[source]¶
Bases:
CoordinateTransform
Creates a coordinate transformation to move a certain axis to the z-axis. The rotation is done along the great circle. The x-axis/y-axis of the new coordinate frame is created to be similar to the old x and y directions. This variant of AxisToZ assumes that the new Z-axis has its main component along the Z-direction
- class xrayutilities.math.transforms.CoordinateTransform(v1, v2, v3)[source]¶
Bases:
Transform
Create a Transformation object which transforms a point into a new coordinate frame. The new frame is determined by the three vectors v1/norm(v1), v2/norm(v2) and v3/norm(v3), which need to be orthogonal!
- class xrayutilities.math.transforms.Transform(matrix)[source]¶
Bases:
object
- __call__(args, rank=1)[source]¶
transforms a vector, matrix or tensor of rank 4 (e.g. elasticity tensor)
- Parameters:
- argslist or array-like
object to transform, list or numpy array of shape (…, n) (…, n, n), (…, n, n, n, n) where n is the size of the transformation matrix.
- rankint
rank of the supplied object. allowed values are 1, 2, and 4
- property imatrix¶
- inverse(args, rank=1)[source]¶
performs inverse transformation a vector, matrix or tensor of rank 4
- Parameters:
- argslist or array-like
object to transform, list or numpy array of shape (…, n) (…, n, n), (…, n, n, n, n) where n is the size of the transformation matrix.
- rankint
rank of the supplied object. allowed values are 1, 2, and 4
- xrayutilities.math.transforms.VecAngle((v1.v2)/(norm(v1)*norm(v2)))[source]¶
alpha = acos((v1.v2)/(norm(v1)*norm(v2)))
- Parameters:
- v1, v2list or array-like
input vector(s), either one vector or an array of vectors with shape (n, 3)
- deg: bool, optional
True: return result in degree, False: in radiants (default: False)
- Returns:
- float or ndarray
the angle included by the two vectors v1 and v2, either a single float or an array with shape (n, )
- xrayutilities.math.transforms.VecCross(v1, v2, out=None)[source]¶
Calculate the vector cross product.
- Parameters:
- v1, v2list or array-like
input vector(s), either one vector or an array of vectors with shape (n, 3)
- outlist or array-like, optional
output vector
- Returns:
- ndarray
cross product either of shape (3, ) or (n, 3)
- xrayutilities.math.transforms.VecDot(v1, v2)[source]¶
Calculate the vector dot product.
- Parameters:
- v1, v2list or array-like
input vector(s), either one vector or an array of vectors with shape (n, 3)
- Returns:
- float or ndarray
innter product of the vectors, either a single float or (n, )
- xrayutilities.math.transforms.VecNorm(v)[source]¶
Calculate the norm of a vector.
- Parameters:
- vlist or array-like
input vector(s), either one vector or an array of vectors with shape (n, 3)
- Returns:
- float or ndarray
vector norm, either a single float or shape (n, )
- xrayutilities.math.transforms.VecUnit(v)[source]¶
Calculate the unit vector of v.
- Parameters:
- vlist or array-like
input vector(s), either one vector or an array of vectors with shape (n, 3)
- Returns:
- ndarray
unit vector of v, either shape (3, ) or (n, 3)
- xrayutilities.math.transforms.XRotation(alpha, deg=True)[source]¶
Returns a transform that represents a rotation about the x-axis by an angle alpha. If deg=True the angle is assumed to be in degree, otherwise the function expects radiants.
- xrayutilities.math.transforms.YRotation(alpha, deg=True)[source]¶
Returns a transform that represents a rotation about the y-axis by an angle alpha. If deg=True the angle is assumed to be in degree, otherwise the function expects radiants.
- xrayutilities.math.transforms.ZRotation(alpha, deg=True)[source]¶
Returns a transform that represents a rotation about the z-axis by an angle alpha. If deg=True the angle is assumed to be in degree, otherwise the function expects radiants.
- xrayutilities.math.transforms.distance(x, y, z, point, vec)[source]¶
calculate the distance between the point (x, y, z) and the line defined by the point and vector vec
- Parameters:
- xfloat or ndarray
x coordinate(s) of the point(s)
- yfloat or ndarray
y coordinate(s) of the point(s)
- zfloat or ndarray
z coordinate(s) of the point(s)
- pointtuple, list or ndarray
3D point on the line to which the distance should be calculated
- vectuple, list or ndarray
3D vector defining the propergation direction of the line
- xrayutilities.math.transforms.getSyntax(vec)[source]¶
returns vector direction in the syntax ‘x+’ ‘z-’ or equivalents therefore works only for principle vectors of the coordinate system like e.g. [1, 0, 0] or [0, 2, 0]
- Parameters:
- veclist or array-like
vector of length 3
- Returns:
- str
vector string following the synthax [xyz][+-]
- xrayutilities.math.transforms.getVector(string)[source]¶
returns unit vector along a rotation axis given in the syntax ‘x+’ ‘z-’ or equivalents
- Parameters:
- string: str
vector string following the synthax [xyz][+-]
- Returns:
- ndarray
vector along the given direction
- xrayutilities.math.transforms.mycross(vec, mat)[source]¶
function implements the cross-product of a vector with each column of a matrix
- xrayutilities.math.transforms.rotarb(vec, axis, ang, deg=True)[source]¶
function implements the rotation around an arbitrary axis by an angle ang positive rotation is anti-clockwise when looking from positive end of axis vector
- Parameters:
- veclist or array-like
vector to rotate
- axislist or array-like
rotation axis
- angfloat
rotation angle in degree (deg=True) or in rad (deg=False)
- degbool
determines the input format of ang (default: True)
- Returns:
- rotvecrotated vector as numpy.array
Examples
>>> rotarb([1, 0, 0], [0, 0, 1], 90) array([6.123234e-17, 1.000000e+00, 0.000000e+00])