Visual Servoing Platform version 3.7.0
Loading...
Searching...
No Matches
vpMatrix_mul.cpp
1/*
2 * ViSP, open source Visual Servoing Platform software.
3 * Copyright (C) 2005 - 2024 by Inria. All rights reserved.
4 *
5 * This software is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or
8 * (at your option) any later version.
9 * See the file LICENSE.txt at the root directory of this source
10 * distribution for additional information about the GNU GPL.
11 *
12 * For using ViSP with software that can not be combined with the GNU
13 * GPL, please contact Inria about acquiring a ViSP Professional
14 * Edition License.
15 *
16 * See https://visp.inria.fr for more information.
17 *
18 * This software was developed at:
19 * Inria Rennes - Bretagne Atlantique
20 * Campus Universitaire de Beaulieu
21 * 35042 Rennes Cedex
22 * France
23 *
24 * If you have questions regarding the use of this file, please contact
25 * Inria at visp@inria.fr
26 *
27 * This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
28 * WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
29 *
30 * Description:
31 * BLAS subroutines.
32 */
33
34#include <visp3/core/vpConfig.h>
35#include <visp3/core/vpMatrix.h>
36
37#ifndef DOXYGEN_SHOULD_SKIP_THIS
38
39#if defined(VISP_HAVE_LAPACK)
41#ifdef VISP_HAVE_MKL
42#include <mkl.h>
43typedef MKL_INT integer;
44
45void vpMatrix::blas_dgemm(char trans_a, char trans_b, unsigned int M_, unsigned int N_, unsigned int K_, double alpha,
46 double *a_data, unsigned int lda_, double *b_data, unsigned int ldb_, double beta,
47 double *c_data, unsigned int ldc_)
48{
49 MKL_INT M = static_cast<MKL_INT>(M_);
50 MKL_INT N = static_cast<MKL_INT>(N_);
51 MKL_INT K = static_cast<MKL_INT>(K_);
52 MKL_INT lda = static_cast<MKL_INT>(lda_);
53 MKL_INT ldb = static_cast<MKL_INT>(ldb_);
54 MKL_INT ldc = static_cast<MKL_INT>(ldc_);
55
56 dgemm(&trans_a, &trans_b, &M, &N, &K, &alpha, a_data, &lda, b_data, &ldb, &beta, c_data, &ldc);
57}
58
59void vpMatrix::blas_dgemv(char trans, unsigned int M_, unsigned int N_, double alpha, double *a_data, unsigned int lda_,
60 double *x_data, int incx_, double beta, double *y_data, int incy_)
61{
62 MKL_INT M = static_cast<MKL_INT>(M_);
63 MKL_INT N = static_cast<MKL_INT>(N_);
64 MKL_INT lda = static_cast<MKL_INT>(lda_);
65 MKL_INT incx = static_cast<MKL_INT>(incx_);
66 MKL_INT incy = static_cast<MKL_INT>(incy_);
67
68 dgemv(&trans, &M, &N, &alpha, a_data, &lda, x_data, &incx, &beta, y_data, &incy);
69}
70
71#elif defined(VISP_HAVE_GSL)
72
73#include <gsl/gsl_blas.h>
74#include <gsl/gsl_matrix.h>
75#include <gsl/gsl_vector.h>
76
88 void vpMatrix::blas_dgemm(char trans_a, char trans_b, unsigned int M, unsigned int N, unsigned int K, double alpha,
89 double *A_data, unsigned int lda, double *B_data, unsigned int ldb, double beta,
90 double *C_data, unsigned int ldc)
91{
92 CBLAS_TRANSPOSE_t TransA = (trans_a == 'n' || trans_a == 'N') ? CblasNoTrans : CblasTrans;
93 CBLAS_TRANSPOSE_t TransB = (trans_b == 'n' || trans_b == 'N') ? CblasNoTrans : CblasTrans;
94
95 unsigned int A_rows = (TransA == CblasNoTrans) ? M : K;
96 unsigned int A_cols = (TransA == CblasNoTrans) ? K : M;
97 unsigned int B_rows = (TransB == CblasNoTrans) ? K : N;
98 unsigned int B_cols = (TransB == CblasNoTrans) ? N : K;
99
100 gsl_matrix_view A = gsl_matrix_view_array_with_tda(A_data, A_rows, A_cols, lda);
101 gsl_matrix_view B = gsl_matrix_view_array_with_tda(B_data, B_rows, B_cols, ldb);
102 gsl_matrix_view C = gsl_matrix_view_array_with_tda(C_data, M, N, ldc);
103
104 gsl_blas_dgemm(TransA, TransB, alpha, &A.matrix, &B.matrix, beta, &C.matrix);
105}
106
112void vpMatrix::blas_dgemv(char trans, unsigned int M, unsigned int N, double alpha, double *A_data, unsigned int lda,
113 double *x_data, int incx, double beta, double *y_data, int incy)
114{
115 CBLAS_TRANSPOSE_t Trans = (trans == 'n' || trans == 'N') ? CblasNoTrans : CblasTrans;
116
117 unsigned int A_rows = (Trans == CblasNoTrans) ? M : N;
118 unsigned int A_cols = (Trans == CblasNoTrans) ? N : M;
119
120 gsl_matrix_view A = gsl_matrix_view_array_with_tda(A_data, A_rows, A_cols, lda);
121 gsl_vector_view x = gsl_vector_view_array_with_stride(x_data, incx, N);
122 gsl_vector_view y = gsl_vector_view_array_with_stride(y_data, incy, M);
123
124 gsl_blas_dgemv(Trans, alpha, &A.matrix, &x.vector, beta, &y.vector);
125}
126
127#else
128#ifdef VISP_HAVE_LAPACK_BUILT_IN
129typedef long int integer;
130#else
131typedef int integer;
132#endif
133
134extern "C" void dgemm_(char *transa, char *transb, integer *M, integer *N, integer *K, double *alpha, double *a,
135 integer *lda, double *b, integer *ldb, double *beta, double *c, integer *ldc);
136
137extern "C" void dgemv_(char *trans, integer *M, integer *N, double *alpha, double *a, integer *lda, double *x,
138 integer *incx, double *beta, double *y, integer *incy);
139
140void vpMatrix::blas_dgemm(char trans_a, char trans_b, unsigned int M_, unsigned int N_, unsigned int K_, double alpha,
141 double *a_data, unsigned int lda_, double *b_data, unsigned int ldb_, double beta,
142 double *c_data, unsigned int ldc_)
143{
144 integer M = static_cast<integer>(M_);
145 integer K = static_cast<integer>(K_);
146 integer N = static_cast<integer>(N_);
147 integer lda = static_cast<integer>(lda_);
148 integer ldb = static_cast<integer>(ldb_);
149 integer ldc = static_cast<integer>(ldc_);
150
151 dgemm_(&trans_a, &trans_b, &M, &N, &K, &alpha, a_data, &lda, b_data, &ldb, &beta, c_data, &ldc);
152}
153
154void vpMatrix::blas_dgemv(char trans, unsigned int M_, unsigned int N_, double alpha, double *a_data, unsigned int lda_,
155 double *x_data, int incx_, double beta, double *y_data, int incy_)
156{
157 integer M = static_cast<integer>(M_);
158 integer N = static_cast<integer>(N_);
159 integer lda = static_cast<integer>(lda_);
160 integer incx = static_cast<integer>(incx_);
161 integer incy = static_cast<integer>(incy_);
162
163 dgemv_(&trans, &M, &N, &alpha, a_data, &lda, x_data, &incx, &beta, y_data, &incy);
164}
165#endif
166END_VISP_NAMESPACE
167#else
168// Work around to avoid warning LNK4221: This object file does not define any
169// previously undefined public symbols
170void dummy_vpMatrix_blas() { }
171#endif
172
173#endif // #ifndef DOXYGEN_SHOULD_SKIP_THIS