Added conditional inclusion of config.h to source files
[gromacs.git] / src / gmxlib / rbin.c
blobd50a191c9f397d55031649805e44360763b92210
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.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * GROningen Mixture of Alchemy and Childrens' Stories
36 /* This file is completely threadsafe - keep it that way! */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
41 #include "typedefs.h"
42 #include "main.h"
43 #include "network.h"
44 #include "rbin.h"
45 #include "smalloc.h"
47 t_bin *mk_bin(void)
49 t_bin *b;
51 snew(b,1);
53 return b;
56 void reset_bin(t_bin *b)
58 b->nreal = 0;
61 int add_binr(FILE *log,t_bin *b,int nr,real r[])
63 #define MULT 4
64 int i,rest,index;
65 double *rbuf;
67 if (b->nreal+nr > b->maxreal) {
68 #ifdef DEBUG
69 fprintf(log,"Before: maxreal=%d, nr=%d, nreal=%d\n",
70 b->maxreal,nr,b->nreal);
71 #endif
72 b->maxreal=b->nreal+nr;
73 rest=b->maxreal % MULT;
74 if (rest != 0)
75 b->maxreal+=MULT-rest;
76 srenew(b->rbuf,b->maxreal);
77 #ifdef DEBUG
78 fprintf(log,"After: maxreal=%d, nr=%d, nreal=%d\n",
79 b->maxreal,nr,b->nreal+nr);
80 #endif
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(FILE *log,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 #ifdef DEBUG
101 fprintf(log,"Before: maxreal=%d, nr=%d, nreal=%d\n",
102 b->maxreal,nr,b->nreal);
103 #endif
104 b->maxreal=b->nreal+nr;
105 rest=b->maxreal % MULT;
106 if (rest != 0)
107 b->maxreal+=MULT-rest;
108 srenew(b->rbuf,b->maxreal);
109 #ifdef DEBUG
110 fprintf(log,"After: maxreal=%d, nr=%d, nreal=%d\n",
111 b->maxreal,nr,b->nreal+nr);
112 #endif
114 /* Copy pointer */
115 rbuf=b->rbuf+b->nreal;
116 for(i=0; (i<nr); i++)
117 rbuf[i]=r[i];
119 index=b->nreal;
120 b->nreal+=nr;
122 return index;
125 void sum_bin(t_bin *b,t_commrec *cr)
127 int i;
129 for(i=b->nreal; (i<b->maxreal); i++)
130 b->rbuf[i]=0;
131 gmx_sumd(b->maxreal,b->rbuf,cr);
134 void extract_binr(t_bin *b,int index,int nr,real r[])
136 int i;
137 double *rbuf;
139 rbuf = b->rbuf+index;
140 for(i=0; (i<nr); i++)
141 r[i]=rbuf[i];
144 void extract_bind(t_bin *b,int index,int nr,double r[])
146 int i;
147 double *rbuf;
149 rbuf = b->rbuf+index;
150 for(i=0; (i<nr); i++)
151 r[i]=rbuf[i];
154 #ifdef DEBUGRBIN
155 int main(int argc,char *argv[])
157 t_commrec *cr;
158 t_bin *rb;
159 double *r;
160 rvec *v;
161 int k,i,ni,mi,n,m;
163 cr=init_par(&argc,argv);
164 n=atoi(argv[1]);
165 m=atoi(argv[2]);
166 fprintf(stdlog,"n=%d\n",n);
167 rb=mk_bin();
168 snew(r,n);
169 snew(v,m);
171 for(k=0; (k < 3); k++) {
172 fprintf(stdlog,"\nk=%d\n",k);
173 reset_bin(rb);
175 for(i=0; (i<n); i++)
176 r[i]=i+k;
177 for(i=0; (i<m); i++) {
178 v[i][XX]=4*i+k;
179 v[i][YY]=4*i+k+1;
180 v[i][ZZ]=4*i+k+2;
183 ni=add_bind(stdlog,rb,n,r);
184 mi=add_binr(stdlog,rb,DIM*m,v[0]);
186 sum_bin(rb,cr);
188 extract_bind(rb,ni,n,r);
189 extract_binr(rb,mi,DIM*m,v[0]);
191 for(i=0; (i<n); i++)
192 fprintf(stdlog,"r[%d] = %e\n",i,r[i]);
193 for(i=0; (i<m); i++)
194 fprintf(stdlog,"v[%d] = (%e,%e,%e)\n",i,v[i][XX],v[i][YY],v[i][ZZ]);
196 fclose(stdlog);
198 return 0;
200 #endif