Actual source code: fmdot.F90
1: !
2: !
3: ! Fortran kernel for the MDot() vector routine
4: !
5: #include <petsc/finclude/petscsys.h>
6: !
7: pure subroutine FortranMDot4(x,y1,y2,y3,y4,n,sum1,sum2,sum3,sum4)
8: use, intrinsic :: ISO_C_binding
9: implicit none (type, external)
10: PetscScalar, intent(inout) :: sum1,sum2,sum3,sum4
11: PetscScalar, intent(in) :: x(*),y1(*),y2(*),y3(*),y4(*)
12: PetscInt, intent(in) :: n
14: PetscInt :: i
16: PETSC_AssertAlignx(16,x(1))
17: PETSC_AssertAlignx(16,y1(1))
18: PETSC_AssertAlignx(16,y2(1))
19: PETSC_AssertAlignx(16,y3(1))
20: PETSC_AssertAlignx(16,y4(1))
22: do i=1,n
23: sum1 = sum1 + x(i)*PetscConj(y1(i))
24: sum2 = sum2 + x(i)*PetscConj(y2(i))
25: sum3 = sum3 + x(i)*PetscConj(y3(i))
26: sum4 = sum4 + x(i)*PetscConj(y4(i))
27: end do
28: end subroutine FortranMDot4
30: pure subroutine FortranMDot3(x,y1,y2,y3,n,sum1,sum2,sum3)
31: use, intrinsic :: ISO_C_binding
32: implicit none (type, external)
33: PetscScalar, intent(inout) :: sum1,sum2,sum3
34: PetscScalar, intent(in) :: x(*),y1(*),y2(*),y3(*)
35: PetscInt, intent(in) :: n
37: PetscInt :: i
39: PETSC_AssertAlignx(16,x(1))
40: PETSC_AssertAlignx(16,y1(1))
41: PETSC_AssertAlignx(16,y2(1))
42: PETSC_AssertAlignx(16,y3(1))
44: do i=1,n
45: sum1 = sum1 + x(i)*PetscConj(y1(i))
46: sum2 = sum2 + x(i)*PetscConj(y2(i))
47: sum3 = sum3 + x(i)*PetscConj(y3(i))
48: end do
49: end subroutine FortranMDot3
51: pure subroutine FortranMDot2(x,y1,y2,n,sum1,sum2)
52: use, intrinsic :: ISO_C_binding
53: implicit none (type, external)
54: PetscScalar, intent(inout) :: sum1,sum2
55: PetscScalar, intent(in) :: x(*),y1(*),y2(*)
56: PetscInt, intent(in) :: n
58: PetscInt :: i
60: PETSC_AssertAlignx(16,x(1))
61: PETSC_AssertAlignx(16,y1(1))
62: PETSC_AssertAlignx(16,y2(1))
64: do i=1,n
65: sum1 = sum1 + x(i)*PetscConj(y1(i))
66: sum2 = sum2 + x(i)*PetscConj(y2(i))
67: end do
68: end subroutine FortranMDot2
70: pure subroutine FortranMDot1(x,y1,n,sum1)
71: use, intrinsic :: ISO_C_binding
72: implicit none (type, external)
73: PetscScalar, intent(inout) :: sum1
74: PetscScalar, intent(in) :: x(*),y1(*)
75: PetscInt, intent(in) :: n
77: PetscInt :: i
79: PETSC_AssertAlignx(16,x(1))
80: PETSC_AssertAlignx(16,y1(1))
82: do i=1,n
83: sum1 = sum1 + x(i)*PetscConj(y1(i))
84: end do
86: end subroutine FortranMDot1