xref: /petsc/src/mat/ftn-kernels/sgemv.F90 (revision 0113e7197f2fe148981cb35224c0458dc22f510d)
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>
6*0113e719SMartin Diehl
7c96caaccSSatish Balaysubroutine MSGemv(bs,ncols,A,x,y)
8c96caaccSSatish Balay  implicit none
9c96caaccSSatish Balay  PetscInt          bs,ncols
10c96caaccSSatish Balay  MatScalar        A(bs,ncols)
11c96caaccSSatish Balay  PetscScalar      x(ncols),y(bs)
12c96caaccSSatish Balay
13c96caaccSSatish Balay  PetscInt         i,j
14c96caaccSSatish Balay
15*0113e719SMartin Diehl  do j=1,bs
16c96caaccSSatish Balay    y(j) = 0.0d0
17*0113e719SMartin Diehl  end do
18c96caaccSSatish Balay
19*0113e719SMartin Diehl  do i=1,ncols
20*0113e719SMartin Diehl    do j=1,bs
21c96caaccSSatish Balay      y(j) = y(j) + A(j,i)*x(i)
22*0113e719SMartin Diehl    end do
23*0113e719SMartin Diehl  end do
24c96caaccSSatish Balay
25*0113e719SMartin Diehlend subroutine MSGemv
26c96caaccSSatish Balay
27c96caaccSSatish Balaysubroutine MSGemvp(bs,ncols,A,x,y)
28c96caaccSSatish Balay  implicit none
29c96caaccSSatish Balay  PetscInt          bs,ncols
30c96caaccSSatish Balay  MatScalar        A(bs,ncols)
31c96caaccSSatish Balay  PetscScalar      x(ncols),y(bs)
32c96caaccSSatish Balay
33c96caaccSSatish Balay  PetscInt         i, j
34c96caaccSSatish Balay
35*0113e719SMartin Diehl  do i=1,ncols
36*0113e719SMartin Diehl    do j=1,bs
37c96caaccSSatish Balay      y(j) = y(j) + A(j,i)*x(i)
38*0113e719SMartin Diehl    end do
39*0113e719SMartin Diehl  end do
40c96caaccSSatish Balay
41*0113e719SMartin Diehlend subroutine MSGemvp
42c96caaccSSatish Balay
43c96caaccSSatish Balaysubroutine MSGemvm(bs,ncols,A,x,y)
44c96caaccSSatish Balay  implicit none
45c96caaccSSatish Balay  PetscInt          bs,ncols
46c96caaccSSatish Balay  MatScalar        A(bs,ncols)
47c96caaccSSatish Balay  PetscScalar      x(ncols),y(bs)
48c96caaccSSatish Balay
49c96caaccSSatish Balay  PetscInt         i, j
50c96caaccSSatish Balay
51*0113e719SMartin Diehl  do i=1,ncols
52*0113e719SMartin Diehl    do j=1,bs
53c96caaccSSatish Balay      y(j) = y(j) - A(j,i)*x(i)
54*0113e719SMartin Diehl    end do
55*0113e719SMartin Diehl  end do
56c96caaccSSatish Balay
57*0113e719SMartin Diehlend subroutine MSGemvm
58c96caaccSSatish Balay
59c96caaccSSatish Balaysubroutine MSGemvt(bs,ncols,A,x,y)
60c96caaccSSatish Balay  implicit none
61c96caaccSSatish Balay  PetscInt          bs,ncols
62c96caaccSSatish Balay  MatScalar        A(bs,ncols)
63c96caaccSSatish Balay  PetscScalar      x(bs),y(ncols)
64c96caaccSSatish Balay
65c96caaccSSatish Balay  PetscInt          i,j
66c96caaccSSatish Balay  PetscScalar      sum
67*0113e719SMartin Diehl  do  i=1,ncols
68c96caaccSSatish Balay    sum = y(i)
69*0113e719SMartin Diehl    do  j=1,bs
70c96caaccSSatish Balay      sum = sum + A(j,i)*x(j)
71*0113e719SMartin Diehl    end do
72c96caaccSSatish Balay    y(i) = sum
73*0113e719SMartin Diehl  end do
74c96caaccSSatish Balay
75*0113e719SMartin Diehlend subroutine MSGemvt
76c96caaccSSatish Balay
77c96caaccSSatish Balaysubroutine MSGemm(bs,A,B,C)
78c96caaccSSatish Balay  implicit none
79c96caaccSSatish Balay  PetscInt    bs
80c96caaccSSatish Balay  MatScalar   A(bs,bs),B(bs,bs),C(bs,bs)
81c96caaccSSatish Balay  PetscScalar sum
82c96caaccSSatish Balay  PetscInt    i,j,k
83c96caaccSSatish Balay
84*0113e719SMartin Diehl  do i=1,bs
85*0113e719SMartin Diehl    do j=1,bs
86c96caaccSSatish Balay      sum = A(i,j)
87*0113e719SMartin Diehl      do k=1,bs
88c96caaccSSatish Balay        sum = sum - B(i,k)*C(k,j)
89*0113e719SMartin Diehl      end do
90c96caaccSSatish Balay      A(i,j) = sum
91*0113e719SMartin Diehl    end do
92*0113e719SMartin Diehl  end do
93c96caaccSSatish Balay
94*0113e719SMartin Diehlend subroutine MSGemm
95c96caaccSSatish Balay
96c96caaccSSatish Balaysubroutine MSGemmi(bs,A,C,B)
97c96caaccSSatish Balay  implicit none
98c96caaccSSatish Balay  PetscInt    bs
99c96caaccSSatish Balay  MatScalar   A(bs,bs),B(bs,bs),C(bs,bs)
100c96caaccSSatish Balay  PetscScalar sum
101c96caaccSSatish Balay
102c96caaccSSatish Balay  PetscInt    i,j,k
103c96caaccSSatish Balay
104*0113e719SMartin Diehl  do i=1,bs
105*0113e719SMartin Diehl    do j=1,bs
106c96caaccSSatish Balay      sum = 0.0d0
107*0113e719SMartin Diehl      do  k=1,bs
108c96caaccSSatish Balay        sum = sum + B(i,k)*C(k,j)
109*0113e719SMartin Diehl      end do
110c96caaccSSatish Balay      A(i,j) = sum
111*0113e719SMartin Diehl    end do
112*0113e719SMartin Diehl  end do
113c96caaccSSatish Balay
114*0113e719SMartin Diehlend subroutine MSGemmi
115