Prepared t_mdatoms for using vector
[gromacs.git] / src / gromacs / mdlib / rbin.cpp
blobf180b6c865c850d8718100bf738d170981ac37f0
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2010,2014,2015, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
40 #include "rbin.h"
42 #include "gromacs/gmxlib/network.h"
43 #include "gromacs/utility/smalloc.h"
45 t_bin *mk_bin(void)
47 t_bin *b;
49 snew(b, 1);
51 return b;
54 void destroy_bin(t_bin *b)
56 if (b->maxreal > 0)
58 sfree(b->rbuf);
61 sfree(b);
64 void reset_bin(t_bin *b)
66 b->nreal = 0;
69 int add_binr(t_bin *b, int nr, real r[])
71 #define MULT 4
72 int i, rest, index;
73 double *rbuf;
75 if (b->nreal+nr > b->maxreal)
77 b->maxreal = b->nreal+nr;
78 rest = b->maxreal % MULT;
79 if (rest != 0)
81 b->maxreal += MULT-rest;
83 srenew(b->rbuf, b->maxreal);
85 /* Copy pointer */
86 rbuf = b->rbuf+b->nreal;
88 #if (defined __ICC && __ICC >= 1500 || defined __ICL && __ICL >= 1500) && defined __MIC__
89 #pragma novector /* Work-around for incorrect vectorization */
90 #endif
91 for (i = 0; (i < nr); i++)
93 rbuf[i] = r[i];
96 index = b->nreal;
97 b->nreal += nr;
99 return index;
102 int add_bind(t_bin *b, int nr, double r[])
104 #define MULT 4
105 int i, rest, index;
106 double *rbuf;
108 if (b->nreal+nr > b->maxreal)
110 b->maxreal = b->nreal+nr;
111 rest = b->maxreal % MULT;
112 if (rest != 0)
114 b->maxreal += MULT-rest;
116 srenew(b->rbuf, b->maxreal);
118 /* Copy pointer */
119 rbuf = b->rbuf+b->nreal;
120 for (i = 0; (i < nr); i++)
122 rbuf[i] = r[i];
125 index = b->nreal;
126 b->nreal += nr;
128 return index;
131 void sum_bin(t_bin *b, t_commrec *cr)
133 int i;
135 for (i = b->nreal; (i < b->maxreal); i++)
137 b->rbuf[i] = 0;
139 gmx_sumd(b->maxreal, b->rbuf, cr);
142 void extract_binr(t_bin *b, int index, int nr, real r[])
144 int i;
145 double *rbuf;
147 rbuf = b->rbuf+index;
148 for (i = 0; (i < nr); i++)
150 r[i] = rbuf[i];
154 void extract_bind(t_bin *b, int index, int nr, double r[])
156 int i;
157 double *rbuf;
159 rbuf = b->rbuf+index;
160 for (i = 0; (i < nr); i++)
162 r[i] = rbuf[i];
166 #ifdef DEBUGRBIN
167 int main(int argc, char *argv[])
169 t_commrec *cr;
170 t_bin *rb;
171 double *r;
172 rvec *v;
173 int k, i, ni, mi, n, m;
175 cr = init_par(&argc, argv);
176 n = std::strtol(argv[1], NULL, 10);
177 m = std::strtol(argv[2], NULL, 10);
178 fprintf(stdlog, "n=%d\n", n);
179 rb = mk_bin();
180 snew(r, n);
181 snew(v, m);
183 for (k = 0; (k < 3); k++)
185 fprintf(stdlog, "\nk=%d\n", k);
186 reset_bin(rb);
188 for (i = 0; (i < n); i++)
190 r[i] = i+k;
192 for (i = 0; (i < m); i++)
194 v[i][XX] = 4*i+k;
195 v[i][YY] = 4*i+k+1;
196 v[i][ZZ] = 4*i+k+2;
199 ni = add_bind(stdlog, rb, n, r);
200 mi = add_binr(stdlog, rb, DIM*m, v[0]);
202 sum_bin(rb, cr);
204 extract_bind(rb, ni, n, r);
205 extract_binr(rb, mi, DIM*m, v[0]);
207 for (i = 0; (i < n); i++)
209 fprintf(stdlog, "r[%d] = %e\n", i, r[i]);
211 for (i = 0; (i < m); i++)
213 fprintf(stdlog, "v[%d] = (%e,%e,%e)\n", i, v[i][XX], v[i][YY], v[i][ZZ]);
216 fflush(stdlog);
218 return 0;
220 #endif