Digital Image Correlation Post Processing \(dicpp\) module¶
Digital Image Correlation (dicpp
)¶
Post processing tools for DIC data of cylindrical shells.
Fitting Data (dicpp.fit_data
)¶
Functions to create best-fit imperfection data.
- dicpp.fit_data.best_fit_cylinder(path, H, R_expected=10.0, best_fit_with_fixed_radius=False, sample_size=None, alpha0=0.5, beta0=0.5, x0=None, y0=None, z0=None, z1=None, clip_box=None, R_min=- 1000000.0, R_max=1000000.0, loadtxt_kwargs={}, verbose=True, output_path=None, ls_kwargs={'ftol': None, 'gtol': None, 'jac': '3-point', 'max_nfev': 1000000, 'method': 'trf', 'xtol': 1e-08})¶
Fit a best cylinder for a given set of measured data
The coordinate transformation which must be performed in order to adjust the raw data to the finite element coordinate system is illustrated below:
This transformation can be represented in matrix form as:
\[\begin{split}\left\{\begin{matrix} x_c\\y_c\\z_c \end{matrix}\right\} = [R_z][R_y][R_x] \left\{\begin{matrix} x_i + x_0\\ y_i + y_0\\ z_i + z_0 \end{matrix}\right\} + \left\{\begin{matrix} 0\\ 0\\ z_1 \end{matrix}\right\}\end{split}\]with:
\[\begin{split}[R_x]=\begin{bmatrix} 1 & 0 & 0 \\ 0 & \cos{\alpha} & -\sin{\alpha} \\ 0 & \sin{\alpha} & \cos{\alpha} \end{bmatrix}\end{split}\]\[\begin{split}[R_y]=\begin{bmatrix} \cos{\beta} & 0 & \sin{\beta} \\ 0 & 1 & 0 \\ -\sin{\beta} & 0 & \cos{\beta} \end{bmatrix}\end{split}\]\[\begin{split}[R_z]=\begin{bmatrix} cos{\gamma} & -sin{\gamma} & 0 \\ sin{\gamma} & cos{\gamma} & 0 \\ 0 & 0 & 1 \end{bmatrix}\end{split}\]For the best-fit cylinder one can assume \(\gamma=0\), such that there are six unknown variables:
the three components of the first translation \(\Delta x_0\), \(\Delta y_0\) and \(\Delta z_0\)
the rotation angles \(\alpha\) and \(\beta\)
the second axial translation \(\Delta z_1\)
The six unknowns are found in a non-linear least-squares optimization problem (solved with
scipy.optimize.least_squares
), where the measured data is transformed to the reference coordinate system and then compared with a reference cylinder in order to compute the residual error.Since the measured data may have an unknown radius \(R\), the solution also involves finding this radius.
- Parameters
- pathstr or np.ndarray
The path of the file containing the data. Can be a full path using
r"C:\Temp\inputfile.txt"
, for example. The input file must have 3 columns “\(x\) \(y\) \(z\)” expressed in Cartesian coordinates.This input can also be a
np.ndarray
object, with \(x\), \(y\), \(z\) in each corresponding column.- Hfloat
The nominal height of the cylinder.
- R_expectedfloat, optional
The nominal radius of the cylinder, used as a first guess to find the best-fit radius (
R_best_fit
). Note that, if not specified, more iterations may be required.- best_fit_with_fixed_radiusbool, optional
If
True
, a constant value given byR_expected
is used to calculate the best fit.- sample_sizeNone or int, optional
If the input file containing the measured data is too big it may be convenient to use only a sample of it in order to calculate the best fit.
- alpha0, beta0, x0, y0, z0, z1: float, optional
Initial guess for alpha, beta, x0, y0, z0, z1.
- clip_boxNone or sequence, optional
Clip input points into [xmin, xmax, ymin, ymax, zmin, zmax].
- verbosebool, optional
If
True
activates log messages.- output_pathstr or None, optional
Save
out
to a pickle dump file.- loadtxt_kwargsdict, optional
Keyword arguments passed to
np.loadtxt
.- ls_kwargsdict, optional
Keyword arguments passed to
scipy.optimize.least_squares
.
- Returns
- outdict
A Python dictionary with the entries:
out['R_best_fit']
floatThe best-fit radius of the input sample.
out['x0']
floatTranslation in x
out['y0']
floatTranslation in y
out['z0']
floatTranslation in z
out['z1']
floatSecond translation in z, after rotation
out['alpharad']
np.ndarray\(\alpha\) angle to rotate input data.
out['betarad']
np.ndarray\(\beta\) angle to rotate input data.
out['input_pts']
np.ndarrayThe input points in a \(3 \times N\) 2-D array.
out['clip_box_mask']
np.ndarrayThe mask after applying the clip_box in a \(N\) 1-D boolean array.
out['output_pts']
np.ndarrayThe transformed points in a \(3 \times N\) 2-D array.
Examples
General usage
For a given cylinder with expected radius and height of
R_expected
andH
:from dicpp.fit_data import best_fit_cylinder out = best_fit_cylinder(path, H=H, R_expected=R_expected) R_best_fit = out['R_best_fit']
Using the transformation data
For a given input data with \(x, y, z\) positions in each line:
x, y, z = np.loadtxt('input_file.txt', unpack=True)
the transformation can be obtained with:
Rx = calc_Rx(alpharad) Ry = calc_Ry(betarad) xnew, ynew, znew = (Ry @ Rx).dot(np.vstack((x + x0, y + y0, z + z0))) znew += z1
and the inverse transformation:
znew -= z1 x, y, z = (Rx.T @ Ry.T).dot(np.vstack((xnew, ynew, znew))) - np.array([x0, y0, y0])
- dicpp.fit_data.best_fit_elliptic_cylinder(path, H, a_expected=10.0, b_expected=10.0, best_fit_with_fixed_a=False, alpha0=0.5, beta0=0.5, gamma0=0.0, x0=None, y0=None, z0=None, z1=None, clip_box=None, a_min=- 1000000.0, a_max=1000000.0, b_min=- 1000000.0, b_max=1000000.0, verbose=True, output_path=None, loadtxt_kwargs={}, ls_kwargs={'ftol': None, 'gtol': None, 'jac': '3-point', 'max_nfev': 1000000, 'method': 'trf', 'xtol': 1e-08})¶
Fit a best cylinder for a given set of measured data
The coordinate transformation which must be performed in order to adjust the raw data to the finite element coordinate system is illustrated below:
This transformation can be represented in matrix form as:
\[\begin{split}\left\{\begin{matrix} x_c\\y_c\\z_c \end{matrix}\right\} = [R_z][R_y][R_x] \left\{\begin{matrix} x_i + x_0\\ y_i + y_0\\ z_i + z_0 \end{matrix}\right\} + \left\{\begin{matrix} 0\\ 0\\ z_1 \end{matrix}\right\}\end{split}\]with:
\[\begin{split}[R_x]=\begin{bmatrix} 1 & 0 & 0 \\ 0 & \cos{\alpha} & -\sin{\alpha} \\ 0 & \sin{\alpha} & \cos{\alpha} \end{bmatrix}\end{split}\]\[\begin{split}[R_y]=\begin{bmatrix} \cos{\beta} & 0 & \sin{\beta} \\ 0 & 1 & 0 \\ -\sin{\beta} & 0 & \cos{\beta} \end{bmatrix}\end{split}\]\[\begin{split}[R_z]=\begin{bmatrix} cos{\gamma} & -sin{\gamma} & 0 \\ sin{\gamma} & cos{\gamma} & 0 \\ 0 & 0 & 1 \end{bmatrix}\end{split}\]There are seven unknown variables:
the three components of the first translation \(\Delta x_0\), \(\Delta y_0\) and \(\Delta z_0\)
the rotation angles \(\alpha\), \(\beta\) and \(\gamma\)
the second axial translation \(\Delta z_1\)
The six unknowns are found in a non-linear least-squares optimization problem (solved with
scipy.optimize.least_squares
), where the measured data is transformed to the reference coordinate system and then compared with a reference cylinder in order to compute the residual error.Since the measured data may have an unknown minor and major radii \(a,b\), the solution also involves finding these radii.
- Parameters
- pathstr or np.ndarray
The path of the file containing the data. Can be a full path using
r"C:\Temp\inputfile.txt"
, for example. The input file must have 3 columns “\(x\) \(y\) \(z\)” expressed in Cartesian coordinates.This input can also be a
np.ndarray
object, with \(x\), \(y\), \(z\) in each corresponding column.- Hfloat
The nominal height of the cylinder.
- a_expected, a_min, a_maxfloat, optional
The major radius of the elliptic cylinder, used as a first guess to find the best-fit major radius (
a_best_fit
). Note that if not specified more iterations may be required.- b_expected, b_min, b_maxfloat, optional
The minor radius of the elliptic cylinder, used as a first guess to find the best-fit minor radius (
b_best_fit
). Note that if not specified more iterations may be required.- best_fit_with_fixed_abool, optional
If
True
, a constant value given bya_expected
is used to calculate the best fit.- alpha0, beta0, gamma0, x0, y0 ,z0, z1: float, optional
Initial guess for alpha, beta, gamma, x0, y0, z0, z1.
- clip_boxNone or sequence, optional
Clip input points into [xmin, xmax, ymin, ymax, zmin, zmax]
- verbosebool, optional
If
True
activates log messages.- output_pathstr or None, optional
Save
out
to a pickle dump file.- loadtxt_kwargsdict, optional
Keyword arguments passed to
np.loadtxt
- ls_kwargsdict, optional
Keyword arguments passed to
scipy.optimize.least_squares
.
- Returns
- outdict
A Python dictionary with the entries:
out['a_best_fit']
,out['b_best_fit']
: floatThe best-fit radii of the input sample.
out['x0']
floatTranslation in x
out['y0']
floatTranslation in y
out['z0']
floatTranslation in z
out['alpharad']
np.ndarray\(\alpha\) angle to rotate input data.
out['betarad']
np.ndarray\(\beta\) angle to rotate input data.
out['gammarad']
np.ndarray\(\gamma\) angle to rotate input data.
out['input_pts']
np.ndarrayThe input points in a \(3 \times N\) 2-D array.
out['clip_box_mask']
np.ndarrayThe mask after applying the clip_box in a \(N\) 1-D boolean array.
out['output_pts']
np.ndarrayThe transformed points in a \(3 \times N\) 2-D array.
Examples
General usage
For a given elliptic cylinder with expected radii and height of
a_expected
,b_expected
andH
:from dicpp.fit_data import best_fit_elliptic_cylinder out = best_fit_elliptic_cylinder(path, H=H, a_expected=a_expected, b_expected=b_expected) a_best_fit = out['a_best_fit'] b_best_fit = out['b_best_fit']
Using the transformation data
For a given input data with \(x, y, z\) positions in each line:
x, y, z = np.loadtxt('input_file.txt', unpack=True)
the transformation can be obtained with:
Rx = calc_Rx(alpharad) Ry = calc_Ry(betarad) Rz = calc_Rz(gammarad) xnew, ynew, znew = (Rz @ Ry @ Rx).dot(np.vstack((x + x0, y + y0, z + z0)))
and the inverse transformation:
x, y, z = (Rx.T @ Ry.T @ Rz.T).dot(np.vstack((xnew, ynew, znew))) - np.array([x0, y0, y0])
- dicpp.fit_data.calc_Rx(a)¶
Rotation matrix around X for a 3D vector
- dicpp.fit_data.calc_Ry(b)¶
Rotation matrix around Y for a 3D vector
- dicpp.fit_data.calc_Rz(c)¶
Rotation matrix around Z for a 3D vector
- dicpp.fit_data.calc_c0(path, m0=50, n0=50, funcnum=2, fem_meridian_bot2top=True, rotatedeg=None, filter_m0=None, filter_n0=None, sample_size=None, loadtxt_kwargs={})¶
Find the coefficients that best fit the \(w_0\) imperfection
The measured data will be fit using one of the following functions, selected using the
funcnum
parameter:Half-Sine Function
\[w_0 = \sum_{i=1}^{m_0}{ \sum_{j=0}^{n_0}{ {c_0}_{ij}^a sin{b_z} sin{b_\theta} +{c_0}_{ij}^b sin{b_z} cos{b_\theta} }}\]Half-Cosine Function (default)
\[w_0 = \sum_{i=0}^{m_0}{ \sum_{j=0}^{n_0}{ {c_0}_{ij}^a cos{b_z} sin{b_\theta} +{c_0}_{ij}^b cos{b_z} cos{b_\theta} }}\]Complete Fourier Series
\[w_0 = \sum_{i=0}^{m_0}{ \sum_{j=0}^{n_0}{ {c_0}_{ij}^a sin{b_z} sin{b_\theta} +{c_0}_{ij}^b sin{b_z} cos{b_\theta} +{c_0}_{ij}^c cos{b_z} sin{b_\theta} +{c_0}_{ij}^d cos{b_z} cos{b_\theta} }}\]where:
\[ \begin{align}\begin{aligned}b_z = i \pi \frac z H_{points}\\b_\theta = j \theta\end{aligned}\end{align} \]where \(H_{points}\) represents the difference between the maximum and the minimum \(z\) values in the imperfection file.
The approximation can be written in matrix form as:
\[w_0 = [g] \{c_0\}\]where \([g]\) carries the base functions and \({c_0}\) the respective amplitudes. The solution consists on finding the best \({c_0}\) that minimizes the least-square error between the measured imperfection pattern and the \(w_0\) function.
- Parameters
- pathstr or np.ndarray
The path of the file containing the data. Can be a full path using
r"C:\Temp\inputfile.txt"
, for example. The input file must have 3 columns “\(\theta\) \(z\) \(imp\)” expressed in Cartesian coordinates.This input can also be a
np.ndarray
object, with \(\theta\), \(z\), \(imp\) in each corresponding column.- m0int
Number of terms along the meridian (\(z\)).
- n0int
Number of terms along the circumference (\(\theta\)).
- funcnumint, optional
As explained above, selects the base functions used for the approximation.
- fem_meridian_bot2topbool, optional
A boolean indicating if the finite element has the \(x\) axis starting at the bottom or at the top.
- rotatedegfloat or None, optional
Rotation angle in degrees telling how much the imperfection pattern should be rotated about the \(X_3\) (or \(Z\)) axis.
- filter_m0list, optional
The values of
m0
that should be filtered (seefilter_c0()
).- filter_n0list, optional
The values of
n0
that should be filtered (seefilter_c0()
).- sample_sizeNone or int, optional
An in specifying how many points of the imperfection file should be used. If
None
is used all points file will be used in the computations.- loadtxt_kwargsdict, optional
Keyword arguments passed to
np.loadtxt
- Returns
- outnp.ndarray
A 1-D array with the best-fit coefficients, in the following order:
if funcnum==1: w0 = 0 for j in range(n0): sinjt = sin(j*t) cosjt = cos(j*t) for i in range(1, m0+1): sinix = sin(i*pi*x) row = 2*(i-1) + 2*j*m0 w0 += c0[row+0]*sinix*sinjt w0 += c0[row+1]*sinix*cosjt elif funcnum==2: w0 = 0 for j in range(n0): sinjt = sin(j*t) cosjt = cos(j*t) for i in range(m0): cosix = cos(i*pi*x) row = 2*i + 2*j*m0 w0 += c0[row+0]*cosix*sinjt w0 += c0[row+1]*cosix*cosjt elif funcnum==3: w0 = 0 for j in range(n0): sinjt = sin(j*t) cosjt = cos(j*t) for i in range(m0): sinix = sin(i*pi*x) cosix = cos(i*pi*x) row = 4*i + 4*j*m0 w0 += c0[row+0]*sinix*sinjt w0 += c0[row+1]*sinix*cosjt w0 += c0[row+2]*cosix*sinjt w0 += c0[row+3]*cosix*cosjt
Notes
If a similar imperfection pattern is expected along the meridian and along the circumference, the analyst can use an optimized relation between
m0
andn0
in order to achieve a higher accuracy for a given computational cost, as proposed by Castro et al. (2014):\[n_0 = m_0 \frac{\pi(R_{bot}+R_{top})}{2H}\]
- dicpp.fit_data.fa(m0, n0, zs_norm, thetas, funcnum=2)¶
Calculates the matrix with the base functions for \(w_0\)
The calculated matrix is directly used to calculate the \(w_0\) displacement field, when the corresponding coefficients \(c_0\) are known, through:
a = fa(m0, n0, zs_norm, thetas, funcnum) w0 = a.dot(c0)
- Parameters
- m0int
The number of terms along the meridian.
- n0int
The number of terms along the circumference.
- zs_normnp.ndarray
The normalized \(z\) coordinates (from 0. to 1.) used to compute the base functions.
- thetasnp.ndarray
The angles in radians representing the circumferential positions.
- funcnumint, optional
The function used for the approximation (see function
calc_c0()
)
- dicpp.fit_data.filter_c0(m0, n0, c0, filter_m0, filter_n0, funcnum=2)¶
Apply filter to the imperfection coefficients \(\{c_0\}\)
A filter consists on removing some frequencies that are known to be related to rigid body modes or spurious measurement noise. The frequencies to be removed should be passed through inputs
filter_m0
andfilter_n0
.- Parameters
- m0int
The number of terms along the meridian.
- n0int
The number of terms along the circumference.
- c0np.ndarray
The coefficients of the imperfection pattern.
- filter_m0list
The values of
m0
that should be filtered.- filter_n0list
The values of
n0
that should be filtered.- funcnumint, optional
The function used for the approximation (see function
calc_c0()
)
- Returns
- c0_filterednp.ndarray
The filtered coefficients of the imperfection pattern.
- dicpp.fit_data.fw0(m0, n0, c0, xs_norm, ts, funcnum=2)¶
Calculates the imperfection field \(w_0\) for a given input
- Parameters
- m0int
The number of terms along the meridian.
- n0int
The number of terms along the circumference.
- c0np.ndarray
The coefficients of the imperfection pattern.
- xs_normnp.ndarray
The meridian coordinate (\(x\)) normalized to be between
0.
and1.
.- tsnp.ndarray
The angles in radians representing the circumferential coordinate (\(\theta\)).
- funcnumint, optional
The function used for the approximation (see function
calc_c0()
)
- Returns
- w0snp.ndarray
An array with the same shape of
xs_norm
containing the calculated imperfections.
Notes
The inputs
xs_norm
andts
must be of the same size.The inputs must satisfy
c0.shape[0] == size*m0*n0
, where:size=2
iffuncnum==1 or funcnum==2
size=4
iffuncnum==3
Interpolate (dicpp.interpolate
)¶
This module includes some interpolation utilities.
- dicpp.interpolate.interp_theta_z_imp(data, destiny, alphadeg, H_measured, H_model, R_bottom, stretch_H=False, z_offset_bot=None, rotatedeg=0.0, ncp=5, power_parameter=2, ignore_bot_h=None, ignore_top_h=None)¶
Interpolates a data set in the \(\theta, z, imp\) format
This function uses the inverse-weighted algorithm (
inv_weighted()
).- Parameters
- datanumpy.ndarray, shape (N, 3)
The data or an array containing the imperfection file in the \((\theta, Z, imp)\) format.
- destinynumpy.ndarray, shape (M, 3)
The destiny coordinates \((x, y, z)\) where the values will be interpolated to.
- alphadegfloat
The cone semi-vertex angle in degrees.
- H_measuredfloat
The total height of the measured test specimen, including eventual resin rings at the edges.
- H_modelfloat
The total height of the new model, including eventual resin rings at the edges.
- R_bottomfloat
The radius of the model taken at the bottom edge.
- stretch_Hbool, optional
Tells if the height of the measured points, which is usually smaller than the height of the test specimen, should be stretched to fill the whole test specimen. If not, the points will be placed in the middle or using the offset given by
z_offset_bot
and the area not covered by the measured points will be interpolated using the closest available points (the imperfection pattern will look like there was an extrusion close to the edges).- z_offset_botfloat, optional
The offset that should be used from the bottom of the measured points to the bottom of the test specimen.
- rotatedegfloat, optional
Rotation angle in degrees telling how much the imperfection pattern should be rotated about the \(X_3\) (or \(Z\)) axis.
- ncpint, optional
Number of closest points used in the inverse-weighted interpolation.
- power_parameterfloat, optional
Power of inverse weighted interpolation function.
- ignore_bot_hNone or float, optional
Nodes close to the bottom edge are ignored according to this meridional distance.
- ignore_top_hNone or float, optional
Nodes close to the top edge are ignored according to this meridional distance.
- Returns
- ansnumpy.ndarray
An array with M elements containing the interpolated values.
- dicpp.interpolate.inv_weighted(values, source, destiny, ncp=5, power_parameter=2)¶
Interpolates the values taken at one group of points into another using an inverse-weighted algorithm
In the inverse-weighted algorithm a number of \(n_{CP}\) measured points of the input parameter
source
that are closest to a given node in the input parameterdestiny
are found and thevalues
are interpolated as follows:\[{w_0}_{node} = \frac{\sum_{i}^{n_{CP}}{{w_0}_i\frac{1}{w_i}}} {\sum_{i}^{n_{CP}}{\frac{1}{w_i}}}\]where \(w_i\) is the inverse weight of each measured point, calculated as:
\[w_i = \left[(x_{node}-x_i)^2+(y_{node}-y_i)^2+(z_{node}-z_i)^2 \right]^p\]with \(p\) being a power parameter that when increased will increase the relative influence of a closest point.
- Parameters
- values: np.ndarray, shape (N)
Values to be interpolated
- sourcenumpy.ndarray, shape (N, ndim)
Source points corresponding to “values”.
- destinynumpy.ndarray, shape (M, ndim)
The new coordinates where the values will be interpolated to.
- ncpint, optional
Number of closest points used in the inverse-weighted interpolation.
- power_parameterfloat, optional
Power of inverse weighted interpolation function.
- Returns
- dist, data_newnumpy.ndarray, numpy.ndarray
Two 1-D arrays with the distances and interpolated values. The size of this array is
destiny.shape[0]
.
Fitting Data (dicpp.stitch
)¶
Functions to stitch best-fit imperfection data.
- dicpp.stitch.stitch(bf1, bf2, pos_deg1, pos_deg2, height_ref, probe_deg, probe_R, probe_zarray, probe_dist_deg=10, init_deg_min=- 10, init_deg_max=10, init_deg_step=11, init_z_min=- 20, init_z_max=20, init_z_step=11, opt_var_z_min=- 10.0, opt_var_z_max=10.0, opt_var_deg_min=- 10.0, opt_var_deg_max=10.0, error_dr1dr2=None, inv_weighted_kwargs={'ncp': 10, 'power_parameter': 1.7}, ls_kwargs={'diff_step': 0.05, 'ftol': None, 'gtol': None, 'jac': '3-point', 'max_nfev': 1000, 'method': 'trf', 'xtol': 0.0001})¶
Calculate \(\Delta \theta_2\) and \(\Delta z_2\) that stitches the second best fit results to the first.
bf2
must be afterbf1
in a positive rotation aboutz
.- Parameters
- bf1, bf2dict
Results returned by the best fit functions. If the results are for a best-fit elliptic cylinder, the imperfections are remapped to a cylinder of radius``R=out[‘a_best_fit’]``.
- pos_deg1, pos_deg_2float
Nominal circumferential position in degrees of the adjacent best fit results.
pos_deg1
must be such that the \(x\) axis becomes normal tobf1
.- height_reffloat
Reference height.
- probe_degfloat
Position in degrees of the probing line.
- probe_Rfloat
Radial position of the probing line.
- probe_zarray: array-like
Array with the z positions of the probing line.
- probe_dist_degfloat, optional
Angular distance in degrees from the probing line. Points beyond this distance are trimmed out.
- init_deg_min, init_deg_max, init_z_min, init_z_maxfloat, optional
Interval to search for the best initial point for the stitching optimization.
- init_deg_step, init_z_stepint, optional
Number of points in the interval to search for the best initial point for the stitching optimization.
- opt_var_z_min, opt_var_z_max, opt_var_deg_min, opt_var_deg_maxfloat, optional
Variation from initial point to be used in the stitching optimization.
- error_dr1dr2function, optional
Error function in the format
f(dr1, dr2)
used in the stitching optimization. The default functions is(dr1 - dr2)**2
.- inv_weighted_kwargsdict, optional
Keyword arguments paased to
inv_weighted()
.- ls_kwargsdict, optional
Keyword arguments passed to
scipy.optimize.least_squares
.
- Returns
- outdict
A Python dictionary with the entries:
out['delta_deg']
floatThe \(\Delta \theta_2\) offset that stitches the second best-fit result to the first.
out['delta_z']
floatThe \(\Delta z_2\) offset that stitches the second best-fit result to the first.
out['dr1']
array-likeRadial imperfection of the first best-fit result at the probing line.
out['dr2']
array-likeRadial imperfection of the second best-fit result at the probing line.
License¶
Copyright (c) 2021-, Saullo G. P. Castro (S.G.P.Castro@tudelft.nl)
All rights reserved.
Redistribution and use in source and binary forms, with or without modification,
are permitted provided that the following conditions are met:
Redistributions of source code must retain the above copyright notice, this
list of conditions and the following disclaimer.
Redistributions in binary form must reproduce the above copyright notice, this
list of conditions and the following disclaimer in the documentation and/or
other materials provided
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR
ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
(INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
(INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.