Skip to content

Files

This commit does not belong to any branch on this repository, and may belong to a fork outside of the repository.

Latest commit

138b53c · May 16, 2025

History

History
354 lines (236 loc) · 11.3 KB

File metadata and controls

354 lines (236 loc) · 11.3 KB

dgemv

Perform one of the matrix-vector operations y = α*A*x + β*y or y = α*A**T*x + β*y.

Usage

var dgemv = require( '@stdlib/blas/base/dgemv' );

dgemv( ord, trans, M, N, α, A, LDA, x, sx, β, y, sy )

Performs one of the matrix-vector operations y = α*A*x + β*y or y = α*A**T*x + β*y, where α and β are scalars, x and y are vectors, and A is an M by N matrix.

var Float64Array = require( '@stdlib/array/float64' );

var A = new Float64Array( [ 1.0, 2.0, 3.0, 4.0, 5.0, 6.0 ] );
var x = new Float64Array( [ 1.0, 1.0, 1.0 ] );
var y = new Float64Array( [ 1.0, 1.0 ] );

dgemv( 'row-major', 'no-transpose', 2, 3, 1.0, A, 3, x, 1, 1.0, y, 1 );
// y => <Float64Array>[ 7.0, 16.0 ]

The function has the following parameters:

  • ord: storage layout.
  • trans: specifies whether A should be transposed, conjugate-transposed, or not transposed.
  • M: number of rows in the matrix A.
  • N: number of columns in the matrix A.
  • α: scalar constant.
  • A: input matrix stored in linear memory as a Float64Array.
  • lda: stride of the first dimension of A (leading dimension of A).
  • x: input Float64Array.
  • sx: index increment for x.
  • β: scalar constant.
  • y: output Float64Array.
  • sy: index increment for y.

The stride parameters determine how operations are performed. For example, to iterate over every other element in x and y,

var Float64Array = require( '@stdlib/array/float64' );

var A = new Float64Array( [ 1.0, 2.0, 3.0, 4.0 ] );
var x = new Float64Array( [ 1.0, 0.0, 1.0, 0.0 ] );
var y = new Float64Array( [ 1.0, 0.0, 1.0, 0.0 ] );

dgemv( 'row-major', 'no-transpose', 2, 2, 1.0, A, 2, x, 2, 1.0, y, 2 );
// y => <Float64Array>[ 4.0, 0.0, 8.0, 0.0 ]

Note that indexing is relative to the first index. To introduce an offset, use typed array views.

var Float64Array = require( '@stdlib/array/float64' );

// Initial arrays...
var x0 = new Float64Array( [ 0.0, 1.0, 1.0 ] );
var y0 = new Float64Array( [ 0.0, 1.0, 1.0 ] );
var A = new Float64Array( [ 1.0, 2.0, 3.0, 4.0, 5.0, 6.0 ] );

// Create offset views...
var x1 = new Float64Array( x0.buffer, x0.BYTES_PER_ELEMENT*1 ); // start at 2nd element
var y1 = new Float64Array( y0.buffer, y0.BYTES_PER_ELEMENT*1 ); // start at 2nd element

dgemv( 'row-major', 'no-transpose', 2, 2, 1.0, A, 2, x1, -1, 1.0, y1, -1 );
// y0 => <Float64Array>[ 0.0, 8.0, 4.0 ]

dgemv.ndarray( trans, M, N, α, A, sa1, sa2, oa, x, sx, ox, β, y, sy, oy )

Performs one of the matrix-vector operations y = α*A*x + β*y or y = α*A**T*x + β*y, using alternative indexing semantics and where α and β are scalars, x and y are vectors, and A is an M by N matrix.

var Float64Array = require( '@stdlib/array/float64' );

var A = new Float64Array( [ 1.0, 2.0, 3.0, 4.0, 5.0, 6.0 ] );
var x = new Float64Array( [ 1.0, 1.0, 1.0 ] );
var y = new Float64Array( [ 1.0, 1.0 ] );

dgemv.ndarray( 'no-transpose', 2, 3, 1.0, A, 3, 1, 0, x, 1, 0, 1.0, y, 1, 0 );
// y => <Float64Array>[ 7.0, 16.0 ]

The function has the following additional parameters:

  • sa1: stride of the first dimension of A.
  • sa2: stride of the second dimension of A.
  • oa: starting index for A.
  • ox: starting index for x.
  • oy: starting index for y.

While typed array views mandate a view offset based on the underlying buffer, the offset parameters support indexing semantics based on starting indices. For example,

var Float64Array = require( '@stdlib/array/float64' );

var A = new Float64Array( [ 1.0, 2.0, 3.0, 4.0, 5.0, 6.0 ] );
var x = new Float64Array( [ 0.0, 1.0, 2.0, 3.0 ] );
var y = new Float64Array( [ 7.0, 8.0, 9.0, 10.0 ] );

dgemv.ndarray( 'no-transpose', 2, 3, 1.0, A, 3, 1, 0, x, 1, 1, 1.0, y, -2, 2 );
// y => <Float64Array>[ 39, 8, 23, 10 ]

Notes

  • dgemv() corresponds to the BLAS level 2 function dgemv.

Examples

var discreteUniform = require( '@stdlib/random/array/discrete-uniform' );
var dgemv = require( '@stdlib/blas/base/dgemv' );

var opts = {
    'dtype': 'float64'
};

var M = 3;
var N = 3;

var A = discreteUniform( M*N, 0, 255, opts );
var x = discreteUniform( N, 0, 255, opts );
var y = discreteUniform( M, 0, 255, opts );

