F06 Chapter Contents
F06 Chapter Introduction (PDF version)
NAG Library Manual

NAG Library Chapter Introduction

F06 – Linear Algebra Support Routines

+ Contents

1  Scope of the Chapter

This chapter is concerned with basic linear algebra routines which perform elementary algebraic operations involving scalars, vectors and matrices. It includes routines which conform to the specifications of the BLAS (Basic Linear Algebra Subprograms).

2  Background to the Problems

A number of the routines in this chapter meet the specification of the Basic Linear Algebra Subprograms (BLAS) as described in Lawson et al. (1979), Dodson et al. (1991), Dongarra et al. (1988) and Dongarra et al. (1990). The first reference describes a set of routines concerned with operations on scalars and vectors: these will be referred to here as the Level-0 and the Level-1 BLAS; the second reference describes a set of routines concerned with operations on sparse vectors: these will be referred to here as the Level-1 Sparse BLAS; the third reference describes a set of routines concerned with matrix-vector operations: these will be referred to here as the Level-2 BLAS; and the fourth reference describes a set of routines concerned with matrix-matrix operations: these will be referred to here as the Level-3 BLAS.
More generally we refer to the scalar routines in the chapter as Level-0 routines, to the vector routines as Level-1 routines, to the matrix-vector and matrix routines as Level-2 routines, and to the matrix-matrix routines as Level-3 routines. The terminology reflects the number of operations involved. For example, a Level-2 routine involves On2 operations for an n×n matrix.

2.1  The Use of BLAS Names

Many of the routines in other chapters of the Library call the routines in this chapter, and in particular a number of the BLAS are called. These routines are usually called by the BLAS name and so, for correct operation of the Library, it is essential that you do not attempt to link your own versions of these routines. If you are in any doubt about how to avoid this, please consult your computer centre or the NAG Response Centre.
The BLAS names are used in order to make use of efficient implementations of the routines when these exist. Such implementations are stringently tested before being used, to ensure that they correctly meet the specification of the BLAS, and that they return the desired accuracy (see, for example, Dodson et al. (1991), Dongarra et al. (1988) and Dongarra et al. (1990)).

2.2  Background Information

Most of the routines in this chapter implement straightforward scalar, vector and matrix operations that need no further explanation beyond a statement of the purpose of the routine. In this section we give some additional background information to those few cases where additional explanation may be necessary. A sub-section is devoted to each topic.

2.2.1  Real plane rotations

There are a number of routines in the chapter concerned with setting up and applying plane rotations. This section discusses the real case and the next section looks at the complex case. For further background information see Golub and Van Loan (1996).
A plane rotation matrix for the i,j plane, Rij, is an orthogonal matrix that is different from the unit matrix only in the elements rii, rjj, rij and rji. If we put
R= rii rij rji rjj , (1)
then, in the real case, it is usual to choose Rij so that
R= -c s -s c ,   c=cosθ,   s=sinθ.
An exception is routine F06FPF which applies the so-called symmetric rotation for which
R= c -s s -c . (2)
The application of plane rotations is straightforward and needs no further elaboration, so further comment is made only on the construction of plane rotations.
The most common use of plane rotations is to choose c and s so that for given a and b,
-c s -s c a b = d 0 . (3)
In such an application the matrix R is often termed a Givens rotation matrix. There are two approaches to the construction of real Givens rotations in Chapter F06.
The BLAS routine F06AAF (DROTG), see Lawson et al. (1979) and Dodson and Grimes (1982), computes c, s and d as
d=σa2+b21/2,
c= a/d, d0, 1, d=0,     s= b/d, d0, 0, d=0, (4)
where σ= signa, a>b signb, ab .
The value z defined as
z= s, s<c  or  c=0 1/c, 0<cs (5)
is also computed and this enables c and s to be reconstructed from the single value z as
c= 0, z=1 1-z21/2, z<1 1/z, z>1     s= 1, z=1 z, z<1 1-c21/2, z>1 .
The other Chapter F06 routines for constructing Givens rotations are based on the computation of the tangent, t=tanθ. t is computed as
t= 0, b=0 b/a, ba.flmax,b0 signb/a.flmax, b>a.flmax signb.flmax, b0,a=0 (6)
where flmax=1/flmin and flmin is the small positive value returned by X02AMF. The values of c and s are then computed or reconstructed via t as
c= 1/1+t21/2, epst1/eps 1, t<eps 1/t, t>1/eps     s= c.t, epst1/eps t, t<eps signt, t>1/eps (7)
where eps is the machine precision. Note that c is always non-negative in this scheme and that the same expressions are used in the initial computation of c and s from a and b as in any subsequent recovery of c and s via t. This is the approach used by many of the NAG Library routines that require plane rotations. d is computed simply as
d=c.a+s.b.
You need not be too concerned with the above detail, since routines are provided for setting up, recovering and applying such rotations.
Another use of plane rotations is to choose c and s so that for given x, y and z 
-c s -s c x y y z c -s s -c = a 0 0 b . (8)
In such an application the matrix R is often termed a Jacobi rotation matrix. The routine that generates a Jacobi rotation (F06BEF) first computes the tangent t and then computes c and s via t as described above for the Givens rotation.

