File: sdsdot.f

package info (click to toggle)
lapack 3.12.1-6
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 78,912 kB
  • sloc: fortran: 622,840; ansic: 217,704; f90: 6,041; makefile: 5,100; sh: 326; python: 270; xml: 182
file content (163 lines) | stat: -rw-r--r-- 4,113 bytes parent folder | download
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
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
*> \brief \b SDSDOT
*
*  =========== DOCUMENTATION ===========
*
* Online html documentation available at
*            http://www.netlib.org/lapack/explore-html/
*
*  Definition:
*  ===========
*
*       REAL FUNCTION SDSDOT(N,SB,SX,INCX,SY,INCY)
*
*       .. Scalar Arguments ..
*       REAL SB
*       INTEGER INCX,INCY,N
*       ..
*       .. Array Arguments ..
*       REAL SX(*),SY(*)
*       ..
*
*> \par Purpose:
*  =============
*>
*> \verbatim
*>
*>   Compute the inner product of two vectors with extended
*>   precision accumulation.
*>
*>   Returns S.P. result with dot product accumulated in D.P.
*>   SDSDOT = SB + sum for I = 0 to N-1 of SX(LX+I*INCX)*SY(LY+I*INCY),
*>   where LX = 1 if INCX .GE. 0, else LX = 1+(1-N)*INCX, and LY is
*>   defined in a similar way using INCY.
*> \endverbatim
*
*  Arguments:
*  ==========
*
*> \param[in] N
*> \verbatim
*>          N is INTEGER
*>          number of elements in input vector(s)
*> \endverbatim
*>
*> \param[in] SB
*> \verbatim
*>          SB is REAL
*>          single precision scalar to be added to inner product
*> \endverbatim
*>
*> \param[in] SX
*> \verbatim
*>          SX is REAL array, dimension ( 1 + ( N - 1 )*abs( INCX ) )
*>          single precision vector with N elements
*> \endverbatim
*>
*> \param[in] INCX
*> \verbatim
*>          INCX is INTEGER
*>          storage spacing between elements of SX
*> \endverbatim
*>
*> \param[in] SY
*> \verbatim
*>          SY is REAL array, dimension ( 1 + ( N - 1 )*abs( INCX ) )
*>          single precision vector with N elements
*> \endverbatim
*>
*> \param[in] INCY
*> \verbatim
*>          INCY is INTEGER
*>          storage spacing between elements of SY
*> \endverbatim
*
*  Authors:
*  ========
*
*> \author Lawson, C. L., (JPL), Hanson, R. J., (SNLA),
*> \author Kincaid, D. R., (U. of Texas), Krogh, F. T., (JPL)
*
*> \author Univ. of Tennessee
*> \author Univ. of California Berkeley
*> \author Univ. of Colorado Denver
*> \author NAG Ltd.
*
*> \ingroup dot
*
*> \par Further Details:
*  =====================
*>
*> \verbatim
*>
*>    REFERENCES
*>
*>    C. L. Lawson, R. J. Hanson, D. R. Kincaid and F. T.
*>    Krogh, Basic linear algebra subprograms for Fortran
*>    usage, Algorithm No. 539, Transactions on Mathematical
*>    Software 5, 3 (September 1979), pp. 308-323.
*>
*>    REVISION HISTORY  (YYMMDD)
*>
*>    791001  DATE WRITTEN
*>    890531  Changed all specific intrinsics to generic.  (WRB)
*>    890831  Modified array declarations.  (WRB)
*>    890831  REVISION DATE from Version 3.2
*>    891214  Prologue converted to Version 4.0 format.  (BAB)
*>    920310  Corrected definition of LX in DESCRIPTION.  (WRB)
*>    920501  Reformatted the REFERENCES section.  (WRB)
*>    070118  Reformat to LAPACK coding style
*> \endverbatim
*>
*  =====================================================================
      REAL FUNCTION SDSDOT(N,SB,SX,INCX,SY,INCY)
*
*  -- Reference BLAS level1 routine --
*  -- Reference BLAS is a software package provided by Univ. of Tennessee,    --
*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
*
*     .. Scalar Arguments ..
      REAL SB
      INTEGER INCX,INCY,N
*     ..
*     .. Array Arguments ..
      REAL SX(*),SY(*)
*     .. Local Scalars ..
      DOUBLE PRECISION DSDOT
      INTEGER I,KX,KY,NS
*     ..
*     .. Intrinsic Functions ..
      INTRINSIC DBLE
*     ..
      DSDOT = SB
      IF (N.LE.0) THEN
         SDSDOT = REAL(DSDOT)
         RETURN
      END IF
      IF (INCX.EQ.INCY .AND. INCX.GT.0) THEN
*
*     Code for equal and positive increments.
*
         NS = N*INCX
         DO I = 1,NS,INCX
            DSDOT = DSDOT + DBLE(SX(I))*DBLE(SY(I))
         END DO
      ELSE
*
*     Code for unequal or nonpositive increments.
*
         KX = 1
         KY = 1
         IF (INCX.LT.0) KX = 1 + (1-N)*INCX
         IF (INCY.LT.0) KY = 1 + (1-N)*INCY
         DO I = 1,N
            DSDOT = DSDOT + DBLE(SX(KX))*DBLE(SY(KY))
            KX = KX + INCX
            KY = KY + INCY
         END DO
      END IF
      SDSDOT = REAL(DSDOT)
      RETURN
*
*     End of SDSDOT
*
      END