4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
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
34 * GROningen Mixture of Alchemy and Childrens' Stories
36 /* This file is completely threadsafe - keep it that way! */
55 void move_rvecs(FILE *log
,bool bForward
,bool bSum
,
56 int left
,int right
,rvec vecs
[],rvec buf
[],
57 int shift
,t_nsborder
*nsb
,t_nrnb
*nrnb
)
59 int i
,j
,j0
=137,j1
=391;
61 #define next ((cur+1) % nsb->nnodes)
62 #define prev ((cur-1+nsb->nnodes) % nsb->nnodes)
65 cur
=(nsb
->nodeid
+nsb
->shift
) % nsb
->nnodes
;
70 for(i
=0; (i
<shift
); i
++) {
74 j1
=j0
+nsb
->homenr
[prev
];
78 j1
=j0
+nsb
->homenr
[next
];
80 for(j
=j0
; (j
<j1
); j
++) {
84 /* Forward pulse around the ring, to increasing NODE number */
87 gmx_tx_rx_real(right
,vecs
[nsb
->index
[cur
]], nsb
->homenr
[cur
]*DIM
,
88 left
,buf
[nsb
->index
[prev
]],nsb
->homenr
[prev
]*DIM
);
90 gmx_tx_rx_real(right
,vecs
[nsb
->index
[cur
]], nsb
->homenr
[cur
]*DIM
,
91 left
, vecs
[nsb
->index
[prev
]],nsb
->homenr
[prev
]*DIM
);
92 /* Wait for communication to end */
96 /* Backward pulse around the ring, to decreasing NODE number */
99 gmx_tx_rx_real(left
, vecs
[nsb
->index
[cur
]], nsb
->homenr
[cur
]*DIM
,
100 right
,buf
[nsb
->index
[next
]],nsb
->homenr
[next
]*DIM
);
102 gmx_tx_rx_real(left
, vecs
[nsb
->index
[cur
]], nsb
->homenr
[cur
]*DIM
,
103 right
,vecs
[nsb
->index
[next
]],nsb
->homenr
[next
]*DIM
);
104 /* Wait for communication to end */
105 gmx_wait(left
,right
);
108 /* Actual summation */
110 for(j
=j0
; (j
<j1
); j
++) {
111 rvec_inc(vecs
[j
],buf
[j
]);
121 inc_nrnb(nrnb
,eNR_FSUM
,nsum
);
126 void move_x(FILE *log
,
127 int left
,int right
,rvec x
[],t_nsborder
*nsb
,
130 move_rvecs(log
,FALSE
,FALSE
,left
,right
,x
,NULL
,nsb
->shift
,nsb
,nrnb
);
131 move_rvecs(log
,TRUE
, FALSE
,left
,right
,x
,NULL
,nsb
->bshift
,nsb
,nrnb
);
136 void move_f(FILE *log
,
137 int left
,int right
,rvec f
[],rvec fadd
[],
138 t_nsborder
*nsb
,t_nrnb
*nrnb
)
140 move_rvecs(log
,TRUE
, TRUE
,left
,right
,f
,fadd
,nsb
->shift
,nsb
,nrnb
);
141 move_rvecs(log
,FALSE
,TRUE
,left
,right
,f
,fadd
,nsb
->bshift
,nsb
,nrnb
);
146 void move_cgcm(FILE *log
,t_commrec
*cr
,rvec cg_cm
[],int nload
[])
150 #define next ((cur+1) % cr->nnodes)
152 for(i
=0; (i
<cr
->nnodes
-1); i
++) {
153 start
= (cur
== 0) ? 0 : nload
[cur
-1];
154 nr
= nload
[cur
] - start
;
155 gmx_tx(cr
->left
, cg_cm
[start
], nr
*sizeof(cg_cm
[0]));
157 fprintf(log
,"move_cgcm: TX start=%d, nr=%d\n",start
,nr
);
159 start
= (next
== 0) ? 0 : nload
[next
-1];
160 nr
= nload
[next
] - start
;
161 gmx_rx(cr
->right
,cg_cm
[start
], nr
*sizeof(cg_cm
[0]));
163 fprintf(log
,"move_cgcm: RX start=%d, nr=%d\n",start
,nr
);
165 gmx_tx_wait(cr
->left
);
166 gmx_rx_wait(cr
->right
);