2.2.2  Complex plane rotations

In the complex case a plane rotation matrix for the i,j plane, Rij is a unitary matrix and, analogously to the real case, it is usual to choose Rij so that
R= -c- s- -s c ,   c2+s2=1, (9)
where a- denotes the complex conjugate of a.
The BLAS (see Lawson et al. (1979)) do not contain a routine for the generation of complex rotations, and so the routines in Chapter F06 are all based upon computing c and s via t=b/a in an analogous manner to the real case. R can be chosen to have either c real, or s real and there are routines for both cases.
When c is real then it is non-negative and the transformation
-c s- -s c a b = d 0 (10)
is such that if a is real then d is also real.
When s is real then the transformation
-c- s -s c a b = d 0 (11)
is such that if b is real then d is also real.

2.2.3  Elementary real (Householder) reflections

There are a number of routines in the chapter concerned with setting up and applying Householder transformations. This section discusses the real case and the next section looks at the complex case. For further background information see Golub and Van Loan (1996).
A real elementary reflector, P, is a matrix of the form
P=I-μuuT,   μuTu=2, (12)
where μ is a scalar and u is a vector, and P is both symmetric and orthogonal. In the routines in Chapter F06, u is expressed in the form
u= ζ z ,   ζ​ a scalar (13)
because in many applications ζ and z are not contiguous elements. The usual use of elementary reflectors is to choose μ and u so that for given α and x 
P α x = β 0 ,   α​ and ​β​ scalars. (14)
Such a transformation is often termed a Householder transformation. There are two choices of μ and u available in Chapter F06.
The first form of the Householder transformation is compatible with that used by LINPACK (see Dongarra et al. (1979)) and has
μ=1/ζ. (15)
This choice makes ζ satisfy
1ζ2.
The second form, and the form used by many of the NAG Library routines, has
μ=1 (16)
which makes
1ζ2.
In both cases the special setting
ζ=0 (17)
is used by the routines to flag the case where P=I.
Note that while there are routines to apply an elementary reflector to a vector, there are no routines available in Chapter F06 to apply an elementary reflector to a matrix. This is because such transformations can readily and efficiently be achieved by calls to the matrix-vector Level 2 BLAS routines. For example, to form PA for a given matrix
PA = I-μuuTA=A-μuuTA = A-μubT,  b=ATu, (18)
and so we can call a matrix-vector product routine to form b=ATu and then call a rank-one update routine to form A-μubT. Of course, we must skip the transformation when ζ has been set to zero.

2.2.4  Elementary complex (Householder) reflections

