@stdlib/lapack-base-dlacpy
v0.2.1
Published
Copy all or part of a matrix A to another matrix B.
Readme
dlacpy
Copy all or part of a matrix
Ato another matrixB.
Installation
npm install @stdlib/lapack-base-dlacpyUsage
var dlacpy = require( '@stdlib/lapack-base-dlacpy' );dlacpy( order, uplo, M, N, A, LDA, B, LDB )
Copies all or part of a matrix A to another matrix B.
var Float64Array = require( '@stdlib/array-float64' );
var A = new Float64Array( [ 1.0, 2.0, 3.0, 4.0 ] );
var B = new Float64Array( 4 );
dlacpy( 'row-major', 'all', 2, 2, A, 2, B, 2 );
// B => <Float64Array>[ 1.0, 2.0, 3.0, 4.0 ]The function has the following parameters:
- order: storage layout.
- uplo: specifies whether to copy the upper or lower triangular/trapezoidal part of a matrix
A. - M: number of rows in
A. - N: number of columns in
A. - A: input
Float64Array. - LDA: stride of the first dimension of
A(a.k.a., leading dimension of the matrixA). - B: output
Float64Array. - LDB: stride of the first dimension of
B(a.k.a., leading dimension of the matrixB).
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 A0 = new Float64Array( [ 1.0, 2.0, 3.0, 4.0, 5.0 ] );
var B0 = new Float64Array( 5 );
// Create offset views...
var A1 = new Float64Array( A0.buffer, A0.BYTES_PER_ELEMENT*1 ); // start at 2nd element
var B1 = new Float64Array( B0.buffer, B0.BYTES_PER_ELEMENT*1 ); // start at 2nd element
dlacpy( 'row-major', 'all', 2, 2, A1, 2, B1, 2 );
// B0 => <Float64Array>[ 0.0, 2.0, 3.0, 4.0, 5.0 ]dlacpy.ndarray( uplo, M, N, A, sa1, sa2, oa, B, sb1, sb2, ob )
Copies all or part of a matrix A to another matrix B using alternative indexing semantics.
var Float64Array = require( '@stdlib/array-float64' );
var A = new Float64Array( [ 1.0, 2.0, 3.0, 4.0 ] );
var B = new Float64Array( [ 0.0, 0.0, 0.0, 0.0 ] );
dlacpy.ndarray( 'all', 2, 2, A, 2, 1, 0, B, 2, 1, 0 );
// B => <Float64Array>[ 1.0, 2.0, 3.0, 4.0 ]The function has the following parameters:
- uplo: specifies whether to copy the upper or lower triangular/trapezoidal part of a matrix
A. - M: number of rows in
A. - N: number of columns in
A. - A: input
Float64Array. - sa1: stride of the first dimension of
A. - sa2: stride of the second dimension of
A. - oa: starting index for
A. - B: output
Float64Array. - sb1: stride of the first dimension of
B. - sb2: stride of the second dimension of
B. - ob: starting index for
B.
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 ] );
var B = new Float64Array( [ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 ] );
dlacpy.ndarray( 'all', 2, 2, A, 2, 1, 1, B, 2, 1, 2 );
// B => <Float64Array>[ 0.0, 0.0, 2.0, 3.0, 4.0, 5.0 ]Notes
Examples
var ndarray2array = require( '@stdlib/ndarray-base-to-array' );
var uniform = require( '@stdlib/random-array-discrete-uniform' );
var numel = require( '@stdlib/ndarray-base-numel' );
var shape2strides = require( '@stdlib/ndarray-base-shape2strides' );
var dlacpy = require( '@stdlib/lapack-base-dlacpy' );
var shape = [ 5, 8 ];
var order = 'row-major';
var strides = shape2strides( shape, order );
var N = numel( shape );
var A = uniform( N, -10, 10, {
'dtype': 'float64'
});
console.log( ndarray2array( A, shape, strides, 0, order ) );
var B = uniform( N, -10, 10, {
'dtype': 'float64'
});
console.log( ndarray2array( B, shape, strides, 0, order ) );
dlacpy( order, 'all', shape[ 0 ], shape[ 1 ], A, strides[ 0 ], B, strides[ 0 ] );
console.log( ndarray2array( B, shape, strides, 0, order ) );C APIs
Usage
#include "stdlib/lapack/base/dlacpy.h"c_dlacpy( layout, uplo, M, N, *A, LDA, *B, LDB )
Copies all or part of a matrix A to another matrix B.
#include "stdlib/lapack/base/shared.h"
const double A[] = { 1.0, 2.0, 3.0, 4.0 };
double B[] = { 0.0, 0.0, 0.0, 0.0 };
c_dlacpy( LAPACK_ROW_MAJOR, LAPACK_UPPER_TRIANGLE, 2, 2, A, 2, B, 2 );The function accepts the following arguments:
- order:
[in] LAPACK_LAYOUTstorage layout. - uplo:
[in] intspecifies whether to copy the upper or lower triangular/trapezoidal part of a matrixA. - M:
[in] LAPACK_INTnumber of rows inA. - N:
[in] LAPACK_INTnumber of columns inA. - A:
[in] double*input matrix. - LDA:
[in] LAPACK_INTstride of the first dimension ofA(a.k.a., leading dimension of the matrixA). - B:
[out] double*output matrix. - LDB:
[in] LAPACK_INTstride of the first dimension ofB(a.k.a., leading dimension of the matrixB).
LAPACK_INT c_dlacpy( const LAPACK_LAYOUT layout, const int uplo, const LAPACK_INT M, const LAPACK_INT N, const double *A, const LAPACK_INT LDA, double *B, const LAPACK_INT LDB );c_dlacpy_ndarray( uplo, M, N, *A, sa1, sa2, oa, *B, sb1, sb2, ob )
Copies all or part of a matrix A to another matrix B using alternative indexing semantics.
#include "stdlib/lapack/base/shared.h"
const double A[] = { 1.0, 2.0, 3.0, 4.0 };
double B[] = { 0.0, 0.0, 0.0, 0.0 };
c_dlacpy_ndarray( LAPACK_UPPER_TRIANGLE, 2, 2, A, 2, 1, 0, B, 2, 1, 0 );The function accepts the following arguments:
- uplo:
[in] intspecifies whether to copy the upper or lower triangular/trapezoidal part of a matrixA. - M:
[in] LAPACK_INTnumber of rows inA. - N:
[in] LAPACK_INTnumber of columns inA. - A:
[in] double*input matrix. - sa1:
[in] LAPACK_INTstride of the first dimension ofA. - sa2:
[in] LAPACK_INTstride of the second dimension ofA. - oa:
[in] LAPACK_INTstarting index forA. - B:
[out] double*output matrix. - sb1:
[in] LAPACK_INTstride of the first dimension ofB. - sb2:
[in] LAPACK_INTstride of the second dimension ofB. - ob:
[in] LAPACK_INTstarting index forB.
LAPACK_INT c_dlacpy_ndarray( const int uplo, const LAPACK_INT M, const LAPACK_INT N, const double *A, const LAPACK_INT strideA1, const LAPACK_INT strideA2, const LAPACK_INT offsetA, double *B, const LAPACK_INT strideB1, const LAPACK_INT strideB2, const LAPACK_INT offsetB );Examples
#include "stdlib/lapack/base/dlacpy.h"
#include "stdlib/lapack/base/shared.h"
#include <stdio.h>
int main( void ) {
// Define a 3x3 input matrix stored in row-major order:
const double A[ 3*3 ] = {
1.0, 2.0, 3.0,
4.0, 5.0, 6.0,
7.0, 8.0, 9.0
};
// Define a 3x3 output matrix:
double B[ 3*3 ] = {
0.0, 0.0, 0.0,
0.0, 0.0, 0.0,
0.0, 0.0, 0.0
};
// Specify the number of elements along each dimension of `A`:
const int M = 3;
const int N = 3;
// Copy elements from the upper triangle of `A` to `B`:
c_dlacpy( LAPACK_ROW_MAJOR, LAPACK_UPPER_TRIANGLE, M, N, A, N, B, N );
// Print the result:
for ( int i = 0; i < M; i++ ) {
for ( int j = 0; j < N; j++ ) {
printf( "B[ %i, %i ] = %lf\n", i, j, B[ (i*N)+j ] );
}
}
// Copy elements from the lower triangle of `A` to `B` using alternative indexing semantics:
c_dlacpy_ndarray( LAPACK_LOWER_TRIANGLE, M, N, A, N, 1, 0, B, N, 1, 0 );
// Print the result:
for ( int i = 0; i < M; i++ ) {
for ( int j = 0; j < N; j++ ) {
printf( "B[ %i, %i ] = %lf\n", i, j, B[ (i*N)+j ] );
}
}
}Notice
This package is part of stdlib, a standard library for JavaScript and Node.js, with an emphasis on numerical and scientific computing. The library provides a collection of robust, high performance libraries for mathematics, statistics, streams, utilities, and more.
For more information on the project, filing bug reports and feature requests, and guidance on how to develop stdlib, see the main project repository.
Community
License
See LICENSE.
Copyright
Copyright © 2016-2026. The Stdlib Authors.
