Statistiques
| Révision :

root / src / blas / HPL_dger.c @ 9

Historique | Voir | Annoter | Télécharger (8,45 ko)

1 1 equemene
/*
2 1 equemene
 * -- High Performance Computing Linpack Benchmark (HPL)
3 1 equemene
 *    HPL - 2.0 - September 10, 2008
4 1 equemene
 *    Antoine P. Petitet
5 1 equemene
 *    University of Tennessee, Knoxville
6 1 equemene
 *    Innovative Computing Laboratory
7 1 equemene
 *    (C) Copyright 2000-2008 All Rights Reserved
8 1 equemene
 *
9 1 equemene
 * -- Copyright notice and Licensing terms:
10 1 equemene
 *
11 1 equemene
 * Redistribution  and  use in  source and binary forms, with or without
12 1 equemene
 * modification, are  permitted provided  that the following  conditions
13 1 equemene
 * are met:
14 1 equemene
 *
15 1 equemene
 * 1. Redistributions  of  source  code  must retain the above copyright
16 1 equemene
 * notice, this list of conditions and the following disclaimer.
17 1 equemene
 *
18 1 equemene
 * 2. Redistributions in binary form must reproduce  the above copyright
19 1 equemene
 * notice, this list of conditions,  and the following disclaimer in the
20 1 equemene
 * documentation and/or other materials provided with the distribution.
21 1 equemene
 *
22 1 equemene
 * 3. All  advertising  materials  mentioning  features  or  use of this
23 1 equemene
 * software must display the following acknowledgement:
24 1 equemene
 * This  product  includes  software  developed  at  the  University  of
25 1 equemene
 * Tennessee, Knoxville, Innovative Computing Laboratory.
26 1 equemene
 *
27 1 equemene
 * 4. The name of the  University,  the name of the  Laboratory,  or the
28 1 equemene
 * names  of  its  contributors  may  not  be used to endorse or promote
29 1 equemene
 * products  derived   from   this  software  without  specific  written
30 1 equemene
 * permission.
31 1 equemene
 *
32 1 equemene
 * -- Disclaimer:
33 1 equemene
 *
34 1 equemene
 * THIS  SOFTWARE  IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
35 1 equemene
 * ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES,  INCLUDING,  BUT NOT
36 1 equemene
 * LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
37 1 equemene
 * A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE UNIVERSITY
38 1 equemene
 * OR  CONTRIBUTORS  BE  LIABLE FOR ANY  DIRECT,  INDIRECT,  INCIDENTAL,
39 1 equemene
 * SPECIAL,  EXEMPLARY,  OR  CONSEQUENTIAL DAMAGES  (INCLUDING,  BUT NOT
40 1 equemene
 * LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
41 1 equemene
 * DATA OR PROFITS; OR BUSINESS INTERRUPTION)  HOWEVER CAUSED AND ON ANY
42 1 equemene
 * THEORY OF LIABILITY, WHETHER IN CONTRACT,  STRICT LIABILITY,  OR TORT
43 1 equemene
 * (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
44 1 equemene
 * OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
45 1 equemene
 * ---------------------------------------------------------------------
46 1 equemene
 */
47 1 equemene
/*
48 1 equemene
 * Include files
49 1 equemene
 */
50 1 equemene
#include "hpl.h"
51 1 equemene
52 1 equemene
#ifndef HPL_dger
53 1 equemene
54 1 equemene
#ifdef STDC_HEADERS
55 1 equemene
void HPL_dger
56 1 equemene
(
57 1 equemene
   const enum HPL_ORDER             ORDER,
58 1 equemene
   const int                        M,
59 1 equemene
   const int                        N,
60 1 equemene
   const double                     ALPHA,
61 1 equemene
   const double *                   X,
62 1 equemene
   const int                        INCX,
63 1 equemene
   double *                         Y,
64 1 equemene
   const int                        INCY,
65 1 equemene
   double *                         A,
66 1 equemene
   const int                        LDA
67 1 equemene
)
68 1 equemene
#else
69 1 equemene
void HPL_dger
70 1 equemene
( ORDER, M, N, ALPHA, X, INCX, Y, INCY, A, LDA )
71 1 equemene
   const enum HPL_ORDER             ORDER;
