2 use num_types,
only: rp, xp
3 use comm,
only: neko_comm, mpi_extra_precision
4 use mpi_f08,
only: mpi_allreduce, mpi_sum, mpi_in_place
12 subroutine copy_mask(a, b, size, mask, mask_size)
13 integer,
intent(in) :: size, mask_size
14 real(kind=rp),
dimension(size),
intent(inout) :: a
15 real(kind=rp),
dimension(size),
intent(in) :: b
16 integer,
dimension(mask_size),
intent(in) :: mask
20 a(mask(i)) = b(mask(i))
23 end subroutine copy_mask
27 subroutine cadd_mask(a, c, size, mask, mask_size)
28 integer,
intent(in) :: size, mask_size
29 real(kind=rp),
dimension(size),
intent(inout) :: a
30 real(kind=rp),
intent(in) :: c
31 integer,
dimension(mask_size),
intent(in) :: mask
35 a(mask(i)) = a(mask(i)) + c
38 end subroutine cadd_mask
42 subroutine invcol1_mask(a, size, mask, mask_size)
43 integer,
intent(in) :: size, mask_size
44 real(kind=rp),
dimension(size),
intent(inout) :: a
45 integer,
dimension(mask_size),
intent(in) :: mask
49 a(mask(i)) = 1.0_rp / a(mask(i))
52 end subroutine invcol1_mask
56 subroutine cmult_mask(a, c, size, mask, mask_size)
57 integer,
intent(in) :: size, mask_size
58 real(kind=rp),
dimension(size),
intent(inout) :: a
59 real(kind=rp),
intent(in) :: c
60 integer,
dimension(mask_size),
intent(in) :: mask
64 a(mask(i)) = c * a(mask(i))
67 end subroutine cmult_mask
71 subroutine col2_mask(a, b, size, mask, mask_size)
72 integer,
intent(in) :: size, mask_size
73 real(kind=rp),
dimension(size),
intent(inout) :: a
74 real(kind=rp),
dimension(size),
intent(in) :: b
75 integer,
dimension(mask_size),
intent(in) :: mask
79 a(mask(i)) = a(mask(i)) * b(mask(i))
82 end subroutine col2_mask
86 subroutine col3_mask(a, b, c, size, mask, mask_size)
87 integer,
intent(in) :: size, mask_size
88 real(kind=rp),
dimension(size),
intent(inout) :: a
89 real(kind=rp),
dimension(size),
intent(in) :: b, c
90 integer,
dimension(mask_size),
intent(in) :: mask
94 a(mask(i)) = b(mask(i)) * c(mask(i))
97 end subroutine col3_mask
101 subroutine sub3_mask(a, b, c, size, mask, mask_size)
102 integer,
intent(in) :: size, mask_size
103 real(kind=rp),
dimension(size),
intent(inout) :: a
104 real(kind=rp),
dimension(size),
intent(in) :: b, c
105 integer,
dimension(mask_size),
intent(in) :: mask
109 a(mask(i)) = b(mask(i)) - c(mask(i))
112 end subroutine sub3_mask
116 function glsc2_mask(a, b, size, mask, mask_size)
117 integer,
intent(in) :: size, mask_size
118 real(kind=rp),
dimension(size),
intent(in) :: a
119 real(kind=rp),
dimension(size),
intent(in) :: b
120 integer,
dimension(mask_size),
intent(in) :: mask
121 real(kind=rp) :: glsc2_mask
127 tmp = tmp + a(mask(i)) * b(mask(i))
130 call mpi_allreduce(mpi_in_place, tmp, 1, &
131 mpi_extra_precision, mpi_sum, neko_comm, ierr)
133 end function glsc2_mask