Partial commit of the project to remove all static variables.
[gromacs.git] / src / gmxlib / rbin.c
blob65bbe6eb3260a18bbc7828631be0884776abab6b
1 /*
2 * $Id$
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.1
11 * Copyright (c) 1991-2001, University of Groningen, The Netherlands
12 * This program is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU General Public License
14 * as published by the Free Software Foundation; either version 2
15 * of the License, or (at your option) any later version.
17 * If you want to redistribute modifications, please consider that
18 * scientific software is very special. Version control is crucial -
19 * bugs must be traceable. We will be happy to consider code for
20 * inclusion in the official distribution, but derived work must not
21 * be called official GROMACS. Details are found in the README & COPYING
22 * files - if they are missing, get the official version at www.gromacs.org.
24 * To help us fund GROMACS development, we humbly ask that you cite
25 * the papers on the package - you can find them in the top README file.
27 * For more info, check our website at http://www.gromacs.org
29 * And Hey:
30 * Gnomes, ROck Monsters And Chili Sauce
33 /* This file is completely threadsafe - keep it that way! */
35 #include "typedefs.h"
36 #include "main.h"
37 #include "network.h"
38 #include "rbin.h"
39 #include "smalloc.h"
41 t_bin *mk_bin(void)
43 t_bin *b;
45 snew(b,1);
47 return b;
50 void reset_bin(t_bin *b)
52 b->nreal = 0;
55 int add_binr(FILE *log,t_bin *b,int nr,real r[])
57 #define MULT 4
58 int i,rest,index;
59 double *rbuf;
61 if (b->nreal+nr > b->maxreal) {
62 #ifdef DEBUG
63 fprintf(log,"Before: maxreal=%d, nr=%d, nreal=%d\n",
64 b->maxreal,nr,b->nreal);
65 #endif
66 b->maxreal=b->nreal+nr;
67 rest=b->maxreal % MULT;
68 if (rest != 0)
69 b->maxreal+=MULT-rest;
70 srenew(b->rbuf,b->maxreal);
71 #ifdef DEBUG
72 fprintf(log,"After: maxreal=%d, nr=%d, nreal=%d\n",
73 b->maxreal,nr,b->nreal+nr);
74 #endif
76 /* Copy pointer */
77 rbuf=b->rbuf+b->nreal;
78 for(i=0; (i<nr); i++)
79 rbuf[i]=r[i];
81 index=b->nreal;
82 b->nreal+=nr;
84 return index;
87 int add_bind(FILE *log,t_bin *b,int nr,double r[])
89 #define MULT 4
90 int i,rest,index;
91 double *rbuf;
93 if (b->nreal+nr > b->maxreal) {
94 #ifdef DEBUG
95 fprintf(log,"Before: maxreal=%d, nr=%d, nreal=%d\n",
96 b->maxreal,nr,b->nreal);
97 #endif
98 b->maxreal=b->nreal+nr;
99 rest=b->maxreal % MULT;
100 if (rest != 0)
101 b->maxreal+=MULT-rest;
102 srenew(b->rbuf,b->maxreal);
103 #ifdef DEBUG
104 fprintf(log,"After: maxreal=%d, nr=%d, nreal=%d\n",
105 b->maxreal,nr,b->nreal+nr);
106 #endif
108 /* Copy pointer */
109 rbuf=b->rbuf+b->nreal;
110 for(i=0; (i<nr); i++)
111 rbuf[i]=r[i];
113 index=b->nreal;
114 b->nreal+=nr;
116 return index;
119 void sum_bin(t_bin *b,t_commrec *cr)
121 int i;
123 for(i=b->nreal; (i<b->maxreal); i++)
124 b->rbuf[i]=0;
125 gmx_sumd(b->maxreal,b->rbuf,cr);
128 void extract_binr(t_bin *b,int index,int nr,real r[])
130 int i;
131 double *rbuf;
133 rbuf = b->rbuf+index;
134 for(i=0; (i<nr); i++)
135 r[i]=rbuf[i];
138 void extract_bind(t_bin *b,int index,int nr,double r[])
140 int i;
141 double *rbuf;
143 rbuf = b->rbuf+index;
144 for(i=0; (i<nr); i++)
145 r[i]=rbuf[i];
148 #ifdef DEBUGRBIN
149 int main(int argc,char *argv[])
151 t_commrec *cr;
152 t_bin *rb;
153 double *r;
154 rvec *v;
155 int k,i,ni,mi,n,m;
157 cr=init_par(&argc,argv);
158 n=atoi(argv[1]);
159 m=atoi(argv[2]);
160 fprintf(stdlog,"n=%d\n",n);
161 rb=mk_bin();
162 snew(r,n);
163 snew(v,m);
165 for(k=0; (k < 3); k++) {
166 fprintf(stdlog,"\nk=%d\n",k);
167 reset_bin(rb);
169 for(i=0; (i<n); i++)
170 r[i]=i+k;
171 for(i=0; (i<m); i++) {
172 v[i][XX]=4*i+k;
173 v[i][YY]=4*i+k+1;
174 v[i][ZZ]=4*i+k+2;
177 ni=add_bind(stdlog,rb,n,r);
178 mi=add_binr(stdlog,rb,DIM*m,v[0]);
180 sum_bin(rb,cr);
182 extract_bind(rb,ni,n,r);
183 extract_binr(rb,mi,DIM*m,v[0]);
185 for(i=0; (i<n); i++)
186 fprintf(stdlog,"r[%d] = %e\n",i,r[i]);
187 for(i=0; (i<m); i++)
188 fprintf(stdlog,"v[%d] = (%e,%e,%e)\n",i,v[i][XX],v[i][YY],v[i][ZZ]);
190 fclose(stdlog);
192 return 0;
194 #endif