Skip to content

Commit beafb3c

Browse files
committed
simplify csr spmv, remove unused var beta_
1 parent 59d33f0 commit beafb3c

File tree

1 file changed

+16
-52
lines changed

1 file changed

+16
-52
lines changed

src/stdlib_sparse_spmv.fypp

Lines changed: 16 additions & 52 deletions
Original file line numberDiff line numberDiff line change
@@ -45,7 +45,7 @@ contains
4545
${t1}$, intent(in), optional :: alpha
4646
${t1}$, intent(in), optional :: beta
4747
character(1), intent(in), optional :: op
48-
${t1}$ :: alpha_, beta_
48+
${t1}$ :: alpha_
4949
character(1) :: op_
5050
integer(ilp) :: k, ik, jk
5151

@@ -57,6 +57,7 @@ contains
5757
else
5858
vec_y = zero_${s1}$
5959
endif
60+
6061
associate( data => matrix%data, index => matrix%index, storage => matrix%storage, nnz => matrix%nnz )
6162
select case(op_)
6263
case(sparse_op_none)
@@ -132,7 +133,7 @@ contains
132133
${t1}$, intent(in), optional :: alpha
133134
${t1}$, intent(in), optional :: beta
134135
character(1), intent(in), optional :: op
135-
${t1}$ :: alpha_, beta_
136+
${t1}$ :: alpha_
136137
character(1) :: op_
137138
integer(ilp) :: i, j
138139
#:if rank == 1
@@ -144,8 +145,11 @@ contains
144145
op_ = sparse_op_none; if(present(op)) op_ = op
145146
alpha_ = one_${k1}$
146147
if(present(alpha)) alpha_ = alpha
147-
beta_ = zero_${k1}$
148-
if(present(beta)) beta_ = beta
148+
if(present(beta)) then
149+
vec_y = beta * vec_y
150+
else
151+
vec_y = zero_${s1}$
152+
endif
149153

150154
associate( data => matrix%data, col => matrix%col, rowptr => matrix%rowptr, &
151155
& nnz => matrix%nnz, nrows => matrix%nrows, ncols => matrix%ncols, storage => matrix%storage )
@@ -156,19 +160,10 @@ contains
156160
do j = rowptr(i), rowptr(i+1)-1
157161
aux = aux + data(j) * vec_x(${rksfx2(rank-1)}$col(j))
158162
end do
159-
if(present(beta)) then
160-
vec_y(${rksfx2(rank-1)}$i) = beta_ * vec_y(${rksfx2(rank-1)}$i) + alpha_ * aux
161-
else
162-
vec_y(${rksfx2(rank-1)}$i) = alpha_ * aux
163-
end if
163+
vec_y(${rksfx2(rank-1)}$i) = vec_y(${rksfx2(rank-1)}$i) + alpha_ * aux
164164
end do
165165

166166
else if( storage == sparse_full .and. op_==sparse_op_transpose ) then
167-
if(present(beta)) then
168-
vec_y = beta * vec_y
169-
else
170-
vec_y = zero_${s1}$
171-
endif
172167
do i = 1, nrows
173168
aux = alpha_ * vec_x(${rksfx2(rank-1)}$i)
174169
do j = rowptr(i), rowptr(i+1)-1
@@ -185,12 +180,7 @@ contains
185180
vec_y(${rksfx2(rank-1)}$col(j)) = vec_y(${rksfx2(rank-1)}$col(j)) + data(j) * aux2
186181
end do
187182
aux = alpha_ * aux + data(j) * aux2
188-
189-
if(present(beta)) then
190-
vec_y(${rksfx2(rank-1)}$i) = beta_ * vec_y(${rksfx2(rank-1)}$i) + aux
191-
else
192-
vec_y(${rksfx2(rank-1)}$i) = aux
193-
end if
183+
vec_y(${rksfx2(rank-1)}$i) = vec_y(${rksfx2(rank-1)}$i) + aux
194184
end do
195185

