xref: /petsc/src/mat/ftn-kernels/sgemv.F90 (revision ff45ff59e5033a47c65b211fde43097dd9cb5a3c)
1c96caaccSSatish Balay!
2c96caaccSSatish Balay!    Fortran kernel for gemv() BLAS operation. This version supports
3c96caaccSSatish Balay!  matrix array stored in single precision but vectors in double
4c96caaccSSatish Balay!
5c96caaccSSatish Balay#include <petsc/finclude/petscsys.h>
60113e719SMartin Diehl
70ccf82acSMartin Diehlpure subroutine MSGemv(bs,ncols,A,x,y)
80ccf82acSMartin Diehl  implicit none (type, external)
90ccf82acSMartin Diehl  PetscInt, intent(in) :: bs,ncols
100ccf82acSMartin Diehl  MatScalar, intent(in) :: A(bs,ncols)
110ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(ncols)
120ccf82acSMartin Diehl  PetscScalar, intent(out) :: y(bs)
13c96caaccSSatish Balay
14d66e387eSMartin Diehl  PetscInt :: i
15c96caaccSSatish Balay
16*ff45ff59SMartin Diehl  y(1:bs) = 0.0
170113e719SMartin Diehl  do i=1,ncols
18d66e387eSMartin Diehl    y(1:bs) = y(1:bs) + A(1:bs,i)*x(i)
190113e719SMartin Diehl  end do
200113e719SMartin Diehlend subroutine MSGemv
21c96caaccSSatish Balay
220ccf82acSMartin Diehlpure subroutine MSGemvp(bs,ncols,A,x,y)
230ccf82acSMartin Diehl  implicit none (type, external)
240ccf82acSMartin Diehl  PetscInt, intent(in) :: bs,ncols
250ccf82acSMartin Diehl  MatScalar, intent(in) :: A(bs,ncols)
260ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(ncols)
270ccf82acSMartin Diehl  PetscScalar, intent(inout) :: y(bs)
28c96caaccSSatish Balay
29d66e387eSMartin Diehl  PetscInt :: i
30c96caaccSSatish Balay
310113e719SMartin Diehl  do i=1,ncols
32d66e387eSMartin Diehl    y(1:bs) = y(1:bs) + A(1:bs,i)*x(i)
330113e719SMartin Diehl  end do
340113e719SMartin Diehlend subroutine MSGemvp
35c96caaccSSatish Balay
360ccf82acSMartin Diehlpure subroutine MSGemvm(bs,ncols,A,x,y)
370ccf82acSMartin Diehl  implicit none (type, external)
380ccf82acSMartin Diehl  PetscInt, intent(in) :: bs,ncols
390ccf82acSMartin Diehl  MatScalar, intent(in) :: A(bs,ncols)
400ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(ncols)
410ccf82acSMartin Diehl  PetscScalar, intent(inout) :: y(bs)
42c96caaccSSatish Balay
43d66e387eSMartin Diehl  PetscInt :: i
44c96caaccSSatish Balay
450113e719SMartin Diehl  do i=1,ncols
46d66e387eSMartin Diehl    y(1:bs) = y(1:bs) - A(1:bs,i)*x(i)
470113e719SMartin Diehl  end do
480113e719SMartin Diehlend subroutine MSGemvm
49c96caaccSSatish Balay
500ccf82acSMartin Diehlpure subroutine MSGemvt(bs,ncols,A,x,y)
510ccf82acSMartin Diehl  implicit none (type, external)
520ccf82acSMartin Diehl  PetscInt, intent(in) :: bs,ncols
530ccf82acSMartin Diehl  MatScalar, intent(in) :: A(bs,ncols)
540ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(bs)
550ccf82acSMartin Diehl  PetscScalar, intent(inout) :: y(ncols)
56c96caaccSSatish Balay
57d66e387eSMartin Diehl  PetscInt :: i
580ccf82acSMartin Diehl
590113e719SMartin Diehl  do  i=1,ncols
60d66e387eSMartin Diehl    y(i) = y(i) + sum(A(1:bs,i)*x(1:bs))
610113e719SMartin Diehl  end do
620113e719SMartin Diehlend subroutine MSGemvt
63c96caaccSSatish Balay
640ccf82acSMartin Diehlpure subroutine MSGemm(bs,A,B,C)
650ccf82acSMartin Diehl  implicit none (type, external)
660ccf82acSMartin Diehl  PetscInt, intent(in) :: bs
670ccf82acSMartin Diehl  MatScalar, intent(in) :: B(bs,bs),C(bs,bs)
680ccf82acSMartin Diehl  MatScalar, intent(inout) :: A(bs,bs)
690ccf82acSMartin Diehl
70d66e387eSMartin Diehl  PetscInt :: i,j
71c96caaccSSatish Balay
720113e719SMartin Diehl  do i=1,bs
730113e719SMartin Diehl    do j=1,bs
74d66e387eSMartin Diehl      A(i,j) = A(i,j) - sum(B(i,1:bs)*C(1:bs,j))
750113e719SMartin Diehl    end do
760113e719SMartin Diehl  end do
770113e719SMartin Diehlend subroutine MSGemm
78c96caaccSSatish Balay
790ccf82acSMartin Diehlpure subroutine MSGemmi(bs,A,C,B)
800ccf82acSMartin Diehl  implicit none (type, external)
810ccf82acSMartin Diehl  PetscInt, intent(in) :: bs
820ccf82acSMartin Diehl  MatScalar, intent(in) :: B(bs,bs),C(bs,bs)
830ccf82acSMartin Diehl  MatScalar, intent(out) :: A(bs,bs)
84c96caaccSSatish Balay
85d66e387eSMartin Diehl  PetscInt :: i,j
86c96caaccSSatish Balay
870113e719SMartin Diehl  do i=1,bs
880113e719SMartin Diehl    do j=1,bs
89d66e387eSMartin Diehl      A(i,j) = sum(B(i,1:bs)*C(1:bs,j))
900113e719SMartin Diehl    end do
910113e719SMartin Diehl  end do
920113e719SMartin Diehlend subroutine MSGemmi
93