A complex elementary reflector, P, is a matrix of the form
P=I-μuuH,   μuHu=2,   μ​ real,
where uH denotes the complex conjugate of uT, and P is both Hermitian and unitary. For convenience in a number of applications this definition can be generalized slightly by allowing μ to be complex and so defining the generalized elementary reflector as
P=I-μuuH,   μ2uHu=μ+μ- (19)
for which P is still unitary, but is no longer Hermitian.
The Chapter F06 routines choose μ and ζ so that
Reμ=1,   Imζ=0 (20)
and this reduces to (12) with the choice (16) when μ and u are real. This choice is used because μ and u can now be chosen so that in the Householder transformation (14) we can make
Imβ=0
and, as in the real case,
1ζ2.
Rather than returning μ and ζ as separate parameters the Chapter F06 routines return the single complex value θ defined as
θ=ζ+i.Imμ,  i=-1.
Obviously ζ and μ can be recovered as
ζ=Reθ,  μ=1+i.Imθ.
The special setting
θ=0
is used to flag the case where P=I, and
Reθ0,  Imθ0
is used to flag the case where
P= γ 0 0 I ,   γ​ a scalar (21)
and in this case θ actually contains the value of γ. Notice that with both (18) and (21) we merely have to supply θ- rather than θ in order to represent PH.

3  Recommendations on Choice and Use of Available Routines

3.1  Naming Scheme

3.1.1  NAG names

Table 1 shows the naming scheme for the routines in this chapter.
    Level-0 Level-1 Level-2 Level-3
integer Chapter F06 routine F06D_F
‘real’ BLAS routine F06A_F F06E_F F06P_F F06Y_F
‘real’ Chapter F06 routine F06B_F F06F_F F06Q_F
        F06R_F  
‘complex’ BLAS routine F06G_F F06S_F F06Z_F
‘complex’ Chapter F06 routine F06C_F F06H_F F06T_F
        F06U_F  
‘mixed type’ BLAS routine F06J_F
‘mixed type’ Chapter F06 routine F06K_F F06V_F
‘real’ and ‘complex’ LAPACK routines F06W_F F06W_F
Table 1
The heading ‘mixed type’ is for routines where a mixture of data types is involved, such as a routine that returns the real Euclidean length of a complex vector. In future marks of the Library, routines may be included in categories that are currently empty and further categories may be introduced.

3.1.2  BLAS names

Those routines which conform to the specifications of the BLAS may be called either by their NAG names or by their BLAS names.
In many implementations of the NAG Library, references to BLAS names may be linked to an efficient machine-specific implementation of the BLAS, usually provided by the vendor of the machine. Such implementations are stringently tested before being used with the NAG Library, to ensure that they correctly meet the specifications of the BLAS, and that they return the desired accuracy. Use of BLAS names is recommended for efficiency.
References to NAG routine names (beginning F06-) are always linked to the code provided in the NAG Library and may be significantly slower than the equivalent BLAS routine.
The names of the Level-2 and Level-3 BLAS follow a simple scheme (which is similar to that used for LAPACK routines in Chapters F07 and F08). Each name has the structure XYYZZZ, where the components have the following meanings:
the initial letter X indicates the data type (real or complex) and precision:
S real, single precision (in Fortran, REAL)
D real, double precision (in Fortran, DOUBLE PRECISION)
C complex, single precision (in Fortran, COMPLEX)
Z complex, double precision (in Fortran, COMPLEX*16 or DOUBLE COMPLEX)
the second and third letters YY indicate the type of the matrix A (and in some cases its storage scheme):
GE general
GB general band
SY symmetric
SP symmetric (packed storage)
SB symmetric band
HE (complex) Hermitian
HP (complex) Hermitian (packed storage)
HB (complex) Hermitian band
TR triangular
TP triangular (packed storage)
TB triangular band
the remaining 1, 2 or 3 letters ZZZ indicate the computation performed:
MV matrix-vector product
MM matrix-matrix product
R rank-1 update
R2 rank-2 update
RK rank-k update
R2K rank-2k update
SV solve a system of linear equations
SM solve a system of linear equations with a matrix of right-hand sides
Thus the routine DGEMV performs a matrix-vector product involving a real general matrix in double precision; the corresponding routine for a complex general matrix is ZGEMV.
The names of the Level-1 BLAS mostly follow the same convention for the initial letter (S-, C-, D- or Z-), except for a few involving data of mixed type, where the first two characters are precision-dependent.

3.1.3  LAPACK names

