xref: /petsc/src/vec/vec/impls/seq/ftn-kernels/fmdot.F90 (revision edb0e59d3c097acd4a4005a4e51d4daa5c739255)
1c96caaccSSatish Balay!
2c96caaccSSatish Balay!
3c96caaccSSatish Balay!    Fortran kernel for the MDot() vector routine
4c96caaccSSatish Balay!
5c96caaccSSatish Balay#include <petsc/finclude/petscsys.h>
6c96caaccSSatish Balay!
70ccf82acSMartin Diehlpure subroutine FortranMDot4(x, y1, y2, y3, y4, n, sum1, sum2, sum3, sum4)
8*fe66ebccSMartin Diehl  use, intrinsic :: ISO_C_binding
90ccf82acSMartin Diehl  implicit none(type, external)
100ccf82acSMartin Diehl  PetscScalar, intent(inout) :: sum1, sum2, sum3, sum4
110ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(*), y1(*), y2(*), y3(*), y4(*)
120ccf82acSMartin Diehl  PetscInt, intent(in) :: n
130113e719SMartin Diehl
14d66e387eSMartin Diehl  PetscInt :: i
15c96caaccSSatish Balay
16c96caaccSSatish Balay  PETSC_AssertAlignx(16, x(1))
17c96caaccSSatish Balay  PETSC_AssertAlignx(16, y1(1))
18c96caaccSSatish Balay  PETSC_AssertAlignx(16, y2(1))
19c96caaccSSatish Balay  PETSC_AssertAlignx(16, y3(1))
20c96caaccSSatish Balay  PETSC_AssertAlignx(16, y4(1))
21c96caaccSSatish Balay
220113e719SMartin Diehl  do i = 1, n
23c96caaccSSatish Balay    sum1 = sum1 + x(i)*PetscConj(y1(i))
24c96caaccSSatish Balay    sum2 = sum2 + x(i)*PetscConj(y2(i))
25c96caaccSSatish Balay    sum3 = sum3 + x(i)*PetscConj(y3(i))
26c96caaccSSatish Balay    sum4 = sum4 + x(i)*PetscConj(y4(i))
270113e719SMartin Diehl  end do
280113e719SMartin Diehlend subroutine FortranMDot4
29c96caaccSSatish Balay
300ccf82acSMartin Diehlpure subroutine FortranMDot3(x, y1, y2, y3, n, sum1, sum2, sum3)
31*fe66ebccSMartin Diehl  use, intrinsic :: ISO_C_binding
320ccf82acSMartin Diehl  implicit none(type, external)
330ccf82acSMartin Diehl  PetscScalar, intent(inout) :: sum1, sum2, sum3
340ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(*), y1(*), y2(*), y3(*)
350ccf82acSMartin Diehl  PetscInt, intent(in) :: n
360113e719SMartin Diehl
37d66e387eSMartin Diehl  PetscInt :: i
38c96caaccSSatish Balay
39c96caaccSSatish Balay  PETSC_AssertAlignx(16, x(1))
40c96caaccSSatish Balay  PETSC_AssertAlignx(16, y1(1))
41c96caaccSSatish Balay  PETSC_AssertAlignx(16, y2(1))
42c96caaccSSatish Balay  PETSC_AssertAlignx(16, y3(1))
430113e719SMartin Diehl
440113e719SMartin Diehl  do i = 1, n
45c96caaccSSatish Balay    sum1 = sum1 + x(i)*PetscConj(y1(i))
46c96caaccSSatish Balay    sum2 = sum2 + x(i)*PetscConj(y2(i))
47c96caaccSSatish Balay    sum3 = sum3 + x(i)*PetscConj(y3(i))
480113e719SMartin Diehl  end do
490113e719SMartin Diehlend subroutine FortranMDot3
50c96caaccSSatish Balay
510ccf82acSMartin Diehlpure subroutine FortranMDot2(x, y1, y2, n, sum1, sum2)
52*fe66ebccSMartin Diehl  use, intrinsic :: ISO_C_binding
530ccf82acSMartin Diehl  implicit none(type, external)
540ccf82acSMartin Diehl  PetscScalar, intent(inout) :: sum1, sum2
550ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(*), y1(*), y2(*)
560ccf82acSMartin Diehl  PetscInt, intent(in) :: n
570113e719SMartin Diehl
58d66e387eSMartin Diehl  PetscInt :: i
59c96caaccSSatish Balay
60c96caaccSSatish Balay  PETSC_AssertAlignx(16, x(1))
61c96caaccSSatish Balay  PETSC_AssertAlignx(16, y1(1))
62c96caaccSSatish Balay  PETSC_AssertAlignx(16, y2(1))
630113e719SMartin Diehl
640113e719SMartin Diehl  do i = 1, n
65c96caaccSSatish Balay    sum1 = sum1 + x(i)*PetscConj(y1(i))
66c96caaccSSatish Balay    sum2 = sum2 + x(i)*PetscConj(y2(i))
670113e719SMartin Diehl  end do
680113e719SMartin Diehlend subroutine FortranMDot2
69c96caaccSSatish Balay
700ccf82acSMartin Diehlpure subroutine FortranMDot1(x, y1, n, sum1)
71*fe66ebccSMartin Diehl  use, intrinsic :: ISO_C_binding
720ccf82acSMartin Diehl  implicit none(type, external)
730ccf82acSMartin Diehl  PetscScalar, intent(inout) :: sum1
740ccf82acSMartin Diehl  PetscScalar, intent(in) :: x(*), y1(*)
750ccf82acSMartin Diehl  PetscInt, intent(in) :: n
760113e719SMartin Diehl
77d66e387eSMartin Diehl  PetscInt :: i
78c96caaccSSatish Balay
79c96caaccSSatish Balay  PETSC_AssertAlignx(16, x(1))
80c96caaccSSatish Balay  PETSC_AssertAlignx(16, y1(1))
81c96caaccSSatish Balay
820113e719SMartin Diehl  do i = 1, n
830113e719SMartin Diehl    sum1 = sum1 + x(i)*PetscConj(y1(i))
840113e719SMartin Diehl  end do
850113e719SMartin Diehl
860113e719SMartin Diehlend subroutine FortranMDot1
87