Octopus
operate.c
Go to the documentation of this file.
1/*
2 Copyright (C) 2006 X. Andrade
3
4 This program is free software; you can redistribute it and/or modify
5 it under the terms of the GNU General Public License as published by
6 the Free Software Foundation; either version 2, or (at your option)
7 any later version.
8
9 This program is distributed in the hope that it will be useful,
10 but WITHOUT ANY WARRANTY; without even the implied warranty of
11 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 GNU General Public License for more details.
13
14 You should have received a copy of the GNU General Public License
15 along with this program; if not, write to the Free Software
16 Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17 02110-1301, USA.
18
19*/
20
21#include <config.h>
23#include <stddef.h>
24
25#include <assert.h>
26#include <vectors.h>
27
40void FC_FUNC_(doperate_ri_vec,
41 DOPERATE_RI_VEC)(const int *opn, const double *restrict w,
42 const int *opnri, const int *opri,
43 const int *rimap_inv, const int *rimap_inv_max,
44 const double *restrict fi, const int *ldfp,
45 double *restrict fo) {
46 const size_t ldf = ldfp[0];
47
48 /* check whether we got aligned vectors or not */
49 int aligned = 1;
50 aligned = aligned && (((long long)fi) % (8 * VEC_SIZE) == 0);
51 aligned = aligned && (((long long)fo) % (8 * VEC_SIZE) == 0);
52 aligned = aligned && ((1 << ldf) % VEC_SIZE == 0);
53
54 if (aligned) {
55#define ALIGNED
56#include "operate_inc.c"
57#undef ALIGNED
58
59 } else {
60 /* not aligned */
61#include "operate_inc.c"
62 }
64
65/* the same as doperate_ri_vec, but allows giving each an appropriate Fortan
66 interface in which fi and fo are actually complex in Fortran Could be inline,
67 but in that case pgcc will not put it in the symbol table. */
68void FC_FUNC_(zoperate_ri_vec,
69 ZOPERATE_RI_VEC)(const int *opn, const double *restrict w,
70 const int *opnri, const int *opri,
71 const int *rimap_inv, const int *rimap_inv_max,
72 const double *restrict fi, const int *ldfp,
73 double *restrict fo) {
74 FC_FUNC_(doperate_ri_vec, DOPERATE_RI_VEC)
75 (opn, w, opnri, opri, rimap_inv, rimap_inv_max, fi, ldfp, fo);
76}
77
78void FC_FUNC_(dgauss_seidel,
79 DGAUSS_SEIDEL)(const int *opn, const double *restrict w,
80 const int *opnri, const int *opri,
81 const int *rimap_inv, const int *rimap_inv_max,
82 const double *restrict factor, double *pot,
83 const double *restrict rho) {
84
85 const int n = opn[0];
86 const int nri = opnri[0];
88 int l, i, j;
89 const int *index;
90 register double a0;
91 register const double fac = *factor;
92
93 for (l = 0; l < nri; l++) {
94 index = opri + n * l;
95 i = rimap_inv[l];
96
97 for (; i < rimap_inv_max[l]; i++) {
98 a0 = w[0] * pot[i + index[0]];
99 for (j = 1; j < n; j++)
100 a0 += w[j] * pot[i + index[j]];
101 pot[i] += fac * (a0 - rho[i]);
102 }
103 }
104}
const ptrdiff_t nri
Definition: operate_inc.c:10
ptrdiff_t l
Definition: operate_inc.c:12
ptrdiff_t i
Definition: operate_inc.c:12
ptrdiff_t j
Definition: operate_inc.c:12
const int *restrict index
Definition: operate_inc.c:13