Statistiques
| Révision :

root / src / blas / HPL_dtrsv.c @ 7

Historique | Voir | Annoter | Télécharger (16,95 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_dtrsv
53 1 equemene
54 1 equemene
#ifdef HPL_CALL_VSIPL
55 1 equemene
56 1 equemene
#ifdef STDC_HEADERS
57 1 equemene
static void HPL_dtrsvLNN
58 1 equemene
(
59 1 equemene
   const int                  N,
60 1 equemene
   const double               * A,
61 1 equemene
   const int                  LDA,
62 1 equemene
   double                     * X,
63 1 equemene
   const int                  INCX
64 1 equemene
)
65 1 equemene
#else
66 1 equemene
static void HPL_dtrsvLNN( N, A, LDA, X, INCX )
67 1 equemene
   const int                  INCX, LDA, N;
68 1 equemene
   const double               * A;
69 1 equemene
   double                     * X;
70 1 equemene
#endif
71 1 equemene
{
72 1 equemene
   int                        i, iaij, ix, j, jaj, jx, ldap1 = LDA + 1;
73 1 equemene
   register double            t0;
74 1 equemene
75 1 equemene
   for( j = 0, jaj = 0, jx  = 0; j < N; j++, jaj += ldap1, jx += INCX )
76 1 equemene
   {
77 1 equemene
      X[jx] /= A[jaj]; t0 = X[jx];
78 1 equemene
      for( i = j+1,    iaij  = jaj+1, ix  = jx + INCX;
79 1 equemene
           i < N; i++, iaij += 1,     ix += INCX ) { X[ix] -= t0 * A[iaij]; }
80 1 equemene
   }
81 1 equemene
}
82 1 equemene
83 1 equemene
#ifdef STDC_HEADERS
84 1 equemene
static void HPL_dtrsvLNU
85 1 equemene
(
86 1 equemene
   const int                  N,
87 1 equemene
   const double               * A,
88 1 equemene
   const int                  LDA,
89 1 equemene
   double                     * X,
90 1 equemene
   const int                  INCX
91 1 equemene
)
92 1 equemene
#else
93 1 equemene
static void HPL_dtrsvLNU( N, A, LDA, X, INCX )
94 1 equemene
   const int                  INCX, LDA, N;
95 1 equemene
   const double               * A;
96 1 equemene
   double                     * X;
97 1 equemene
#endif
98 1 equemene
{
99 1 equemene
   int                        i, iaij, ix, j, jaj, jx, ldap1 = LDA + 1;
100 1 equemene
   register double            t0;
101 1 equemene
102 1 equemene
   for( j = 0, jaj = 0, jx = 0; j < N; j++, jaj += ldap1, jx += INCX )
103 1 equemene
   {
104 1 equemene
      t0 = X[jx];
105 1 equemene
      for( i = j+1,    iaij  = jaj+1, ix  = jx + INCX;
106 1 equemene
           i < N; i++, iaij += 1,     ix += INCX ) { X[ix] -= t0 * A[iaij]; }
107 1 equemene
   }
108 1 equemene
}
109 1 equemene
110 1 equemene
#ifdef STDC_HEADERS
111 1 equemene
static void HPL_dtrsvLTN
112 1 equemene
(
113 1 equemene
   const int                  N,
114 1 equemene
   const double               * A,
115 1 equemene
   const int                  LDA,
116 1 equemene
   double                     * X,
117 1 equemene
   const int                  INCX
118 1 equemene
)
119 1 equemene
#else
120 1 equemene
static void HPL_dtrsvLTN( N, A, LDA, X, INCX )
121 1 equemene
   const int                  INCX, LDA, N;
122 1 equemene
   const double               * A;
123 1 equemene
   double                     * X;
124 1 equemene
#endif
125 1 equemene
{
126 1 equemene
   int                        i, iaij, ix, j, jaj, jx, ldap1 = LDA + 1;
127 1 equemene
   register double            t0;
128 1 equemene
129 1 equemene
   for( j = N-1,     jaj  = (N-1)*(ldap1), jx  = (N-1)*INCX;
130 1 equemene
        j >= 0; j--, jaj -= ldap1,         jx -= INCX )
131 1 equemene
   {
132 1 equemene
      t0 = X[jx];
133 1 equemene
      for( i = j+1,    iaij  = 1+jaj, ix  = jx + INCX;
134 1 equemene
           i < N; i++, iaij += 1,     ix += INCX ) { t0 -= A[iaij] * X[ix]; }
135 1 equemene
      t0 /= A[jaj]; X[jx] = t0;
136 1 equemene
   }
137 1 equemene
}
138 1 equemene
139 1 equemene
#ifdef STDC_HEADERS
140 1 equemene
static void HPL_dtrsvLTU
141 1 equemene
(
142 1 equemene
   const int                  N,
143 1 equemene
   const double               * A,
144 1 equemene
   const int                  LDA,
145 1 equemene
   double                     * X,
146 1 equemene
   const int                  INCX
147 1 equemene
)
148 1 equemene
#else
149 1 equemene
static void HPL_dtrsvLTU( N, A, LDA, X, INCX )
150 1 equemene
   const int                  INCX, LDA, N;
151 1 equemene
   const double               * A;
152 1 equemene
   double                     * X;
153 1 equemene
#endif
154 1 equemene
{
155 1 equemene
   int                        i, iaij, ix, j, jaj, jx, ldap1 = LDA + 1;
156 1 equemene
   register double            t0;
157 1 equemene
158 1 equemene
   for( j = N-1,     jaj  = (N-1)*(ldap1), jx  = (N-1)*INCX;
159 1 equemene
        j >= 0; j--, jaj -= ldap1,         jx -= INCX )
160 1 equemene
   {
161 1 equemene
      t0 = X[jx];
162 1 equemene
      for( i = j+1,    iaij  = 1+jaj, ix  = jx + INCX;
163 1 equemene
           i < N; i++, iaij += 1,     ix += INCX ) { t0 -= A[iaij] * X[ix]; }
164 1 equemene
      X[jx] = t0;
165 1 equemene
   }
166 1 equemene
}
167 1 equemene
168 1 equemene
169 1 equemene
#ifdef STDC_HEADERS
170 1 equemene
static void HPL_dtrsvUNN
171 1 equemene
(
172 1 equemene
   const int                  N,
173 1 equemene
   const double               * A,
174 1 equemene
   const int                  LDA,
175 1 equemene
   double                     * X,
176 1 equemene
   const int                  INCX
177 1 equemene
)
178 1 equemene
#else
179 1 equemene
static void HPL_dtrsvUNN( N, A, LDA, X, INCX )
180 1 equemene
   const int                  INCX, LDA, N;
181 1 equemene
   const double               * A;
182 1 equemene
   double                     * X;
183 1 equemene
#endif
184 1 equemene
{
185 1 equemene
   int                        i, iaij, ix, j, jaj, jx;
186 1 equemene
   register double            t0;
187 1 equemene
188 1 equemene
   for( j = N-1,     jaj  = (N-1)*LDA, jx  = (N-1)*INCX;
189 1 equemene
        j >= 0; j--, jaj -= LDA,       jx -= INCX )
190 1 equemene
   {
191 1 equemene
      X[jx] /= A[j+jaj]; t0 = X[jx];
192 1 equemene
      for( i = 0, iaij = jaj, ix = 0; i < j; i++, iaij += 1, ix += INCX )
193 1 equemene
      { X[ix] -= t0 * A[iaij]; }
194 1 equemene
   }
195 1 equemene
}
196 1 equemene
197 1 equemene
198 1 equemene
#ifdef STDC_HEADERS
199 1 equemene
static void HPL_dtrsvUNU
200 1 equemene
(
201 1 equemene
   const int                  N,
202 1 equemene
   const double               * A,
203 1 equemene
   const int                  LDA,
204 1 equemene
   double                     * X,
205 1 equemene
   const int                  INCX
206 1 equemene
)
207 1 equemene
#else
208 1 equemene
static void HPL_dtrsvUNU( N, A, LDA, X, INCX )
209 1 equemene
   const int                  INCX, LDA, N;
210 1 equemene
   const double               * A;
211 1 equemene
   double                     * X;
212 1 equemene
#endif
213 1 equemene
{
214 1 equemene
   int                        i, iaij, ix, j, jaj, jx;
215 1 equemene
   register double            t0;
216 1 equemene
217 1 equemene
   for( j = N-1,     jaj  = (N-1)*LDA, jx  = (N-1)*INCX;
218 1 equemene
        j >= 0; j--, jaj -= LDA,       jx -= INCX )
219 1 equemene
   {
220 1 equemene
      t0 = X[jx];
221 1 equemene
      for( i = 0, iaij = jaj, ix = 0; i < j; i++, iaij += 1, ix += INCX )
222 1 equemene
      { X[ix] -= t0 * A[iaij]; }
223 1 equemene
   }
224 1 equemene
}
225 1 equemene
226 1 equemene
227 1 equemene
#ifdef STDC_HEADERS
228 1 equemene
static void HPL_dtrsvUTN
229 1 equemene
(
230 1 equemene
   const int                  N,
231 1 equemene
   const double               * A,
232 1 equemene
   const int                  LDA,
233 1 equemene
   double                     * X,
234 1 equemene
   const int                  INCX
235 1 equemene
)
236 1 equemene
#else
237 1 equemene
static void HPL_dtrsvUTN( N, A, LDA, X, INCX )
238 1 equemene
   const int                  INCX, LDA, N;
239 1 equemene
   const double               * A;
240 1 equemene
   double                     * X;
241 1 equemene
#endif
242 1 equemene
{
243 1 equemene
   int                        i, iaij, ix, j, jaj, jx;
244 1 equemene
   register double            t0;
245 1 equemene
246 1 equemene
   for( j = 0, jaj = 0,jx = 0; j < N; j++, jaj += LDA, jx += INCX )
247 1 equemene
   {
248 1 equemene
      t0 = X[jx];
249 1 equemene
      for( i = 0, iaij = jaj, ix = 0; i < j; i++, iaij += 1, ix += INCX )
250 1 equemene
      { t0 -= A[iaij] * X[ix]; }
251 1 equemene
      t0 /= A[iaij]; X[jx] = t0;
252 1 equemene
   }
253 1 equemene
}
254 1 equemene
255 1 equemene
#ifdef STDC_HEADERS
256 1 equemene
static void HPL_dtrsvUTU
257 1 equemene
(
258 1 equemene
   const int                  N,
259 1 equemene
   const double               * A,
260 1 equemene
   const int                  LDA,
261 1 equemene
   double                     * X,
262 1 equemene
   const int                  INCX
263 1 equemene
)
264 1 equemene
#else
265 1 equemene
static void HPL_dtrsvUTU( N, A, LDA, X, INCX )
266 1 equemene
   const int                  INCX, LDA, N;
267 1 equemene
   const double               * A;
268 1 equemene
   double                     * X;
269 1 equemene
#endif
270 1 equemene
{
271 1 equemene
   int                        i, iaij, ix, j, jaj, jx;
272 1 equemene
   register double            t0;
273 1 equemene
274 1 equemene
   for( j = 0, jaj = 0, jx = 0; j < N; j++, jaj += LDA, jx += INCX )
275 1 equemene
   {
276 1 equemene
      t0 = X[jx];
277 1 equemene
      for( i = 0, iaij = jaj, ix = 0; i < j; i++, iaij += 1, ix += INCX )
278 1 equemene
      { t0 -= A[iaij] * X[ix]; }
279 1 equemene
      X[jx] = t0;
280 1 equemene
   }
281 1 equemene
}
282 1 equemene
283 1 equemene
#ifdef STDC_HEADERS
284 1 equemene
static void HPL_dtrsv0
285 1 equemene
(
286 1 equemene
   const enum HPL_UPLO        UPLO,
287 1 equemene
   const enum HPL_TRANS       TRANS,
288 1 equemene
   const enum HPL_DIAG        DIAG,
289 1 equemene
   const int                  N,
290 1 equemene
   const double               * A,
291 1 equemene
   const int                  LDA,
292 1 equemene
   double                     * X,
293 1 equemene
   const int                  INCX
294 1 equemene
)
295 1 equemene
#else
296 1 equemene
static void HPL_dtrsv0( UPLO, TRANS, DIAG, N, A, LDA, X, INCX )
297 1 equemene
   const enum HPL_UPLO        UPLO;
298 1 equemene
   const enum HPL_TRANS       TRANS;
299 1 equemene
   const enum HPL_DIAG        DIAG;
300 1 equemene
   const int                  INCX, LDA, N;
301 1 equemene
   const double               * A;
302 1 equemene
   double                     * X;
303 1 equemene
#endif
304 1 equemene
{
305 1 equemene
   if( N == 0 ) return;
306 1 equemene
307 1 equemene
   if( UPLO == HplUpper )
308 1 equemene
   {
309 1 equemene
      if( TRANS == HplNoTrans )
310 1 equemene
      {
311 1 equemene
         if( DIAG == HplNonUnit ) { HPL_dtrsvUNN( N,    A, LDA, X, INCX ); }
312 1 equemene
         else                     { HPL_dtrsvUNU( N,    A, LDA, X, INCX ); }
313 1 equemene
      }
314 1 equemene
      else
315 1 equemene
      {
316 1 equemene
         if( DIAG == HplNonUnit ) { HPL_dtrsvUTN( N,    A, LDA, X, INCX ); }
317 1 equemene
         else                     { HPL_dtrsvUTU( N,    A, LDA, X, INCX ); }
318 1 equemene
      }
319 1 equemene
   }
320 1 equemene
   else
321 1 equemene
   {
322 1 equemene
      if( TRANS == HplNoTrans )
323 1 equemene
      {
324 1 equemene
         if( DIAG == HplNonUnit ) { HPL_dtrsvLNN( N,    A, LDA, X, INCX ); }
325 1 equemene
         else                     { HPL_dtrsvLNU( N,    A, LDA, X, INCX ); }
326 1 equemene
      }
327 1 equemene
      else
328 1 equemene
      {
329 1 equemene
         if( DIAG == HplNonUnit ) { HPL_dtrsvLTN( N,    A, LDA, X, INCX ); }
330 1 equemene
         else                     { HPL_dtrsvLTU( N,    A, LDA, X, INCX ); }
331 1 equemene
      }
332 1 equemene
   }
333 1 equemene
}
334 1 equemene
335 1 equemene
#endif
336 1 equemene
337 1 equemene
#ifdef STDC_HEADERS
338 1 equemene
void HPL_dtrsv
339 1 equemene
(
340 1 equemene
   const enum HPL_ORDER             ORDER,
341 1 equemene
   const enum HPL_UPLO              UPLO,
342 1 equemene
   const enum HPL_TRANS             TRANS,
343 1 equemene
   const enum HPL_DIAG              DIAG,
344 1 equemene
   const int                        N,
345 1 equemene
   const double *                   A,
346 1 equemene
   const int                        LDA,
347 1 equemene
   double *                         X,
348 1 equemene
   const int                        INCX
349 1 equemene
)
350 1 equemene
#else
351 1 equemene
void HPL_dtrsv
352 1 equemene
( ORDER, UPLO, TRANS, DIAG, N, A, LDA, X, INCX )
353 1 equemene
   const enum HPL_ORDER             ORDER;
354 1 equemene
   const enum HPL_UPLO              UPLO;
355 1 equemene
   const enum HPL_TRANS             TRANS;
356 1 equemene
   const enum HPL_DIAG              DIAG;
357 1 equemene
   const int                        N;
358 1 equemene
   const double *                   A;
359 1 equemene
   const int                        LDA;
360 1 equemene
   double *                         X;
361 1 equemene
   const int                        INCX;
362 1 equemene
#endif
363 1 equemene
{
364 1 equemene
/*
365 1 equemene
 * Purpose
366 1 equemene
 * =======
367 1 equemene
 *
368 1 equemene
 * HPL_dtrsv solves one of the systems of equations
369 1 equemene
 *
370 1 equemene
 *     A * x = b,   or   A^T * x = b,
371 1 equemene
 *
372 1 equemene
 * where b and x are n-element vectors and  A  is an n by n non-unit, or
373 1 equemene
 * unit, upper or lower triangular matrix.
374 1 equemene
 *
375 1 equemene
 * No test for  singularity  or  near-singularity  is included  in  this
376 1 equemene
 * routine. Such tests must be performed before calling this routine.
377 1 equemene
 *
378 1 equemene
 * Arguments
379 1 equemene
 * =========
380 1 equemene
 *
381 1 equemene
 * ORDER   (local input)                 const enum HPL_ORDER
382 1 equemene
 *         On entry, ORDER  specifies the storage format of the operands
383 1 equemene
 *         as follows:
384 1 equemene
 *            ORDER = HplRowMajor,
385 1 equemene
 *            ORDER = HplColumnMajor.
386 1 equemene
 *
387 1 equemene
 * UPLO    (local input)                 const enum HPL_UPLO
388 1 equemene
 *         On  entry,   UPLO   specifies  whether  the  upper  or  lower
389 1 equemene
 *         triangular  part  of the array  A  is to be referenced.  When
390 1 equemene
 *         UPLO==HplUpper, only  the upper triangular part of A is to be
391 1 equemene
 *         referenced, otherwise only the lower triangular part of A is
392 1 equemene
 *         to be referenced.
393 1 equemene
 *
394 1 equemene
 * TRANS   (local input)                 const enum HPL_TRANS
395 1 equemene
 *         On entry,  TRANS  specifies  the equations  to  be  solved as
396 1 equemene
 *         follows:
397 1 equemene
 *            TRANS==HplNoTrans     A   * x = b,
398 1 equemene
 *            TRANS==HplTrans       A^T * x = b.
399 1 equemene
 *
400 1 equemene
 * DIAG    (local input)                 const enum HPL_DIAG
401 1 equemene
 *         On entry,  DIAG  specifies  whether  A  is unit triangular or
402 1 equemene
 *         not. When DIAG==HplUnit,  A is assumed to be unit triangular,
403 1 equemene
 *         and otherwise, A is not assumed to be unit triangular.
404 1 equemene
 *
405 1 equemene
 * N       (local input)                 const int
406 1 equemene
 *         On entry, N specifies the order of the matrix A. N must be at
407 1 equemene
 *         least zero.
408 1 equemene
 *
409 1 equemene
 * A       (local input)                 const double *
410 1 equemene
 *         On entry,  A  points  to an array of size equal to or greater
411 1 equemene
 *         than LDA * n. Before entry with  UPLO==HplUpper,  the leading
412 1 equemene
 *         n by n upper triangular  part of the array A must contain the
413 1 equemene
 *         upper triangular  matrix and the  strictly  lower  triangular
414 1 equemene
 *         part of A is not referenced.  When  UPLO==HplLower  on entry,
415 1 equemene
 *         the  leading n by n lower triangular part of the array A must
416 1 equemene
 *         contain the lower triangular matrix  and  the  strictly upper
417 1 equemene
 *         triangular part of A is not referenced.
418 1 equemene
 *
419 1 equemene
 *         Note  that  when  DIAG==HplUnit,  the diagonal elements of  A
420 1 equemene
 *         not referenced  either,  but are assumed to be unity.
421 1 equemene
 *
422 1 equemene
 * LDA     (local input)                 const int
423 1 equemene
 *         On entry,  LDA  specifies  the  leading  dimension  of  A  as
424 1 equemene
 *         declared  in  the  calling  (sub) program.  LDA  must  be  at
425 1 equemene
 *         least MAX(1,n).
426 1 equemene
 *
427 1 equemene
 * X       (local input/output)          double *
428 1 equemene
 *         On entry,  X  is an incremented array of dimension  at  least
429 1 equemene
 *         ( 1 + ( n - 1 ) * abs( INCX ) )  that  contains the vector x.
430 1 equemene
 *         Before entry,  the  incremented array  X  must contain  the n
431 1 equemene
 *         element right-hand side vector b. On exit,  X  is overwritten
432 1 equemene
 *         with the solution vector x.
433 1 equemene
 *
434 1 equemene
 * INCX    (local input)                 const int
435 1 equemene
 *         On entry, INCX specifies the increment for the elements of X.
436 1 equemene
 *         INCX must not be zero.
437 1 equemene
 *
438 1 equemene
 * ---------------------------------------------------------------------
439 1 equemene
 */
440 1 equemene
#ifdef HPL_CALL_CBLAS
441 1 equemene
   cblas_dtrsv( ORDER, UPLO, TRANS, DIAG, N, A, LDA, X, INCX );
442 1 equemene
#endif
443 1 equemene
#ifdef HPL_CALL_VSIPL
444 1 equemene
   if( ORDER == HplColumnMajor )
445 1 equemene
   {
446 1 equemene
      HPL_dtrsv0( UPLO, TRANS, DIAG, N, A, LDA, X, INCX );
447 1 equemene
   }
448 1 equemene
   else
449 1 equemene
   {
450 1 equemene
      HPL_dtrsv0( ( UPLO  == HplUpper   ? HplLower : HplUpper   ),
451 1 equemene
                  ( TRANS == HplNoTrans ? HplTrans : HplNoTrans ),
452 1 equemene
                  DIAG, N, A, LDA, X, INCX );
453 1 equemene
   }
454 1 equemene
#endif
455 1 equemene
#ifdef HPL_CALL_FBLAS
456 1 equemene
#ifdef StringSunStyle
457 1 equemene
#ifdef HPL_USE_F77_INTEGER_DEF
458 1 equemene
   F77_INTEGER               IONE = 1;
459 1 equemene
#else
460 1 equemene
   int                       IONE = 1;
461 1 equemene
#endif
462 1 equemene
#endif
463 1 equemene
#ifdef StringStructVal
464 1 equemene
   F77_CHAR                  fuplo, ftran, fdiag;
465 1 equemene
#endif
466 1 equemene
#ifdef StringStructPtr
467 1 equemene
   F77_CHAR                  fuplo, ftran, fdiag;
468 1 equemene
#endif
469 1 equemene
#ifdef StringCrayStyle
470 1 equemene
   F77_CHAR                  fuplo, ftran, fdiag;
471 1 equemene
#endif
472 1 equemene
473 1 equemene
#ifdef HPL_USE_F77_INTEGER_DEF
474 1 equemene
   const F77_INTEGER         F77N = N, F77lda = LDA, F77incx = INCX;
475 1 equemene
#else
476 1 equemene
#define F77N              N
477 1 equemene
#define F77lda            LDA
478 1 equemene
#define F77incx           INCX
479 1 equemene
#endif
480 1 equemene
   char                      cuplo, ctran, cdiag;
481 1 equemene
482 1 equemene
   if( ORDER == HplColumnMajor )
483 1 equemene
   {
484 1 equemene
      cuplo = ( UPLO  == HplUpper   ? 'U' : 'L' );
485 1 equemene
      ctran = ( TRANS == HplNoTrans ? 'N' : 'T' );
486 1 equemene
   }
487 1 equemene
   else
488 1 equemene
   {
489 1 equemene
      cuplo = ( UPLO  == HplUpper   ? 'L' : 'U' );
490 1 equemene
      ctran = ( TRANS == HplNoTrans ? 'T' : 'N' );
491 1 equemene
   }
492 1 equemene
   cdiag = ( DIAG == HplNonUnit ? 'N' : 'U' );
493 1 equemene
494 1 equemene
#ifdef StringSunStyle
495 1 equemene
   F77dtrsv( &cuplo, &ctran, &cdiag, &F77N, A, &F77lda, X, &F77incx,
496 1 equemene
             IONE, IONE, IONE );
497 1 equemene
#endif
498 1 equemene
#ifdef StringCrayStyle
499 1 equemene
   ftran = HPL_C2F_CHAR( ctran ); fdiag = HPL_C2F_CHAR( cdiag );
500 1 equemene
   fuplo = HPL_C2F_CHAR( cuplo );
501 1 equemene
   F77dtrsv( fuplo,  ftran,  fdiag,  &F77N, A, &F77lda, X, &F77incx );
502 1 equemene
#endif
503 1 equemene
#ifdef StringStructVal
504 1 equemene
   fuplo.len = 1; fuplo.cp = &cuplo; ftran.len = 1; ftran.cp = &ctran;
505 1 equemene
   fdiag.len = 1; fdiag.cp = &cdiag;
506 1 equemene
   F77dtrsv( fuplo,  ftran,  fdiag,  &F77N, A, &F77lda, X, &F77incx );
507 1 equemene
#endif
508 1 equemene
#ifdef StringStructPtr
509 1 equemene
   fuplo.len = 1; fuplo.cp = &cuplo; ftran.len = 1; ftran.cp = &ctran;
510 1 equemene
   fdiag.len = 1; fdiag.cp = &cdiag;
511 1 equemene
   F77dtrsv( &fuplo, &ftran, &fdiag, &F77N, A, &F77lda, X, &F77incx );
512 1 equemene
#endif
513 1 equemene
514 1 equemene
#endif
515 1 equemene
/*
516 1 equemene
 * End of HPL_dtrsv
517 1 equemene
 */
518 1 equemene
}
519 1 equemene
520 1 equemene
#endif