xref: /petsc/src/vec/vec/impls/seq/ftn-kernels/fmaxpy.F90 (revision 34c645fd3b0199e05bec2fcc32d3597bfeb7f4f2)
1!
2!
3!    Fortran kernel for the MAXPY() vector routine
4!
5#include <petsc/finclude/petscsys.h>
6!
7
8      Subroutine FortranMAXPY4(x, a0, a1, a2, a3, y0, y1, y2, y3, n)
9      implicit none
10      PetscScalar  a0,a1,a2,a3
11      PetscScalar  x(*),y0(*)
12      PetscScalar y1(*),y2(*),y3(*)
13      PetscInt n, i
14
15      PETSC_AssertAlignx(16,x(1))
16      PETSC_AssertAlignx(16,y0(1))
17      PETSC_AssertAlignx(16,y1(1))
18      PETSC_AssertAlignx(16,Y2(1))
19      PETSC_AssertAlignx(16,Y3(1))
20      do i=1,n
21          x(i)  = x(i) + (a0*y0(i) + a1*y1(i) + a2*y2(i) + a3*y3(i))
22      enddo
23      end
24
25      subroutine FortranMAXPY3(x,a0,a1,a2,y0,y1,y2,n)
26      implicit none
27      PetscScalar  a0,a1,a2,x(*)
28      PetscScalar y0(*),y1(*),y2(*)
29      PetscInt n
30      PetscInt i
31      PETSC_AssertAlignx(16,x(1))
32      PETSC_AssertAlignx(16,y0(1))
33      PETSC_AssertAlignx(16,y1(1))
34      PETSC_AssertAlignx(16,y2(1))
35      do 10,i=1,n
36         x(i) = x(i) + (a0*y0(i) + a1*y1(i) + a2*y2(i))
37  10  continue
38      end
39
40      Subroutine FortranMAXPY2(x, a0, a1, y0, y1, n)
41      implicit none
42      PetscScalar  a0,a1,x(*)
43      PetscScalar  y0(*),y1(*)
44      PetscInt n, i
45      PETSC_AssertAlignx(16,x(1))
46      PETSC_AssertAlignx(16,y0(1))
47      PETSC_AssertAlignx(16,y1(1))
48      do i=1,n
49          x(i)  = x(i) + (a0*y0(i) + a1*y1(i))
50      enddo
51      end
52