aboutsummaryrefslogtreecommitdiffhomepage
path: root/bench/btl/libs/C_BLAS/C_BLAS_interface.hh
blob: e637c0e5e58b39f7593aa89fb2a1bbe73bf35470 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
//=====================================================
// File   :  C_BLAS_interface.hh
// Author :  L. Plagne <laurent.plagne@edf.fr)>
// Copyright (C) EDF R&D,  lun sep 30 14:23:28 CEST 2002
//=====================================================
//
// This program is free software; you can redistribute it and/or
// modify it under the terms of the GNU General Public License
// as published by the Free Software Foundation; either version 2
// of the License, or (at your option) any later version.
//
// This program is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
// GNU General Public License for more details.
// You should have received a copy of the GNU General Public License
// along with this program; if not, write to the Free Software
// Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA  02111-1307, USA.
//
#ifndef C_BLAS_PRODUIT_MATRICE_VECTEUR_HH
#define C_BLAS_PRODUIT_MATRICE_VECTEUR_HH

#include "f77_interface.hh"

extern "C"
{
#include "cblas.h"
}

template<class real>
class C_BLAS_interface : public f77_interface_base<real>
{
public :

  typedef typename f77_interface_base<real>::gene_matrix gene_matrix;
  typedef typename f77_interface_base<real>::gene_vector gene_vector;

  static inline std::string name( void )
  {
    return "C_BLAS";
  }

  static  inline void matrix_vector_product(gene_matrix & A, gene_vector & B, gene_vector & X, int N)
  {
    cblas_dgemv(CblasColMajor,CblasNoTrans,N,N,1.0,A,N,B,1,0.0,X,1);
  }

  static  inline void atv_product(gene_matrix & A, gene_vector & B, gene_vector & X, int N)
  {
    cblas_dgemv(CblasColMajor,CblasTrans,N,N,1.0,A,N,B,1,0.0,X,1);
  }

  static  inline void matrix_matrix_product(gene_matrix & A, gene_matrix & B, gene_matrix & X, int N)
  {
    cblas_dgemm(CblasColMajor,CblasNoTrans,CblasNoTrans,N,N,N,1.0,A,N,B,N,0.0,X,N);
  }

  static  inline void transposed_matrix_matrix_product(gene_matrix & A, gene_matrix & B, gene_matrix & X, int N)
  {
    cblas_dgemm(CblasColMajor,CblasTrans,CblasTrans,N,N,N,1.0,A,N,B,N,0.0,X,N);
  }

  static  inline void ata_product(gene_matrix & A, gene_matrix & X, int N)
  {
    cblas_dgemm(CblasColMajor,CblasTrans,CblasNoTrans,N,N,N,1.0,A,N,A,N,0.0,X,N);
  }

  static  inline void aat_product(gene_matrix & A, gene_matrix & X, int N)
  {
    cblas_dgemm(CblasColMajor,CblasNoTrans,CblasTrans,N,N,N,1.0,A,N,A,N,0.0,X,N);
  }

  static  inline void axpy(real coef, const gene_vector & X, gene_vector & Y, int N)
  {
    cblas_daxpy(N,coef,X,1,Y,1);
  }

};

template<>
class C_BLAS_interface<float> : public f77_interface_base<float>
{
public :

  static inline std::string name( void )
  {
    return "C_BLAS";
  }

  static  inline void matrix_vector_product(gene_matrix & A, gene_vector & B, gene_vector & X, int N)
  {
    cblas_sgemv(CblasColMajor,CblasNoTrans,N,N,1.0,A,N,B,1,0.0,X,1);
  }

  static  inline void atv_product(gene_matrix & A, gene_vector & B, gene_vector & X, int N)
  {
    cblas_sgemv(CblasColMajor,CblasTrans,N,N,1.0,A,N,B,1,0.0,X,1);
  }

  static  inline void matrix_matrix_product(gene_matrix & A, gene_matrix & B, gene_matrix & X, int N)
  {
    cblas_sgemm(CblasColMajor,CblasNoTrans,CblasNoTrans,N,N,N,1.0,A,N,B,N,0.0,X,N);
  }

  static  inline void transposed_matrix_matrix_product(gene_matrix & A, gene_matrix & B, gene_matrix & X, int N)
  {
    cblas_sgemm(CblasColMajor,CblasNoTrans,CblasNoTrans,N,N,N,1.0,A,N,B,N,0.0,X,N);
  }

  static  inline void ata_product(gene_matrix & A, gene_matrix & X, int N)
  {
    cblas_sgemm(CblasColMajor,CblasTrans,CblasNoTrans,N,N,N,1.0,A,N,A,N,0.0,X,N);
  }

  static  inline void aat_product(gene_matrix & A, gene_matrix & X, int N)
  {
    cblas_sgemm(CblasColMajor,CblasNoTrans,CblasTrans,N,N,N,1.0,A,N,A,N,0.0,X,N);
  }

  static  inline void axpy(float coef, const gene_vector & X, gene_vector & Y, int N)
  {
    cblas_saxpy(N,coef,X,1,Y,1);
  }

};


#endif