72 1 equemene
   const int                        M;
73 1 equemene
   const int                        N;
74 1 equemene
   const double                     ALPHA;
75 1 equemene
   const double *                   X;
76 1 equemene
   const int                        INCX;
77 1 equemene
   double *                         Y;
78 1 equemene
   const int                        INCY;
79 1 equemene
   double *                         A;
80 1 equemene
   const int                        LDA;
81 1 equemene
#endif
82 1 equemene
{
83 1 equemene
/*
84 1 equemene
 * Purpose
85 1 equemene
 * =======
86 1 equemene
 *
87 1 equemene
 * HPL_dger performs the rank 1 operation
88 1 equemene
 *
89 1 equemene
 *     A := alpha * x * y^T + A,
90 1 equemene
 *
91 1 equemene
 * where alpha is a scalar,  x is an m-element vector, y is an n-element
92 1 equemene
 * vector and A is an m by n matrix.
93 1 equemene
 *
94 1 equemene
 * Arguments
95 1 equemene
 * =========
96 1 equemene
 *
97 1 equemene
 * ORDER   (local input)                 const enum HPL_ORDER
98 1 equemene
 *         On entry, ORDER  specifies the storage format of the operands
99 1 equemene
 *         as follows:
100 1 equemene
 *            ORDER = HplRowMajor,
101 1 equemene
 *            ORDER = HplColumnMajor.
102 1 equemene
 *
103 1 equemene
 * M       (local input)                 const int
104 1 equemene
 *         On entry,  M  specifies  the number of rows of  the matrix A.
105 1 equemene
 *         M must be at least zero.
106 1 equemene
 *
107 1 equemene
 * N       (local input)                 const int
108 1 equemene
 *         On entry, N  specifies the number of columns of the matrix A.
109 1 equemene
 *         N must be at least zero.
110 1 equemene
 *
111 1 equemene
 * ALPHA   (local input)                 const double
112 1 equemene
 *         On entry, ALPHA specifies the scalar alpha.   When  ALPHA  is
113 1 equemene
 *         supplied as zero then  X and Y  need not be set on input.
114 1 equemene
 *
115 1 equemene
 * X       (local input)                 const double *
116 1 equemene
 *         On entry,  X  is an incremented array of dimension  at  least
117 1 equemene
 *         ( 1 + ( m - 1 ) * abs( INCX ) )  that  contains the vector x.
118 1 equemene
 *
119 1 equemene
 * INCX    (local input)                 const int
120 1 equemene
 *         On entry, INCX specifies the increment for the elements of X.
121 1 equemene
 *         INCX must not be zero.
122 1 equemene
 *
123 1 equemene
 * Y       (local input)                 double *
124 1 equemene
 *         On entry,  Y  is an incremented array of dimension  at  least
125 1 equemene
 *         ( 1 + ( n - 1 ) * abs( INCY ) )  that  contains the vector y.
126 1 equemene
 *
127 1 equemene
 * INCY    (local input)                 const int
128 1 equemene
 *         On entry, INCY specifies the increment for the elements of Y.
129 1 equemene
 *         INCY must not be zero.
130 1 equemene
 *
131 1 equemene
 * A       (local input/output)          double *
132 1 equemene
 *         On entry,  A  points  to an array of size equal to or greater
133 1 equemene
 *         than LDA * n.  Before  entry, the leading m by n part  of the
134 1 equemene
 *         array  A  must contain the matrix coefficients. On exit, A is
135 1 equemene
 *         overwritten by the updated matrix.
136 1 equemene
 *
137 1 equemene
 * LDA     (local input)                 const int
138 1 equemene
 *         On entry,  LDA  specifies  the  leading  dimension  of  A  as
139 1 equemene
 *         declared  in  the  calling  (sub) program.  LDA  must  be  at
140 1 equemene
 *         least MAX(1,m).
141 1 equemene
 *
142 1 equemene
 * ---------------------------------------------------------------------
143 1 equemene
 */
144 1 equemene
#ifdef HPL_CALL_CBLAS
145 1 equemene
   cblas_dger( ORDER, M, N, ALPHA, X, INCX, Y, INCY, A, LDA );
146 1 equemene
#endif
147 1 equemene
#ifdef HPL_CALL_VSIPL
148 1 equemene
   register double           t0;
149 1 equemene
   int                       i, iaij, ix, iy, j, jaj, jx, jy;
150 1 equemene
151 1 equemene
   if( ( M == 0 ) || ( N == 0 ) || ( ALPHA == HPL_rzero ) ) return;
152 1 equemene
153 1 equemene
   if( ORDER == HplColumnMajor )
154 1 equemene
   {
155 1 equemene
      for( j = 0, jaj = 0, jy = 0; j < N; j++, jaj += LDA, jy += INCY )
156 1 equemene
      {
157 1 equemene
         t0 = ALPHA * Y[jy];
158 1 equemene
         for( i = 0, iaij = jaj, ix = 0; i < M; i++, iaij += 1, ix += INCX )
159 1 equemene
         { A[iaij] += X[ix] * t0; }
160 1 equemene
      }
161 1 equemene
   }
162 1 equemene
   else
163 1 equemene
   {
164 1 equemene
      for( j = 0, jaj = 0, jx = 0; j < M; j++, jaj += LDA, jx += INCX )
165 1 equemene
      {
166 1 equemene
         t0 = ALPHA * X[jx];
167 1 equemene
         for( i = 0, iaij = jaj, iy = 0; i < N; i++, iaij += 1, iy += INCY )
168 1 equemene
         { A[iaij] += Y[iy] * t0; }
169 1 equemene
      }
170 1 equemene
   }
171 1 equemene
#endif
172 9 equemene
173 1 equemene
#ifdef HPL_CALL_FBLAS
174 1 equemene
   double                    alpha = ALPHA;
175 1 equemene
#ifdef HPL_USE_F77_INTEGER_DEF
176 1 equemene
   const F77_INTEGER         F77M    = M,   F77N    = N,
177 1 equemene
                             F77lda  = LDA, F77incx = INCX, F77incy = INCY;
178 1 equemene
#else
179 1 equemene
#define F77M                 M
180 1 equemene
#define F77N                 N
181 1 equemene
#define F77lda               LDA
182 1 equemene
#define F77incx              INCX
183 1 equemene
#define F77incy              INCY
184 1 equemene
#endif
185 1 equemene
186 1 equemene
   if( ORDER == HplColumnMajor )
187 1 equemene
   {  F77dger( &F77M, &F77N, &alpha, X, &F77incx, Y, &F77incy, A, &F77lda ); }
188 1 equemene
   else
189 1 equemene
   {  F77dger( &F77N, &F77M, &alpha, Y, &F77incy, X, &F77incx, A, &F77lda ); }
190 1 equemene
#endif
191 9 equemene
192 9 equemene
#ifdef HPL_CALL_CUBLAS
193 9 equemene
   double                    alpha = ALPHA;
194 9 equemene
195 9 equemene
#define CUBLASM                 M
196 9 equemene
#define CUBLASN                 N
197 9 equemene
#define CUBLASlda               LDA
198 9 equemene
#define CUBLASincx              INCX
199 9 equemene
#define CUBLASincy              INCY
200 9 equemene
201 9 equemene
   if( ORDER == HplColumnMajor )
202 9 equemene
   {
203 9 equemene
     CUBLAS_DGER( &CUBLASM, &CUBLASN, &alpha,
204 9 equemene
                  X, &CUBLASincx, Y, &CUBLASincy, A, &CUBLASlda );
205 9 equemene
   }
206 9 equemene
   else
207 9 equemene
   {
208 9 equemene
     CUBLAS_DGER( &CUBLASN, &CUBLASM, &alpha,
209 9 equemene
                  Y, &CUBLASincy, X, &CUBLASincx, A, &CUBLASlda ); }
210 9 equemene
#endif
211 1 equemene
/*
212 1 equemene
 * End of HPL_dger
213 1 equemene
 */
214 1 equemene
}
215 1 equemene
216 1 equemene
#endif