There are some LAPACK routines in this chapter that have BLAS-like functionalty. Four are equivalent to BLAS routines but for matrices stored in Rectangular Full Packed (RFP) format. The naming convention for these is as above with the addition of the matrix types:
HF (complex) Hermitian (RFP storage)
TF triangular (RFP storage)
SF symmetric (RFP storage)
There are an additonal two that compute norms of RFP matrices. These have second and third letters LA (signifying LAPACK), fourth letter N (signifying norm), and fifth and sixth letter signifying matrix type as above. For example ZLANHF computes the norm of a Hermitian matrix in RFP format.

3.2  The Level-0 Scalar Routines

The Level-0 routines perform operations on scalars or on vectors or matrices of order 2.

3.3  The Level-1 Vector Routines

The Level-1 routines perform operations either on a single vector or on a pair of vectors.

3.4  The Level-2 Matrix-vector and Matrix Routines

The Level-2 routines perform operations involving either a matrix on its own, or a matrix and one or more vectors.

3.5  The Level-3 Matrix-matrix Routines

The Level-3 routines perform operations involving matrix-matrix products.

3.6  Vector Arguments

Vector arguments (except in the Level-1 Sparse BLAS) are represented by a one-dimensional array, immediately followed by an increment argument whose name consists of the three characters INC followed by the name of the array. For example, a vector x is represented by the two arguments X and INCX. The length of the vector, n say, is passed as a separate argument, N.
The increment argument is the spacing (stride) in the array between the elements of the vector. For instance, if INCX=2, then the elements of x are in locations x1,x3,,x2n-1 of the array X and the intermediate locations x2,x4,,x2n-2 are not referenced.
When INCX>0, the vector element xi is in the array element X1+i-1×INCX. When INCX0, the elements are stored in the reverse order so that the vector element xi is in the array element X1-n-i×INCX and hence, in particular, the element xn is in X1. The declared length of the array X in the calling subroutine must be at least 1+N-1×INCX.
Negative increments are permitted only for:
Zero increments are formally permitted for Level-1 routines with more than one argument (in which case the element X1 is accessed repeatedly), but their use is strongly discouraged since the effect may be implementation-dependent. There is usually an alternative routine in this chapter, with a simplified parameter list, to achieve the required purpose. Zero increments are not permitted in the Level-2 BLAS.
In the Level-1 Sparse BLAS, each routine operates on two vectors x and y. The vector x is stored as a compressed sparse vector, and is represented by the three arguments NZ, X and INDX; NZ is the number of ‘interesting’ (usually nonzero) elements of x, and INDX is a one-dimensional index array such that
xINDXk=Xk,  k=1,2,,NZ.
The (mathematical) length of the vector, n say, does not need to be supplied; it is assumed that 1INDXkn. For example, the vector
x=0,4,0,0,1,0,0,0,6,0
could be represented with NZ=3, X=4,1,6, INDX=2,5,9. The second vector y is stored conventionally, and is represented simply by the one-dimensional array Y, with yi in Yi; the increment is assumed to be 1. Only the elements YINDXk are referenced.
Non-positive values of NZ are permitted, in which case the routines return immediately — except that functions set their value to zero before returning. For those routines where Y is an output argument the values in the array INDX must be distinct; violating this condition may yield incorrect results.

3.7  Matrix Arguments and Storage Schemes

In this chapter the following different storage schemes are used for matrices:
These storage schemes are compatible with those used in Chapters F07 and F08. (Different schemes for packed or band storage are used in a few older routines in Chapters F01, F02, F03 and F04.)
Chapter F01 provides some utility routines for conversion between storage schemes.
In the examples, * indicates an array element which need not be set and is not referenced by the routines. The examples illustrate only the relevant leading rows and columns of the arrays; array arguments may of course have additional rows or columns, according to the usual rules for passing array arguments in Fortran.

3.7.1  Conventional storage

Please see Section 3.3.1 in the F07 Chapter Introduction for full details.

3.7.2  Packed storage

Please see Section 3.3.2 in the F07 Chapter Introduction for full details.

3.7.3  Rectangular Full Packed (RFP) storage

Please see Section 3.3.3 in the F07 Chapter Introduction for full details.

