@@ -43,25 +43,25 @@ contains
43
43
if ( norm_sq0 > zero_${s}$ ) then
44
44
45
45
R = B
46
- call A%apply (X, R, alpha= -one_${s}$, beta=one_${s}$) ! R = B - A*X
46
+ call A%matvec (X, R, alpha= -one_${s}$, beta=one_${s}$, op='N' ) ! R = B - A*X
47
47
48
- call M%apply (R,P, alpha= one_${s}$, beta=zero_${s}$) ! P = M^{-1}*R
48
+ call M%matvec (R,P, alpha= one_${s}$, beta=zero_${s}$, op='N' ) ! P = M^{-1}*R
49
49
50
50
tolsq = tol*tol
51
51
52
52
zr1 = zero_${s}$
53
53
zr2 = one_${s}$
54
54
do while ( (iter < maxiter) .AND. (norm_sq > tolsq * norm_sq0) )
55
-
56
- call M%apply (R,S, alpha= one_${s}$, beta=zero_${s}$) ! S = M^{-1}*R
55
+
56
+ call M%matvec (R,S, alpha= one_${s}$, beta=zero_${s}$, op='N' ) ! S = M^{-1}*R
57
57
zr2 = A%inner_product( R, S )
58
58
59
59
if (iter>0) then
60
60
beta = zr2 / zr1
61
61
P = S + beta * P
62
62
end if
63
-
64
- call A%apply (P, Q, alpha= one_${s}$, beta=zero_${s}$) ! Q = A*P
63
+
64
+ call A%matvec (P, Q, alpha= one_${s}$, beta=zero_${s}$, op='N' ) ! Q = A*P
65
65
zv2 = A%inner_product( P, Q )
66
66
67
67
alpha = zr2 / zv2
@@ -134,14 +134,14 @@ contains
134
134
#:else
135
135
call diag(A,diagonal)
136
136
#:endif
137
- M_%apply => precond_jacobi
137
+ M_%matvec => precond_jacobi
138
138
case default
139
- M_%apply => precond_none
139
+ M_%matvec => precond_none
140
140
end select
141
141
where(abs(diagonal)>epsilon(zero_${s}$)) diagonal = one_${s}$/diagonal
142
142
end if
143
143
! matvec for the operator
144
- op%apply => matvec
144
+ op%matvec => matvec
145
145
146
146
! direchlet boundary conditions mask
147
147
if(present(di))then
@@ -176,31 +176,37 @@ contains
176
176
workspace_ => null()
177
177
contains
178
178
179
- subroutine matvec(x,y,alpha,beta)
179
+ subroutine matvec(x,y,alpha,beta,op)
180
+ #:if matrix == "dense"
181
+ use stdlib_linalg_blas, only: gemv
182
+ #:endif
180
183
${t}$, intent(in) :: x(:)
181
184
${t}$, intent(inout) :: y(:)
182
185
${t}$, intent(in) :: alpha
183
186
${t}$, intent(in) :: beta
187
+ character(1), intent(in) :: op
184
188
#:if matrix == "dense"
185
- y = alpha * matmul (A,x) + beta * y
189
+ call gemv(op,m=size(A,1),n=size(A,2), alpha=alpha,a=A,lda=size (A,1),x=x,incx=1, beta=beta,y=y,incy=1)
186
190
#:else
187
- call spmv( A , x, y , alpha, beta )
191
+ call spmv( A , x, y , alpha, beta , op )
188
192
#:endif
189
193
y = merge( zero_${s}$, y, di_ )
190
194
end subroutine
191
195
192
- subroutine precond_none(x,y,alpha,beta)
196
+ subroutine precond_none(x,y,alpha,beta,op )
193
197
${t}$, intent(in) :: x(:)
194
198
${t}$, intent(inout) :: y(:)
195
199
${t}$, intent(in) :: alpha
196
200
${t}$, intent(in) :: beta
201
+ character(1), intent(in) :: op
197
202
y = merge( zero_${s}$, x, di_ )
198
203
end subroutine
199
- subroutine precond_jacobi(x,y,alpha,beta)
204
+ subroutine precond_jacobi(x,y,alpha,beta,op )
200
205
${t}$, intent(in) :: x(:)
201
206
${t}$, intent(inout) :: y(:)
202
207
${t}$, intent(in) :: alpha
203
208
${t}$, intent(in) :: beta
209
+ character(1), intent(in) :: op
204
210
y = merge( zero_${s}$, diagonal * x, di_ ) ! inverted diagonal
205
211
end subroutine
206
212
end subroutine
0 commit comments