4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
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
27 * Green Red Orange Magenta Azure Cyan Skyblue
29 static char *SRCID_network_c
= "$Id$";
35 void def_sync_ring(int pid
,int nprocs
,int left
,int right
)
40 for (i
=0; (i
<nprocs
); i
++) {
42 gmx_txs(right
,record(tag
));
43 gmx_rxs(left
,record(tag
));
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
);
61 void def_wait(int send
,int 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
)
79 int NR
,bufs
,j
,i
,cur
=0;
82 bufs
=nr
*sizeof(buf
[0][0]);
85 bufs
=bufs
-(bufs
% 24);
86 NR
=bufs
/sizeof(double);
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
++)
110 void def_sumf(int nr
,float r
[],t_commrec
*cr
)
113 int NR
,bufs
,j
,i
,cur
=0;
116 bufs
=nr
*sizeof(float);
119 bufs
=bufs
-(bufs
% 24);
120 NR
=bufs
/sizeof(float);
128 for(i
=0; (i
<nr
); 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
++)
143 void def_sumi(int nr
,int r
[],t_commrec
*cr
)
146 int NR
,bufs
,j
,i
,cur
=0;
152 bufs
=bufs
-(bufs
% 24);
161 for(i
=0; (i
<nr
); 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
++)
179 int i860main(int argc
,char *argv
[],FILE *log
,
180 int nprocs
,int pid
,int left
,int right
)
187 for(i
=0; (i
<MAXBUF
); i
++)
194 gmx_sum(MAXBUF
,buf
,&cr
);
197 for(i=0; (i<MAXBUF); i++)
198 fprintf(log,"buf[%2d]=%g\n",i,buf[i]);