3.7.4  Band storage

Please see Section 3.3.4 in the F07 Chapter Introduction for full details.

3.7.5  Unit triangular matrices

Please see Section 3.3.5 in the F07 Chapter Introduction for full details.

3.7.6  Real diagonal elements of complex Hermitian matrices

Please see Section 3.3.6 in the F07 Chapter Introduction for full details.

3.7.7  Spiked matrices

A few routines in this chapter (F06QSF, F06QWF, F06TSF and F06TWF) deal with upper spiked matrices. These are upper triangular matrices with an additional nonzero row or column below the diagonal.
The position of the spike is defined by indices k1 and k2; it is assumed that k1<k2. A row spike has nonzero elements in the k2th row, ak2,k for k=k1,k1+1,,k2-1; a column spike has nonzero elements in the k1th column, ak+1,k1 for k=k1,k1+1,,k2-1. For example, when n=6, k1=2 and k2=5:
Row spike Column spike
a11 a12 a13 a14 a15 a16 a22 a23 a24 a25 a26 a33 a34 a35 a36 a44 a45 a46 a52 a53 a54 a55 a56 a66 a11 a12 a13 a14 a15 a16 a22 a23 a24 a25 a26 a32 a33 a34 a35 a36 a42 a44 a45 a46 a52 a55 a56 a66
The storage scheme adopted by the routines in this chapter is for the upper triangular part of the spiked matrix to be stored conventionally in a two-dimensional array A, with the subdiagonal elements of the spike stored in a separate vector.

3.8  Option Parameters

Many of the routines in this chapter have one or more option parameters, of type CHARACTER. The descriptions in the routine documents refer only to upper-case values (for example UPLO='U' or UPLO='L'); however, in every case, the corresponding lower-case characters may be supplied (with the same meaning). Any other value is illegal.
A longer character string can be passed as the actual parameter, making the calling program more readable, but only the first character is significant. (This is a feature of Fortran.) For example:
 CALL DTRSV('Upper','Transpose','Non-unit',...)
The following option arguments are used in this chapter:

3.8.1  Matrix norms

The option argument NORM specifies different matrix norms whose definitions are given here for reference (for a general m by n matrix A):
If A is symmetric or Hermitian, A1=A.
The argument NORM can also be used to specify the maximum absolute value maxi,jaij (if NORM='M'), but this is not a norm in the strict mathematical sense.

3.9  Error Handling

Routines in this chapter do not use the usual NAG Library error-handling mechanism, involving the parameter IFAIL.
If one of the Level-2 or Level-3 BLAS routines is called with an invalid value of one of its arguments, then an error message is output on the error message unit (see X04AAF), giving the name of the routine and the number of the first invalid argument, and execution of the program is terminated. The following values of arguments are invalid:
Zero values for the matrix dimensions M, N or K are considered valid.
The other routines in this chapter do not report any errors in their arguments. Normally, if called, for example, with an unspecified value for one of the option arguments, or with a negative value of one of the problem dimensions M or N, they simply do nothing and return immediately.

4  Functionality Index

Level 0 (Scalar) operations, 
    complex numbers, 
        apply similarity rotation to 2 by 2 Hermitian matrix F06CHF
        generate a plane rotation, storing the tangent, real cosine F06CAF
        generate a plane rotation, storing the tangent, real sine F06CBF
        quotient of two numbers, with overflow flag F06CLF
        recover cosine and sine from given tangent, real cosine F06CCF
        recover cosine and sine from given tangent, real sine F06CDF
    real numbers, 
        apply similarity rotation to 2 by 2 symmetric matrix F06BHF
        compute (a2 + b2)1 / 2 F06BNF
        compute Euclidean norm from scaled form F06BMF
        eigenvalue of 2 by 2 symmetric matrix F06BPF
        generate a Jacobi plane rotation F06BEF
        generate a plane rotation F06AAF (DROTG)
        generate a plane rotation storing the tangent F06BAF
        quotient of two numbers, with overflow flag F06BLF
        recover cosine and sine from given tangent F06BCF
