Skip to content

Commit c744ebe

Browse files
committed
updating dlarf1l to use firstv scanner properly
1 parent 741907c commit c744ebe

File tree

1 file changed

+19
-13
lines changed

1 file changed

+19
-13
lines changed

SRC/dlarf1l.f

Lines changed: 19 additions & 13 deletions
Original file line numberDiff line numberDiff line change
@@ -147,7 +147,7 @@ SUBROUTINE DLARF1L( SIDE, M, N, V, INCV, TAU, C, LDC, WORK )
147147
* ..
148148
* .. Local Scalars ..
149149
LOGICAL APPLYLEFT
150-
INTEGER I, LASTV, LASTC, J
150+
INTEGER I, LASTV, LASTC, J, FIRSTV
151151
* ..
152152
* .. External Subroutines ..
153153
EXTERNAL DGEMV, DGER
@@ -160,7 +160,7 @@ SUBROUTINE DLARF1L( SIDE, M, N, V, INCV, TAU, C, LDC, WORK )
160160
* .. Executable Statements ..
161161
*
162162
APPLYLEFT = LSAME( SIDE, 'L' )
163-
LASTV = 0
163+
FIRSTV = 1
164164
LASTC = 0
165165
IF( TAU.NE.ZERO ) THEN
166166
! Set up variables for scanning V. LASTV begins pointing to the end
@@ -170,7 +170,12 @@ SUBROUTINE DLARF1L( SIDE, M, N, V, INCV, TAU, C, LDC, WORK )
170170
ELSE
171171
LASTV = N
172172
END IF
173+
I = 1
173174
! Look for the last non-zero row in V.
175+
DO WHILE( LASTV.GT.FIRSTV .AND. V( I ).EQ.ZERO )
176+
FIRSTV = FIRSTV + 1
177+
I = I + INCV
178+
END DO
174179
IF( APPLYLEFT ) THEN
175180
! Scan for the last non-zero column in C(1:lastv,:).
176181
LASTC = ILADLC(LASTV, N, C, LDC)
@@ -190,15 +195,16 @@ SUBROUTINE DLARF1L( SIDE, M, N, V, INCV, TAU, C, LDC, WORK )
190195
IF( LASTV.GT.0 ) THEN
191196
! Check if m = 1. This means v = 1, So we just need to compute
192197
! C := HC = (1-\tau)C.
193-
IF( LASTV.EQ.1 ) THEN
194-
CALL DSCAL(LASTC, ONE - TAU, C, LDC)
198+
IF( LASTV.EQ.FIRSTV ) THEN
199+
CALL DSCAL(LASTC, ONE - TAU, C( FIRSTV, 1), LDC)
195200
ELSE
196201
*
197202
* w(1:lastc,1) := C(1:lastv,1:lastc)**T * v(1:lastv,1)
198203
*
199204
! w(1:lastc,1) := C(1:lastv-1,1:lastc)**T * v(1:lastv-1,1)
200-
CALL DGEMV( 'Transpose', LASTV-1, LASTC, ONE, C(1,1),
201-
$ LDC, V(1), INCV, ZERO, WORK, 1)
205+
CALL DGEMV( 'Transpose', LASTV-FIRSTV, LASTC, ONE,
206+
$ C(FIRSTV,1), LDC, V(I), INCV, ZERO,
207+
$ WORK, 1)
202208
! w(1:lastc,1) += C(lastv,1:lastc)**T * v(lastv,1) = C(lastv,1:lastc)**T
203209
CALL DAXPY(LASTC, ONE, C(LASTV,1), LDC, WORK, 1)
204210
*
@@ -208,8 +214,8 @@ SUBROUTINE DLARF1L( SIDE, M, N, V, INCV, TAU, C, LDC, WORK )
208214
! = C(...) - tau * w(1:lastc,1)**T
209215
CALL DAXPY(LASTC, -TAU, WORK, 1, C(LASTV,1), LDC)
210216
! C(1:lastv-1,1:lastc) := C(...) - tau * v(1:lastv-1,1)*w(1:lastc,1)**T
211-
CALL DGER(LASTV-1, LASTC, -TAU, V(1), INCV, WORK, 1,
212-
$ C(1,1), LDC)
217+
CALL DGER(LASTV-FIRSTV, LASTC, -TAU, V(I), INCV,
218+
$ WORK, 1, C(FIRSTV,1), LDC)
213219
END IF
214220
END IF
215221
ELSE
@@ -219,15 +225,15 @@ SUBROUTINE DLARF1L( SIDE, M, N, V, INCV, TAU, C, LDC, WORK )
219225
IF( LASTV.GT.0 ) THEN
220226
! Check if n = 1. This means v = 1, so we just need to compute
221227
! C := CH = C(1-\tau).
222-
IF( LASTV.EQ.1 ) THEN
228+
IF( LASTV.EQ.FIRSTV ) THEN
223229
CALL DSCAL(LASTC, ONE - TAU, C, 1)
224230
ELSE
225231
*
226232
* w(1:lastc,1) := C(1:lastc,1:lastv) * v(1:lastv,1)
227233
*
228234
! w(1:lastc,1) := C(1:lastc,1:lastv-1) * v(1:lastv-1,1)
229-
CALL DGEMV( 'No transpose', LASTC, LASTV-1, ONE,
230-
$ C(1,1), LDC, V(1), INCV, ZERO, WORK, 1 )
235+
CALL DGEMV( 'No transpose', LASTC, LASTV-FIRSTV,
236+
$ ONE, C(1,FIRSTV), LDC, V(I), INCV, ZERO, WORK, 1 )
231237
! w(1:lastc,1) += C(1:lastc,lastv) * v(lastv,1) = C(1:lastc,lastv)
232238
CALL DAXPY(LASTC, ONE, C(1,LASTV), 1, WORK, 1)
233239
*
@@ -237,8 +243,8 @@ SUBROUTINE DLARF1L( SIDE, M, N, V, INCV, TAU, C, LDC, WORK )
237243
! = C(...) - tau * w(1:lastc,1)
238244
CALL DAXPY(LASTC, -TAU, WORK, 1, C(1,LASTV), 1)
239245
! C(1:lastc,1:lastv-1) := C(...) - tau * w(1:lastc,1) * v(1:lastv-1)**T
240-
CALL DGER( LASTC, LASTV-1, -TAU, WORK, 1, V(1),
241-
$ INCV, C(1,1), LDC )
246+
CALL DGER( LASTC, LASTV-FIRSTV, -TAU, WORK, 1, V(I),
247+
$ INCV, C(1,FIRSTV), LDC )
242248
END IF
243249
END IF
244250
END IF

0 commit comments

Comments
 (0)