Solves least squares problem for QR-decomposed matrix.
Case 1: Matrix - vector operation
IppStatus ippmQRBackSubst_mv_32f(const Ipp32f* pSrc1, int src1Stride1, int src1Stride2, Ipp32f* pBuffer, const Ipp32f* pSrc2, int src2Stride2, Ipp32f* pDst, int dstStride2, int width, int height);
IppStatus ippmQRBackSubst_mv_64f(const Ipp64f* pSrc1, int src1Stride1, int src1Stride2, Ipp64f* pBuffer, const Ipp64f* pSrc2, int src2Stride2, Ipp64f* pDst, int dstStride2, int width, int height);
IppStatus ippmQRBackSubst_mv_32f_P(const Ipp32f** ppSrc1, int src1RoiShift, Ipp32f* pBuffer, const Ipp32f** ppSrc2, int src2RoiShift, Ipp32f** ppDst, int dstRoiShift, int width, int height);
IppStatus ippmQRBackSubst_mv_64f_P(const Ipp64f** ppSrc1, int src1RoiShift, Ipp64f* pBuffer, const Ipp64f** ppSrc2, int src2RoiShift, Ipp64f** ppDst, int dstRoiShift, int width, int height);
Case 2: Matrix - vector array operation
IppStatus ippmQRBackSubst_mva_32f(const Ipp32f* pSrc1, int src1Stride0, int src1Stride1, int src1Stride2, Ipp32f* pBuffer, const Ipp32f* pSrc2, int src2Stride0, int src2Stride2, Ipp32f* pDst, int dstStride0, int dstStride2, int width, int height, int count);
IppStatus ippmQRBackSubst_mva_64f(const Ipp64f* pSrc1, int src1Stride0, int src1Stride1, int src1Stride2, Ipp64f* pBuffer, const Ipp64f* pSrc2, int src2Stride0, int src2Stride2, Ipp64f* pDst, int dstStride0, int dstStride2, int width, int height, int count);
IppStatus ippmQRBackSubst_mva_32f_P(const Ipp32f** ppSrc1, int src1RoiShift, Ipp32f* pBuffer, const Ipp32f** ppSrc2, int src2RoiShift, int src2Stride0, Ipp32f** ppDst, int dstRoiShift, int dstStride0, int width, int height, int count);
IppStatus ippmQRBackSubst_mva_64f_P(const Ipp64f** ppSrc1, int src1RoiShift, Ipp64f* pBuffer, const Ipp64f** ppSrc2, int src2RoiShift, int src2Stride0, Ipp64f** ppDst, int dstRoiShift, int dstStride0, int width, int height, int count);
IppStatus ippmQRBackSubst_mva_32f_L(const Ipp32f** ppSrc1, int src1RoiShift, int src1Stride1, int src1Stride2, Ipp32f* pBuffer, const Ipp32f** ppSrc2, int src2RoiShift, int src2Stride2, Ipp32f** ppDst, int dstRoiShift, int dstStride2, int width, int height, int count);
IppStatus ippmQRBackSubst_mva_64f_L(const Ipp64f** ppSrc1, int src1RoiShift, int src1Stride1, int src1Stride2, Ipp64f* pBuffer, const Ipp64f** ppSrc2, int src2RoiShift, int src2Stride2, Ipp64f** ppDst, int dstRoiShift, int dstStride2, int width, int height, int count);
Case 3: Matrix array - vector array operation
IppStatus ippmQRBackSubst_mava_32f(const Ipp32f* pSrc1, int src1Stride0, int src1Stride1, int src1Stride2, Ipp32f* pBuffer, const Ipp32f* pSrc2, int src2Stride0, int src2Stride2, Ipp32f* pDst, int dstStride0, int dstStride2, int width, int height, int count);
IppStatus ippmQRBackSubst_mava_64f(const Ipp64f* pSrc1, int src1Stride0, int src1Stride1, int src1Stride2, Ipp64f* pBuffer, const Ipp64f* pSrc2, int src2Stride0, int src2Stride2, Ipp64f* pDst, int dstStride0, int dstStride2, int width, int height, int count);
IppStatus ippmQRBackSubst_mava_32f_P(const Ipp32f** ppSrc1, int src1RoiShift, int src1Stride0, Ipp32f* pBuffer, const Ipp32f** ppSrc2, int src2RoiShift, int src2Stride0, Ipp32f** ppDst, int dstRoiShift, int dstStride0, int width, int height, int count);
IppStatus ippmQRBackSubst_mava_64f_P(const Ipp64f** ppSrc1, int src1RoiShift, int src1Stride0, Ipp64f* pBuffer, const Ipp64f** ppSrc2, int src2RoiShift, int src2Stride0, Ipp64f** ppDst, int dstRoiShift, int dstStride0, int width, int height, int count);
IppStatus ippmQRBackSubst_mava_32f_L(const Ipp32f** ppSrc1, int src1RoiShift, int src1Stride1, int src1Stride2, Ipp32f* pBuffer, const Ipp32f** ppSrc2, int src2RoiShift, int src2Stride2, Ipp32f** ppDst, int dstRoiShift, int dstStride2, int width, int height, int count);
IppStatus ippmQRBackSubst_mava_64f_L(const Ipp64f** ppSrc1, int src1RoiShift, int src1Stride1, int src1Stride2, Ipp64f* pBuffer, const Ipp64f** ppSrc2, int src2RoiShift, int src2Stride2, Ipp64f** ppDst, int dstRoiShift, int dstStride2, int width, int height, int count);
pSrc1, ppSrc1 |
Pointer to the source matrix or array of matrices. The matrix(ces) must be a result of calling ippmQRDecomp . |
src1Stride0 |
Stride between the matrices in the source matrix array. |
src1Stride1 |
Stride between the rows in the source matrix(ces). |
src1Stride2 |
Stride between the elements in the source matrix(ces). |
src1RoiShift |
ROI shift in the source matrix(ces). |
pSrc2, ppSrc2 |
Pointer to the source vector or array of vectors. |
src2Stride0 |
Stride between the vectors in the source vector array. |
src2Stride2 |
Stride between the elements in the source vector(s). |
src2RoiShift |
ROI shift in the source vector(s). |
pBuffer |
Pointer to a pre-allocated auxiliary array to be used for internal computations. The number of elements in the array must be at least height. |
pDst, ppDst |
Pointer to the destination matrix or array of matrices. |
dstStride0 |
Stride between the vectors in the destination array. |
dstStride2 |
Stride between the elements of the destination vector(s). |
dstRoiShift |
ROI shift in the destination vector(s). |
width |
Matrix width and destination vector length. |
height |
Matrix height and source vector length. |
count |
Number of matrices and right-hand part vectors in the array. |
The function ippmQRBackSubst is declared in the ippm.h header file. The function solves the least squares problem for an overdetermined system of linear equations
Ax = b,
where
A is the matrix of linear equations system, stored in pSrc1 or ppSrc1,
b is the vector of the right-hand side, stored in pSrc2 or ppSrc2,
x is the unknown vector, stored in pDst or ppDst.
The number of equations in the system height is equal to or greater than the number of unknown variables width, rank (A) = width.
A typical least squares problem is as follows: given a matrix A and a vector b, find the vector x that minimizes the L2-norm
You should call the function ippmQRDecomp to compute the QR decomposition of A before calling the ippmQRBackSubst function (Please see description for the function ippmQRDecomp.)
The following example demonstrates how to use the functions ippmQRDecomp_m_32f and ippmQRBackSubst_mva_32f. For more information, see also examples in Getting Started.
IppStatus QRFactorization_32f(void) {
/* Source matrix with width=4 and height=5 */
Ipp32f pSrc[5*4] = {1, 1, 1, 1,
1, 3, 1, 1,
1,-1, 3, 1,
1, 1, 1, 3,
1, 1, 1, -1 };
int srcStride2 = sizeof(Ipp32f);
int srcStride1 = 4*sizeof(Ipp32f);
/* Solver right-part is 2 vectors with length=5 */
Ipp32f pSrc2[2*5] = { 2, 1, 6, 3, 1,
3, 4, 5, 6, 1 };
int src2Stride2 = sizeof(Ipp32f);
int src2Stride0 = 5*sizeof(Ipp32f);
Ipp32f pDecomp[5*4]; /* Decomposed matrix location */
int decompStride2 = sizeof(Ipp32f);
int decompStride1 = 4*sizeof(Ipp32f);
Ipp32f pDst[2*4]; /* Solver destination location */
int dstStride2 = sizeof(Ipp32f);
int dstStride0 = 4*sizeof(Ipp32f);
int width = 4;
int height = 5;
int count = 2;
Ipp32f pBuffer[5]; /* Buffer location */
IppStatus status = ippmQRDecomp_m_32f((const Ipp32f*)pSrc,
srcStride1, srcStride2, pBuffer,
pDecomp, decompStride1, decompStride2, width, height);
status = ippmQRBackSubst_mva_32f((const Ipp32f*)pDecomp,
decompStride1, decompStride2, pBuffer, pSrc2, src2Stride0,
src2Stride2, pDst, dstStride0, dstStride2, width, height, count);
/*
// It is required for QR decomposition function to check return status
// for catching wrong result in case of invalid input data
*/
if (status == ippStsNoErr) {
status = ippmQRBackSubst_mva_32f((const Ipp32f*)pDecomp,
decompStride1, decompStride2, pBuffer, pSrc2, src2Stride0,
src2Stride2, pDst, dstStride0, dstStride2, width, height, count);
printf_m_Ipp32f("QRDecomp result:", pDecomp, 4, 5, status);
printf_va_Ipp32f("2 destination vectors:", pDst, 4, 2, status);
} else {
printf("Function returns status: %s \n", ippGetStatusString(status));
}
return status;
}
The program above produces the following output:
QRDecomp result:
-2.236068 -2.236068 -3.130495 -2.236068
0.309017 -2.828427 1.414214 -0.000000
0.309017 -0.414214 -1.095445 0.000000
0.309017 -0.000000 -0.130449 -2.828427
0.309017 -0.000000 -0.130449 -0.414214
2 destination vectors:
0.500000 -0.500000 1.500000 0.500000
0.583333 0.333333 1.166667 1.250000
ippStsOk |
Returns no error. |
ippStsNullPtrErr |
Returns an error when one of the input pointers is NULL. |
ippStsSizeErr |
Returns an error when the size of the source matrix is equal to 0. |
ippStsStrideMatrixErr |
Returns an error when the stride value is not positive or not divisible by size of data type. |
ippStsRoiShiftMatrixErr |
RoiShift value is negative or not divisible by size of data type. |
ippStsCountMatrixErr |
Returns an error when the count value is less or equal to zero. |
ippStsSizeMatchMatrixErr |
Returns an error when the sizes of the source matrices are unsuitable. |
Copyright © 2000 - 2010, Intel Corporation. All rights reserved.