Level 1 (Vector) operations, 
    complex vector(s), 
        add scalar times a vector to another vector F06GCF (ZAXPY)
        apply a complex plane rotation F06HPF
        apply an elementary reflection to a vector F06HTF
        apply a real plane rotation F06KPF
        broadcast a scalar into a vector F06HBF
        copy a real vector to a complex vector F06KFF
        copy a vector F06GFF (ZCOPY)
        dot product of two vectors, conjugated F06GBF (ZDOTC)
        dot product of two vectors, unconjugated F06GAF (ZDOTU)
        Euclidean norm of a vector F06JJF (DZNRM2)
        generate an elementary reflection F06HRF
        generate a sequence of plane rotations F06HQF
        index of element of largest absolute value F06JMF (IZAMAX)
        multiply vector by a complex scalar F06GDF (ZSCAL)
        multiply vector by a complex scalar, preserving input vector F06HDF
        multiply vector by a real scalar F06JDF (ZDSCAL)
        multiply vector by a real scalar, preserving input vector F06KDF
        multiply vector by complex diagonal matrix F06HCF
        multiply vector by real diagonal matrix F06KCF
        multiply vector by reciprocal of a real scalar F06KEF
        negate a vector F06HGF
        sum of absolute values of vector-elements F06JKF (DZASUM)
        swap two vectors F06GGF (ZSWAP)
        update Euclidean norm in scaled form F06KJF
    Complex vector(s), 
        apply plane rotation, 
            real cosine, complex sine F06HMF (ZROT)
    integer vector(s), 
        broadcast a scalar into a vector F06DBF
        copy a vector F06DFF
    real vector(s), 
        add scalar times a vector to another vector F06ECF (DAXPY)
        apply an elementary reflection to a vector (Linpack style) F06FUF
        apply an elementary reflection to a vector (NAG style) F06FTF
        apply a symmetric plane rotation to two vectors F06FPF
        apply plane rotation F06EPF (DROT)
        broadcast a scalar into a vector F06FBF
        copy a vector F06EFF (DCOPY)
        cosine of angle between two vectors F06FAF
        dot product of two vectors F06EAF (DDOT)
        elements of largest and smallest absolute value F06FLF
        Euclidean norm of a vector F06EJF (DNRM2)
        generate an elementary reflection (Linpack style) F06FSF
        generate an elementary reflection (NAG style) F06FRF
        generate a sequence of plane rotations F06FQF
        index of element of largest absolute value F06JLF (IDAMAX)
        index of last non-negligible element F06KLF
        multiply vector by a scalar F06EDF (DSCAL)
        multiply vector by a scalar, preserving input vector F06FDF
        multiply vector by diagonal matrix F06FCF
        multiply vector by reciprocal of a scalar F06FEF
        negate a vector F06FGF
        sum of absolute values of vector-elements F06EKF (DASUM)
        swap two vectors F06EGF (DSWAP)
        update Euclidean norm in scaled form F06FJF
        weighted Euclidean norm of a vector F06FKF
