xref: /petsc/src/mat/ftn-kernels/sgemv.F90 (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
1!
2!    Fortran kernel for gemv() BLAS operation. This version supports
3!  matrix array stored in single precision but vectors in double
4!
5#include <petsc/finclude/petscsys.h>
6!
7      subroutine MSGemv(bs,ncols,A,x,y)
8      implicit none
9      PetscInt          bs,ncols
10      MatScalar        A(bs,ncols)
11      PetscScalar      x(ncols),y(bs)
12
13      PetscInt         i,j
14
15      do 5, j=1,bs
16        y(j) = 0.0d0
17 5    continue
18
19      do 10, i=1,ncols
20        do 20, j=1,bs
21          y(j) = y(j) + A(j,i)*x(i)
22 20     continue
23 10   continue
24
25      return
26      end
27
28
29      subroutine MSGemvp(bs,ncols,A,x,y)
30      implicit none
31      PetscInt          bs,ncols
32      MatScalar        A(bs,ncols)
33      PetscScalar      x(ncols),y(bs)
34
35      PetscInt         i, j
36
37      do 10, i=1,ncols
38        do 20, j=1,bs
39          y(j) = y(j) + A(j,i)*x(i)
40 20     continue
41 10   continue
42
43      return
44      end
45
46
47      subroutine MSGemvm(bs,ncols,A,x,y)
48      implicit none
49      PetscInt          bs,ncols
50      MatScalar        A(bs,ncols)
51      PetscScalar      x(ncols),y(bs)
52
53      PetscInt         i, j
54
55      do 10, i=1,ncols
56        do 20, j=1,bs
57          y(j) = y(j) - A(j,i)*x(i)
58 20     continue
59 10   continue
60
61      return
62      end
63
64
65      subroutine MSGemvt(bs,ncols,A,x,y)
66      implicit none
67      PetscInt          bs,ncols
68      MatScalar        A(bs,ncols)
69      PetscScalar      x(bs),y(ncols)
70
71      PetscInt          i,j
72      PetscScalar      sum
73      do 10, i=1,ncols
74        sum = y(i)
75        do 20, j=1,bs
76          sum = sum + A(j,i)*x(j)
77 20     continue
78        y(i) = sum
79 10   continue
80
81      return
82      end
83
84      subroutine MSGemm(bs,A,B,C)
85      implicit none
86      PetscInt    bs
87      MatScalar   A(bs,bs),B(bs,bs),C(bs,bs)
88      PetscScalar sum
89      PetscInt    i,j,k
90
91
92      do 10, i=1,bs
93        do 20, j=1,bs
94          sum = A(i,j)
95          do 30, k=1,bs
96            sum = sum - B(i,k)*C(k,j)
97 30       continue
98          A(i,j) = sum
99 20     continue
100 10   continue
101
102      return
103      end
104
105
106      subroutine MSGemmi(bs,A,C,B)
107      implicit none
108      PetscInt    bs
109      MatScalar   A(bs,bs),B(bs,bs),C(bs,bs)
110      PetscScalar sum
111
112      PetscInt    i,j,k
113
114
115      do 10, i=1,bs
116        do 20, j=1,bs
117          sum = 0.0d0
118          do 30, k=1,bs
119            sum = sum + B(i,k)*C(k,j)
120 30       continue
121          A(i,j) = sum
122 20     continue
123 10   continue
124
125      return
126      end
127