196186
else if( storage == sparse_upper .and. op_/=sparse_op_hermitian )then
@@ -199,26 +189,13 @@ contains
199189
aux2 = alpha_ * vec_x(${rksfx2(rank-1)}$i)
200190
do j = rowptr(i)+1, rowptr(i+1)-1
201191
aux = aux + data(j) * vec_x(${rksfx2(rank-1)}$col(j))
192+
vec_y(${rksfx2(rank-1)}$col(j)) = vec_y(${rksfx2(rank-1)}$col(j)) + data(j) * aux2
202193
end do
203-
if(present(beta)) then
204-
do j = rowptr(i)+1, rowptr(i+1)-1
205-
vec_y(${rksfx2(rank-1)}$col(j)) = beta_ * vec_y(${rksfx2(rank-1)}$col(j)) + data(j) * aux2
206-
end do
207-
else
208-
do j = rowptr(i)+1, rowptr(i+1)-1
209-
vec_y(${rksfx2(rank-1)}$col(j)) = data(j) * aux2
210-
end do
211-
end if
212194
vec_y(${rksfx2(rank-1)}$i) = vec_y(${rksfx2(rank-1)}$i) + alpha_ * aux
213195
end do
214196

215197
#:if t1.startswith('complex')
216198
else if( storage == sparse_full .and. op_==sparse_op_hermitian) then
217-
if(present(beta)) then
218-
vec_y = beta * vec_y
219-
else
220-
vec_y = zero_${s1}$
221-
endif
222199
do i = 1, nrows
223200
aux = alpha_ * vec_x(${rksfx2(rank-1)}$i)
224201
do j = rowptr(i), rowptr(i+1)-1
@@ -235,12 +212,7 @@ contains
235212
vec_y(${rksfx2(rank-1)}$col(j)) = vec_y(${rksfx2(rank-1)}$col(j)) + conjg(data(j)) * aux2
236213
end do
237214
aux = alpha_ * aux + conjg(data(j)) * aux2
238-
239-
if(present(beta)) then
240-
vec_y(${rksfx2(rank-1)}$i) = beta_ * vec_y(${rksfx2(rank-1)}$i) + aux
241-
else
242-
vec_y(${rksfx2(rank-1)}$i) = aux
243-
end if
215+
vec_y(${rksfx2(rank-1)}$i) = vec_y(${rksfx2(rank-1)}$i) + aux
244216
end do
245217

246218
else if( storage == sparse_upper .and. op_==sparse_op_hermitian )then
@@ -249,16 +221,8 @@ contains
249221
aux2 = alpha_ * vec_x(${rksfx2(rank-1)}$i)
250222
do j = rowptr(i)+1, rowptr(i+1)-1
251223
aux = aux + conjg(data(j)) * vec_x(${rksfx2(rank-1)}$col(j))
224+
vec_y(${rksfx2(rank-1)}$col(j)) = vec_y(${rksfx2(rank-1)}$col(j)) + conjg(data(j)) * aux2
252225
end do
253-
if(present(beta)) then
254-
do j = rowptr(i)+1, rowptr(i+1)-1
255-
vec_y(${rksfx2(rank-1)}$col(j)) = beta_ * vec_y(${rksfx2(rank-1)}$col(j)) + conjg(data(j)) * aux2
256-
end do
257-
else
258-
do j = rowptr(i)+1, rowptr(i+1)-1
259-
vec_y(${rksfx2(rank-1)}$col(j)) = conjg(data(j)) * aux2
260-
end do
261-
end if
262226
vec_y(${rksfx2(rank-1)}$i) = vec_y(${rksfx2(rank-1)}$i) + alpha_ * aux
263227
end do
264228
#:endif
@@ -279,7 +243,7 @@ contains
279243
${t1}$, intent(in), optional :: alpha
280244
${t1}$, intent(in), optional :: beta
281245
character(1), intent(in), optional :: op
282-
${t1}$ :: alpha_, beta_
246+
${t1}$ :: alpha_
283247
character(1) :: op_
284248
integer(ilp) :: i, j
285249
#:if rank == 1
@@ -385,7 +349,7 @@ contains
385349
${t1}$, intent(in), optional :: alpha
386350
${t1}$, intent(in), optional :: beta
387351
character(1), intent(in), optional :: op
388-
${t1}$ :: alpha_, beta_
352+
${t1}$ :: alpha_
389353
character(1) :: op_
390354
integer(ilp) :: i, j, k
391355

@@ -452,7 +416,7 @@ contains
452416
${t1}$, intent(in), optional :: alpha
453417
${t1}$, intent(in), optional :: beta
454418
character(1), intent(in), optional :: op
455-
${t1}$ :: alpha_, beta_
419+
${t1}$ :: alpha_
456420
character(1) :: op_
457421
integer(ilp) :: i, nz, rowidx, num_chunks, rm
458422

0 commit comments

Comments
 (0)