dgemv( 'row-major', 'no-transpose', M, N, 1.0, A, N, x, -1, 1.0, y, -1 );
console.log( y );

C APIs

Usage

#include "stdlib/blas/base/dgemv.h"

c_dgemv( order, trans, M, N, alpha, *A, LDA, *X, strideX, beta, *Y, strideY )

Performs one of the matrix-vector operations y = α*A*x + β*y or y = α*A^T*x + β*y, where α and β are scalars, x and y are vectors, and A is an M by N matrix.

#include "stdlib/blas/base/shared.h"

double A[] = { 1.0, 0.0, 0.0, 2.0, 1.0, 0.0, 3.0, 2.0, 1.0 };
const double x[] = { 1.0, 2.0, 3.0 };
const double y[] = { 1.0, 2.0, 3.0 };

c_dgemv( CblasColMajor, CblasNoTrans, 3, 3, 1.0, A, 3, x, 1, 1.0, y, 1 );

The function accepts the following arguments:

  • order: [in] CBLAS_LAYOUT storage layout.
  • trans: [in] CBLAS_TRANSPOSE specifies whether A should be transposed, conjugate-transposed, or not transposed.
  • M: [in] CBLAS_INT number of rows in the matrix A.
  • N: [in] CBLAS_INT number of columns in the matrix A.
  • alpha: [in] double scalar.
  • A: [inout] double* input matrix.
  • LDA: [in] CBLAS_INT stride of the first dimension of A (a.k.a., leading dimension of the matrix A).
  • X: [in] double* first input vector.
  • strideX: [in] CBLAS_INT index increment for X.
  • beta: [in] double scalar.
  • Y: [in] double* second input vector.
  • strideY: [in] CBLAS_INT index increment for Y.
void c_dgemv( const CBLAS_LAYOUT order, const CBLAS_TRANSPOSE trans, const CBLAS_INT M, const CBLAS_INT N, const double alpha, const double *A, const CBLAS_INT LDA, const double *x, const CBLAS_INT strideX, const double beta, double *y, const CBLAS_INT strideY )

c_dgemv_ndarray( trans, M, N, alpha, *A, strideA1, strideA2, offsetA, *X, strideX, offsetX, beta, *Y, strideY, offsetY )

Performs one of the matrix-vector operations y = α*A*x + β*y or y = α*A^T*x + β*y, where α and β are scalars, x and y are vectors, and A is an M by N matrix using indexing alternative semantics.

#include "stdlib/blas/base/shared.h"

double A[] = { 1.0, 0.0, 0.0, 2.0, 1.0, 0.0, 3.0, 2.0, 1.0 };
const double x[] = { 1.0, 2.0, 3.0 };
const double y[] = { 1.0, 2.0, 3.0 };

c_dgemv_ndarray( CblasNoTrans, 3, 3, 1.0, A, 1, 3, 0, x, 1, 0, 1.0, y, 1, 0 );

The function accepts the following arguments:

  • trans: [in] CBLAS_TRANSPOSE specifies whether A should be transposed, conjugate-transposed, or not transposed.
  • M: [in] CBLAS_INT number of rows in the matrix A.
  • N: [in] CBLAS_INT number of columns in the matrix A.
  • alpha: [in] double scalar.
  • A: [inout] double* input matrix.
  • strideA1: [in] CBLAS_INT stride of the first dimension of A.
  • strideA2: [in] CBLAS_INT stride of the second dimension of A.
  • offsetA: [in] CBLAS_INT starting index for A.
  • X: [in] double* first input vector.
  • strideX: [in] CBLAS_INT index increment for X.
  • offsetX: [in] CBLAS_INT starting index for X.
  • beta: [in] double scalar.
  • Y: [in] double* second input vector.
  • strideY: [in] CBLAS_INT index increment for Y.
  • offsetY: [in] CBLAS_INT starting index for Y.
void c_dgemv_ndarray( const CBLAS_TRANSPOSE trans, const CBLAS_INT M, const CBLAS_INT N, const double alpha, const double *A, const CBLAS_INT strideA1, const CBLAS_INT strideA2, const CBLAS_INT offsetA, const double *x, const CBLAS_INT strideX, const CBLAS_INT offsetX, const double beta, double *y, const CBLAS_INT strideY, const CBLAS_INT offsetY )

Examples

#include "stdlib/blas/base/dgemv.h"
#include "stdlib/blas/base/shared.h"
#include <stdio.h>

int main( void ) {
    // Create a strided array:
    const double A[] = { 1.0, 0.0, 0.0, 2.0, 1.0, 0.0, 3.0, 2.0, 1.0 };
    const double x[] = { 1.0, 2.0, 3.0 };
    double y[] = { 1.0, 2.0, 3.0 };

    // Specify the number of elements along each dimension of `A`:
    const int M = 3;
    const int N = 3;

    // Perform the matrix-vector operations `y = α*A*x + β*y`:
    c_dgemv( CblasRowMajor, CblasNoTrans, M, N, 1.0, A, M, x, 1, 1.0, y, 1 );

    // Print the result:
    for ( int i = 0; i < N; i++ ) {
        printf( "y[ %i ] = %lf\n", i, y[ i ] );
    }

    // Perform the symmetric rank 2 operation `A = α*x*y^T + α*y*x^T + A`:
    c_dgemv_ndarray( CblasNoTrans, 3, 3, 1.0, A, 3, 1, 0, x, 1, 0, 1.0, y, 1, 0 );

    // Print the result:
    for ( int i = 0; i < N; i++ ) {
        printf( "y[ %i ] = %lf\n", i, y[ i ] );
    }
}