added Verlet scheme and NxN non-bonded functionality
[gromacs.git] / src / gmxlib / rbin.c
blob152bad932f4188ec4a8bb8ed25796a62ab8126a2
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * GROningen Mixture of Alchemy and Childrens' Stories
35 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include "typedefs.h"
41 #include "main.h"
42 #include "network.h"
43 #include "rbin.h"
44 #include "smalloc.h"
46 t_bin *mk_bin(void)
48 t_bin *b;
50 snew(b,1);
52 return b;
55 void destroy_bin(t_bin *b)
57 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) {
76 b->maxreal=b->nreal+nr;
77 rest=b->maxreal % MULT;
78 if (rest != 0)
79 b->maxreal+=MULT-rest;
80 srenew(b->rbuf,b->maxreal);
82 /* Copy pointer */
83 rbuf=b->rbuf+b->nreal;
84 for(i=0; (i<nr); i++)
85 rbuf[i]=r[i];
87 index=b->nreal;
88 b->nreal+=nr;
90 return index;
93 int add_bind(t_bin *b,int nr,double r[])
95 #define MULT 4
96 int i,rest,index;
97 double *rbuf;
99 if (b->nreal+nr > b->maxreal) {
100 b->maxreal=b->nreal+nr;
101 rest=b->maxreal % MULT;
102 if (rest != 0)
103 b->maxreal+=MULT-rest;
104 srenew(b->rbuf,b->maxreal);
106 /* Copy pointer */
107 rbuf=b->rbuf+b->nreal;
108 for(i=0; (i<nr); i++)
109 rbuf[i]=r[i];
111 index=b->nreal;
112 b->nreal+=nr;
114 return index;
117 void sum_bin(t_bin *b,t_commrec *cr)
119 int i;
121 for(i=b->nreal; (i<b->maxreal); i++)
122 b->rbuf[i]=0;
123 gmx_sumd(b->maxreal,b->rbuf,cr);
126 void extract_binr(t_bin *b,int index,int nr,real r[])
128 int i;
129 double *rbuf;
131 rbuf = b->rbuf+index;
132 for(i=0; (i<nr); i++)
133 r[i]=rbuf[i];
136 void extract_bind(t_bin *b,int index,int nr,double r[])
138 int i;
139 double *rbuf;
141 rbuf = b->rbuf+index;
142 for(i=0; (i<nr); i++)
143 r[i]=rbuf[i];
146 #ifdef DEBUGRBIN
147 int main(int argc,char *argv[])
149 t_commrec *cr;
150 t_bin *rb;
151 double *r;
152 rvec *v;
153 int k,i,ni,mi,n,m;
155 cr=init_par(&argc,argv);
156 n=strtol(argv[1],NULL,10);
157 m=strtol(argv[2],NULL,10);
158 fprintf(stdlog,"n=%d\n",n);
159 rb=mk_bin();
160 snew(r,n);
161 snew(v,m);
163 for(k=0; (k < 3); k++) {
164 fprintf(stdlog,"\nk=%d\n",k);
165 reset_bin(rb);
167 for(i=0; (i<n); i++)
168 r[i]=i+k;
169 for(i=0; (i<m); i++) {
170 v[i][XX]=4*i+k;
171 v[i][YY]=4*i+k+1;
172 v[i][ZZ]=4*i+k+2;
175 ni=add_bind(stdlog,rb,n,r);
176 mi=add_binr(stdlog,rb,DIM*m,v[0]);
178 sum_bin(rb,cr);
180 extract_bind(rb,ni,n,r);
181 extract_binr(rb,mi,DIM*m,v[0]);
183 for(i=0; (i<n); i++)
184 fprintf(stdlog,"r[%d] = %e\n",i,r[i]);
185 for(i=0; (i<m); i++)
186 fprintf(stdlog,"v[%d] = (%e,%e,%e)\n",i,v[i][XX],v[i][YY],v[i][ZZ]);
188 fflush(stdlog);
190 return 0;
192 #endif