?geqrf

计算矩阵的QR分解,即

接口定义

C Interface:

void dgeqrf_(const int *m, const int *n, double *a, const int *lda, double *tau, double *work, const int *lwork, int *info);

void sgeqrf_(const int *m, const int *n, float *a, const int *lda, float *tau, float *work, const int *lwork, int *info);

void cgeqrf_(const int *m, const int *n, float _Complex *a, const int *lda, float _Complex *tau, float _Complex *work, const int *lwork, int *info);

void zgeqrf_(const int *m, const int *n, double _Complex *a, const int *lda, double _Complex *tau, double _Complex *work, const int *lwork, int *info);

Fortran Interface:

DGEQRF(m, n, a, lda, tau, work, lwork, info);

SGEQRF(m, n, a, lda, tau, work, lwork, info);

CGEQRF(m, n, a, lda, tau, work, lwork, info);

ZGEQRF(m, n, a, lda, tau, work, lwork, info);

参数

参数名

类型

描述

输入/输出

m

整型数

矩阵A的行数。

输入

n

整型数

矩阵A的列数。

输入

a

  • 在sgeqrf中为单精度浮点型数组
  • 在dgeqrf中为双精度浮点型数组
  • 在cgeqrf中为单精度复数型数组
  • 在zgeqrf中为双精度复数型数组
  • 调用前保存待分解的矩阵A。
  • 调用后对角线以及对角线以上的部分保存大小为min(m,n)*n的矩阵R(当m>=n时R为上三角矩阵),对角线以下元素和tau共同表示正交矩阵Q(参见说明)。

输入/输出

lda

整型数

A的leading dimension大小,要求lda>=max(1, m)。

输出

tau

  • 在sgeqrf中为单精度浮点型数组
  • 在dgeqrf中为双精度浮点型数组
  • 在cgeqrf中为单精度复数型数组
  • 在zgeqrf中为双精度复数型数组

初等反射的系数,长度为min(m,n)(参见说明)。

输出

work

  • 在sgeqrf中为单精度浮点型数组
  • 在dgeqrf中为双精度浮点型数组
  • 在cgeqrf中为单精度复数型数组
  • 在zgeqrf中为双精度复数型数组

临时存储空间,使用lwork=-1调用后work[0]为最优的lwork值。

输出

lwork

整数型

work数组的长度。

lwork=-1时查询最优work大小,结果保存在work[0]中,否则要求lwork>=n。

输入

info

整数型

执行结果:

  • 等于0:成功。
  • 小于0:第-info个参数值不合法。

输出

分解结果矩阵Q通过一系列初等反射乘积表示:Q=H(1)*H(2)*...*H(k), k=min(m,n)。H(i)=I-tau*v*v'。tau为标量,v为向量且前i-1个元素为0,第i个元素为1,剩余元素保存在a的第i列中(a的下三角部分)。

依赖

#include "klapack.h"

示例

C Interface:

    int m = 6; 
    int n = 4; 
    int lda = 6; 
    int info = 0; 
    double tau[4]; 
    double *work = NULL; 
    double qwork; 
    int lwork = -1; 
    /* 
     * A (6x4, stored in column-major): 
     *  2.229  1.273  0.087  0.035 
     *  8.667  4.317  4.091  3.609 
     *  0.205  7.810  2.553  6.507 
     *  2.758  2.911  8.791  5.051 
     *  8.103  1.396  1.317  4.738 
     *  8.859  3.161  0.808  5.972 
     */ 
    double a[] = {2.229, 8.667, 0.205, 2.758, 8.103, 8.859, 
                    1.273, 4.317, 7.810, 2.911, 1.396, 3.161, 
                    0.087, 4.091, 2.553, 8.791, 1.317, 0.808, 
                    0.035, 3.609, 6.507, 5.051, 4.738, 5.972}; 
    /* Query optimal work size */ 
    dgeqrf_(&m, &n, a, &lda, tau, &qwork, &lwork, &info); 
    if (info != 0) { 
        return ERROR; 
    } 
    lwork = (int)qwork; 
    work = (double *)malloc(sizeof(double) * lwork); 
    /* Calculate QR */ 
    dgeqrf_(&m, &n, a, &lda, tau, work, &lwork, &info); 
    free(work); 
    /* 
     * Output: 
     * tau 
     *   1.1464 1.0946 1.3670 1.7415 
     * A output (stored in column-major) 
     *  -15.2274 -5.8577  -5.1387  -9.0573 
     *   0.4965  -8.2070  -4.5806  -6.5291 
     *   0.0117   0.8600   7.4527   1.3953 
     *   0.1580   0.1986  -0.6634  -3.5188 
     *   0.4642  -0.2130  -0.0183   0.2805 
     *   0.5075  -0.0510   0.1506   0.2642 
     */

Fortran Interface:

        PARAMETER (m = 6) 
        PARAMETER (n = 4) 
        PARAMETER (lda = 6) 
        INTEGER :: info = 0 
        REAL(8) :: tau(4) 
        REAL(8) :: qwork(1) 
        INTEGER :: lwork = -1 
        REAL(8), ALLOCATABLE :: work(:) 
* 
*       A (6x4, stored in column-major): 
*         2.229  1.273  0.087  0.035 
*         8.667  4.317  4.091  3.609 
*         0.205  7.810  2.553  6.507 
*         2.758  2.911  8.791  5.051 
*         8.103  1.396  1.317  4.738 
*         8.859  3.161  0.808  5.972 
* 
        REAL(8) :: a(m, n) 
        DATA a / 2.229, 8.667, 0.205, 2.758, 8.103, 8.859, 
     $           1.273, 4.317, 7.810, 2.911, 1.396, 3.161, 
     $           0.087, 4.091, 2.553, 8.791, 1.317, 0.808, 
     $           0.035, 3.609, 6.507, 5.051, 4.738, 5.972 / 
 
        EXTERNAL DGEQRF 
*       Query optimal work size 
        CALL DGEQRF(m, n, a, lda, tau, qwork, lwork, info) 
        IF (info.NE.0) THEN 
            CALL EXIT(1) 
        END IF 
        lwork = INT(qwork(1)) 
        ALLOCATE(work(lwork)) 
*       Calculate QR 
        CALL DGEQRF(m, n, a, lda, tau, work, lwork, info) 
        DEALLOCATE(work) 
 
*       Output: 
*       tau 
*         1.1464 1.0946 1.3670 1.7415 
*       A output (stored in column-major) 
*        -15.2274 -5.8577  -5.1387  -9.0573 
*         0.4965  -8.2070  -4.5806  -6.5291 
*         0.0117   0.8600   7.4527   1.3953 
*         0.1580   0.1986  -0.6634  -3.5188 
*         0.4642  -0.2130  -0.0183   0.2805 
*         0.5075  -0.0510   0.1506   0.2642