Level 2 (Matrix-vector and matrix) operations, 
    complex matrix and vector(s), 
        apply sequence of plane rotations to a rectangular matrix, 
            complex cosine, real sine F06TYF
            real cosine, complex sine F06TXF
            real cosine and sine F06VXF
        compute a norm or the element of largest absolute value, 
            band matrix F06UBF
            general matrix F06UAF
            Hermitian band matrix F06UEF
            Hermitian matrix F06UCF
            Hermitian matrix, packed form F06UDF
            Hermitian matrix, RFP format F06WNF (ZLANHF)
            Hermitian tridiagonal matrix F06UPF
            Hessenberg matrix F06UMF
            symmetric band matrix F06UHF
            symmetric matrix F06UFF
            symmetric matrix, packed form F06UGF
            trapezoidal matrix F06UJF
            triangular band matrix F06ULF
            triangular matrix, packed form F06UKF
            tridiagonal matrix F06UNF
        compute upper Hessenberg matrix by applying sequence of plane rotations to an upper triangular matrix F06TVF
        compute upper spiked matrix by applying sequence of plane rotations to an upper triangular matrix F06TWF
        matrix initialization F06THF
        matrix-vector product, 
            Hermitian band matrix F06SDF (ZHBMV)
            Hermitian matrix F06SCF (ZHEMV)
            Hermitian packed matrix F06SEF (ZHPMV)
            rectangular band matrix F06SBF (ZGBMV)
            rectangular matrix F06SAF (ZGEMV)
            symmetric matrix F06TAF
            symmetric packed matrix F06TCF
            triangular band matrix F06SGF (ZTBMV)
            triangular matrix F06SFF (ZTRMV)
            triangular packed matrix F06SHF (ZTPMV)
        permute rows or columns of a matrix, 
            permutations represented by an integer array F06VJF
            permutations represented by a real array F06VKF
        QR factorization by sequence of plane rotations, 
            of rank-1 update of upper triangular matrix F06TPF
            of upper triangular matrix augmented by a full row F06TQF
        QR factorization of UZ or RQ factorization of ZU, where U is upper triangular and Z is a sequence of plane rotations F06TTF
        QR or RQ factorization by sequence of plane rotations, 
            of upper Hessenberg matrix F06TRF
            of upper spiked matrix F06TSF
        rank-1 update, 
            Hermitian matrix F06SPF (ZHER)
            Hermitian packed matrix F06SQF (ZHPR)
            rectangular matrix, conjugated vector F06SNF (ZGERC)
            rectangular matrix, unconjugated vector F06SMF (ZGERU)
            symmetric matrix F06TBF
            symmetric packed matrix F06TDF
        rank-2 update, 
            Hermitian matrix F06SRF (ZHER2)
            Hermitian packed matrix F06SSF (ZHPR2)
            matrix copy, rectangular or trapezoidal F06TFF
        solution of a system of equations, 
            triangular band matrix F06SKF (ZTBSV)
            triangular matrix F06SJF (ZTRSV)
            triangular packed matrix F06SLF (ZTPSV)
        unitary similarity transformation of a Hermitian matrix, 
            as sequence of plane rotations F06TMF
    real matrix and vector(s), 
        apply sequence of plane rotations to a rectangular matrix F06QXF
        compute a norm or the element of largest absolute value, 
            band matrix F06RBF
            general matrix F06RAF
            Hessenberg matrix F06RMF
            matrix initialization F06QHF
            symmetric band matrix F06REF
            symmetric matrix F06RCF
            symmetric matrix, packed form F06RDF
            symmetric matrix, RFP format F06WAF (DLANSF)
            symmetric tridiagonal matrix F06RPF
            trapezoidal matrix F06RJF
            triangular band matrix F06RLF
            triangular matrix, packed form F06RKF
            tridiagonal matrix F06RNF
        compute upper Hessenberg matrix by applying sequence of plane rotations to an upper triangular matrix F06QVF
        compute upper spiked matrix by applying sequence of plane rotations to an upper triangular matrix F06QWF
        matrix-vector product, 
            rectangular band matrix F06PBF (DGBMV)
            rectangular matrix F06PAF (DGEMV)
            symmetric band matrix F06PDF (DSBMV)
            symmetric matrix F06PCF (DSYMV)
            symmetric packed matrix F06PEF (DSPMV)
            triangular band matrix F06PGF (DTBMV)
            triangular matrix F06PFF (DTRMV)
            triangular packed matrix F06PHF (DTPMV)
        orthogonal similarity transformation of a symmetric matrix, 
            as sequence of plane rotations F06QMF
        permute rows or columns of a matrix, 
            permutations represented by an integer array F06QJF
            permutations represented by a real array F06QKF
        QR factorization by sequence of plane rotations, 
            of rank-1 update of upper triangular matrix F06QPF
            of upper triangular matrix augmented by a full row F06QQF
        QR factorization of UZ or RQ factorization of ZU, where U is upper triangular and Z is a sequence of plane rotations F06QTF
        QR or RQ factorization by sequence of plane rotations, 
            of upper Hessenberg matrix F06QRF
            of upper spiked matrix F06QSF
        rank-1 update, 
            rectangular matrix F06PMF (DGER)
            symmetric matrix F06PPF (DSYR)
            symmetric packed matrix F06PQF (DSPR)
        rank-2 update, 
            matrix copy, rectangular or trapezoidal F06QFF
            symmetric matrix F06PRF (DSYR2)
            symmetric packed matrix F06PSF (DSPR2)
        solution of a system of equations, 
            triangular band matrix F06PKF (DTBSV)
            triangular matrix F06PJF (DTRSV)
            triangular packed matrix F06PLF (DTPSV)
