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