changed reading hint
[gromacs/adressmacs.git] / src / gmxlib / network.c
blobba17a211e44c3f6ede16fd0450f08119d80cd042
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 2.0
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
16 * Please refer to:
17 * GROMACS: A message-passing parallel molecular dynamics implementation
18 * H.J.C. Berendsen, D. van der Spoel and R. van Drunen
19 * Comp. Phys. Comm. 91, 43-56 (1995)
21 * Also check out our WWW page:
22 * http://md.chem.rug.nl/~gmx
23 * or e-mail to:
24 * gromacs@chem.rug.nl
26 * And Hey:
27 * Green Red Orange Magenta Azure Cyan Skyblue
29 static char *SRCID_network_c = "$Id$";
31 #include "typedefs.h"
32 #include "smalloc.h"
33 #include "network.h"
35 void def_sync_ring(int pid,int nprocs,int left,int right)
37 int i;
38 int tag=0;
40 for (i=0; (i<nprocs); i++) {
41 if (pid == 0) {
42 gmx_txs(right,record(tag));
43 gmx_rxs(left,record(tag));
45 else {
46 gmx_rxs(left,record(tag));
47 gmx_txs(right,record(tag));
52 void def_tx_rx(int send_pid,void *send_buf,int send_bufsize,
53 int rec_pid,void *rec_buf,int rec_bufsize)
55 gmx_tx(send_pid,send_buf,send_bufsize);
56 gmx_rx(rec_pid,rec_buf,rec_bufsize);
57 gmx_tx_wait(send_pid);
58 gmx_rx_wait(rec_pid);
61 void def_wait(int send,int receive)
63 gmx_tx_wait(send);
64 gmx_rx_wait(receive);
67 void def_stat(FILE *fp,char *msg)
69 fprintf(fp,"def_stat: %s (from %s, %d)\n",msg,__FILE__,__LINE__);
72 void def_reset_idle(void)
76 void def_sumd(int nr,double r[],t_commrec *cr)
78 double *buf[2];
79 int NR,bufs,j,i,cur=0;
80 #define next (1-cur)
82 bufs=nr*sizeof(buf[0][0]);
83 #ifdef _amb_
84 bufs+=23;
85 bufs=bufs-(bufs % 24);
86 NR=bufs/sizeof(double);
87 #else
88 NR=nr;
89 #endif
91 snew(buf[0],NR);
92 snew(buf[1],NR);
94 for(i=0; (i<nr); i++)
95 buf[cur][i]=r[i];
96 for(j=0; (j<cr->nprocs-1); j++) {
97 gmx_tx(cr->left,buf[cur],bufs);
98 gmx_rx(cr->right,buf[next],bufs);
99 gmx_tx_wait(cr->left);
100 gmx_rx_wait(cr->right);
101 for(i=0; (i<nr); i++)
102 r[i]+=buf[next][i];
104 cur=next;
106 sfree(buf[1]);
107 sfree(buf[0]);
110 void def_sumf(int nr,float r[],t_commrec *cr)
112 float *buf[2];
113 int NR,bufs,j,i,cur=0;
114 #define next (1-cur)
116 bufs=nr*sizeof(float);
117 #ifdef _amb_
118 bufs+=23;
119 bufs=bufs-(bufs % 24);
120 NR=bufs/sizeof(float);
121 #else
122 NR=nr;
123 #endif
125 snew(buf[0],NR);
126 snew(buf[1],NR);
128 for(i=0; (i<nr); i++)
129 buf[cur][i]=r[i];
130 for(j=0; (j<cr->nprocs-1); j++) {
131 gmx_tx(cr->left,buf[cur],bufs);
132 gmx_rx(cr->right,buf[next],bufs);
133 gmx_wait(cr->left,cr->right);
134 for(i=0; (i<nr); i++)
135 r[i]+=buf[next][i];
137 cur=next;
139 sfree(buf[1]);
140 sfree(buf[0]);
143 void def_sumi(int nr,int r[],t_commrec *cr)
145 int *buf[2];
146 int NR,bufs,j,i,cur=0;
147 #define next (1-cur)
149 bufs=nr*sizeof(int);
150 #ifdef _amb_
151 bufs+=23;
152 bufs=bufs-(bufs % 24);
153 NR=bufs/sizeof(int);
154 #else
155 NR=nr;
156 #endif
158 snew(buf[0],NR);
159 snew(buf[1],NR);
161 for(i=0; (i<nr); i++)
162 buf[cur][i]=r[i];
163 for(j=0; (j<cr->nprocs-1); j++) {
164 gmx_tx(cr->left,buf[cur],bufs);
165 gmx_rx(cr->right,buf[next],bufs);
166 gmx_wait(cr->left,cr->right);
167 for(i=0; (i<nr); i++)
168 r[i]+=buf[next][i];
170 cur=next;
172 sfree(buf[1]);
173 sfree(buf[0]);
176 #ifdef GSUM
177 #include "main.h"
179 int i860main(int argc,char *argv[],FILE *log,
180 int nprocs,int pid,int left,int right)
182 #define MAXBUF 2737
183 real buf[MAXBUF];
184 int i;
185 t_commrec cr;
187 for(i=0; (i<MAXBUF); i++)
188 buf[i]=i;
189 cr.nprocs=nprocs;
190 cr.pid=pid;
191 cr.left=left;
192 cr.right=right;
194 gmx_sum(MAXBUF,buf,&cr);
197 for(i=0; (i<MAXBUF); i++)
198 fprintf(log,"buf[%2d]=%g\n",i,buf[i]);
200 return 0;
202 #endif