Level 3 (Matrix-matrix) operations, 
    complex matrices, 
        matrix-matrix product, 
            one matrix Hermitian F06ZCF (ZHEMM)
            one matrix symmetric F06ZTF (ZSYMM)
            one matrix triangular F06ZFF (ZTRMM)
            two rectangular matrices F06ZAF (ZGEMM)
        rank-2k update, 
            of a Hermitian matrix F06ZRF (ZHER2K)
            of a symmetric matrix F06ZWF (ZSYR2K)
        rank-k update, 
            of a Hermitian matrix F06ZPF (ZHERK)
            of a Hermitian matrix, RFP format F06WQF (ZHFRK)
            of a symmetric matrix F06ZUF (ZSYRK)
        solution of triangular systems of equations F06ZJF (ZTRSM)
        solution of triangular systems of equations, RFP format F06WPF (ZTFSM)
    real matrices, 
        matrix-matrix product, 
            one matrix symmetric F06YCF (DSYMM)
            one matrix triangular F06YFF (DTRMM)
            rectangular matrices F06YAF (DGEMM)
        rank-2k update of a symmetric matrix F06YRF (DSYR2K)
        rank-k update, 
            of a symmetric matrix F06YPF (DSYRK)
            of a symmetric matrix, RFP format F06WCF (DSFRK)
        solution of triangular systems of equations F06YJF (DTRSM)
        solution of triangular systems of equations, RFP format F06WBF (DTFSM)
Sparse level 1 (vector) operations, 
    complex vector(s), 
        add scalar times sparse vector to a full vector F06GTF (ZAXPYI)
        dot product of a sparse and a full vector (conjugated) F06GSF (ZDOTCI)
        dot product of a sparse and a full vector (unconjugated) F06GRF (ZDOTUI)
        gather and set to zero a sparse vector F06GVF (ZGTHRZ)
        gather sparse vector F06GUF (ZGTHR)
        scatter sparse vector F06GWF (ZSCTR)
    real vector(s), 
        add scalar times sparse vector to a full vector F06ETF (DAXPYI)
        apply plane rotation to a sparse and a full vector F06EXF (DROTI)
        dot product of a sparse and a full vector F06ERF (DDOTI)
        gather and set to zero a sparse vector F06EVF (DGTHRZ)
        gather sparse vector F06EUF (DGTHR)
        scatter sparse vector F06EWF (DSCTR)
BLAS Routines

5  Auxiliary Routines Associated with Library Routine Parameters

None.

6  Routines Withdrawn or Scheduled for Withdrawal

None.

7  References

Dodson D S and Grimes R G (1982) Remark on Algorithm 539 ACM Trans. Math. Software 8 403–404
Dodson D S, Grimes R G and Lewis J G (1991) Sparse extensions to the Fortran basic linear algebra subprograms ACM Trans. Math. Software 17 253–263
Dongarra J J, Du Croz J J, Duff I S and Hammarling S (1990) A set of Level 3 basic linear algebra subprograms ACM Trans. Math. Software 16 1–28
Dongarra J J, Du Croz J J, Hammarling S and Hanson R J (1988) An extended set of FORTRAN basic linear algebra subprograms ACM Trans. Math. Software 14 1–32
Dongarra J J, Moler C B, Bunch J R and Stewart G W (1979) LINPACK Users' Guide SIAM, Philadelphia
Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore
Lawson C L, Hanson R J, Kincaid D R and Krogh F T (1979) Basic linear algebra supbrograms for Fortran usage ACM Trans. Math. Software 5 308–325

F06 Chapter Contents
F06 Chapter Introduction (PDF version)
NAG Library Manual

© The Numerical Algorithms Group Ltd, Oxford, UK. 2012