9f91cc217eb1c9e73db1ad109dfccdc2e08d708c
[alexxy/gromacs.git] / src / mdlib / mvxvf.c
1 /*
2  * 
3  *                This source code is part of
4  * 
5  *                 G   R   O   M   A   C   S
6  * 
7  *          GROningen MAchine for Chemical Simulations
8  * 
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  * 
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  * 
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * GROningen Mixture of Alchemy and Childrens' Stories
34  */
35 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <sysstuff.h>
41 #include <string.h>
42 #include "typedefs.h"
43 #include "main.h"
44 #include "mvdata.h"
45 #include "network.h"
46 #include "smalloc.h"
47 #include "gmx_fatal.h"
48 #include "symtab.h"
49 #include "main.h"
50 #include "typedefs.h"
51 #include "vec.h"
52 #include "tgroup.h"
53 #include "nrnb.h"
54 #include "partdec.h"
55
56 void move_rvecs(const t_commrec *cr,gmx_bool bForward,gmx_bool bSum,
57                 int left,int right,rvec vecs[],rvec buf[],
58                 int shift,t_nrnb *nrnb)
59 {
60     int    i,j,j0=137,j1=391;
61     int    cur,nsum;
62     int    *index;
63 #define next ((cur + 1) % cr->nnodes)
64 #define prev ((cur - 1 + cr->nnodes) % cr->nnodes)
65     
66 #define HOMENRI(ind,i) ((ind)[(i)+1] - (ind)[(i)])
67     
68     index = pd_index(cr);
69     
70     if (bSum)
71         cur = (cr->nodeid + pd_shift(cr)) % cr->nnodes;
72     else
73         cur = cr->nodeid;
74     
75     nsum=0;
76     for(i=0; (i<shift); i++) {
77         if (bSum) {
78             if (bForward) {
79                 j0 = index[prev];
80                 j1 = index[prev+1];
81             }
82             else {
83                 j0 = index[next];
84                 j1 = index[next+1];
85             }
86             for(j=j0; (j<j1); j++) {
87                 clear_rvec(buf[j]);
88             }
89         }
90         /* Forward pulse around the ring, to increasing NODE number */
91         if (bForward) {
92             if (bSum)
93                 gmx_tx_rx_real(cr,
94                                GMX_RIGHT,vecs[index[cur ]],HOMENRI(index,cur )*DIM,
95                                GMX_LEFT, buf [index[prev]],HOMENRI(index,prev)*DIM);
96             else
97                 gmx_tx_rx_real(cr,
98                                GMX_RIGHT,vecs[index[cur ]],HOMENRI(index,cur )*DIM,
99                                GMX_LEFT, vecs[index[prev]],HOMENRI(index,prev)*DIM);
100             /* Wait for communication to end */
101             gmx_wait(cr, right,left);
102         }
103         
104         /* Backward pulse around the ring, to decreasing NODE number */
105         else {
106             if (bSum)
107                 gmx_tx_rx_real(cr,
108                                GMX_LEFT, vecs[index[cur ]],HOMENRI(index,cur )*DIM,
109                                GMX_RIGHT,buf [index[next]],HOMENRI(index,next)*DIM);
110             else
111                 gmx_tx_rx_real(cr,
112                                GMX_LEFT, vecs[index[cur ]],HOMENRI(index,cur )*DIM,
113                                GMX_RIGHT,vecs[index[next]],HOMENRI(index,next)*DIM);
114             /* Wait for communication to end */
115             gmx_wait(cr, left,right);
116         }
117         
118         /* Actual summation */
119         if (bSum) {
120             for(j=j0; (j<j1); j++) {
121                 rvec_inc(vecs[j],buf[j]);
122             }
123             nsum+=(j1-j0);
124         }
125         if (bForward) 
126             cur=prev;
127         else
128             cur=next;
129     }  
130     if (nsum > 0)
131         inc_nrnb(nrnb,eNR_FSUM,nsum);
132 #undef next
133 #undef prev
134 }
135
136
137 void move_reals(const t_commrec *cr,gmx_bool bForward,gmx_bool bSum,
138                 int left,int right,real reals[],real buf[],
139                 int shift,t_nrnb *nrnb)
140 {
141     int    i,j,j0=137,j1=391;
142     int    cur,nsum;
143     int    *index;
144 #define next ((cur + 1) % cr->nnodes)
145 #define prev ((cur - 1 + cr->nnodes) % cr->nnodes)
146     
147 #define HOMENRI(ind,i) ((ind)[(i)+1] - (ind)[(i)])
148     
149     index = pd_index(cr);
150     
151     if (bSum)
152     {
153         cur = (cr->nodeid + pd_shift(cr)) % cr->nnodes;
154     }
155     else
156     {
157         cur = cr->nodeid;
158     }
159     
160     nsum=0;
161     for(i=0; (i<shift); i++) 
162     {
163         if (bSum) 
164         {
165             if (bForward) 
166             {
167                 j0 = index[prev];
168                 j1 = index[prev+1];
169             }
170             else 
171             {
172                 j0 = index[next];
173                 j1 = index[next+1];
174             }
175             for(j=j0; (j<j1); j++) 
176             {
177                 buf[j] = 0.0;
178             }
179         }
180         /* Forward pulse around the ring, to increasing NODE number */
181         if (bForward) 
182         {
183             if (bSum)
184             {    gmx_tx_rx_real(cr,
185                                 GMX_RIGHT,reals+index[cur ],HOMENRI(index,cur ),
186                                 GMX_LEFT, buf+index[prev],HOMENRI(index,prev));
187             }
188             else
189             {
190                 gmx_tx_rx_real(cr,
191                                GMX_RIGHT,reals+index[cur ],HOMENRI(index,cur ),
192                                GMX_LEFT, reals+index[prev],HOMENRI(index,prev));
193             }
194             /* Wait for communication to end */
195             gmx_wait(cr, right,left);
196         }
197         else
198         {
199                  /* Backward pulse around the ring, to decreasing NODE number */
200             if (bSum)
201             {    
202                 gmx_tx_rx_real(cr,
203                                GMX_LEFT, reals+index[cur ],HOMENRI(index,cur ),
204                                GMX_RIGHT,buf+index[next],HOMENRI(index,next));
205             }
206             else
207             {
208                 gmx_tx_rx_real(cr,
209                                GMX_LEFT, reals+index[cur ],HOMENRI(index,cur ),
210                                GMX_RIGHT,reals+index[next],HOMENRI(index,next));
211                 /* Wait for communication to end */
212             }
213             gmx_wait(cr, left,right);
214         }
215
216         /* Actual summation */
217         if (bSum) 
218         {
219             for(j=j0; (j<j1); j++) 
220             {
221                 reals[j] += buf[j];
222             }
223             nsum+=(j1-j0);
224         }
225         if (bForward) 
226         {
227             cur=prev;
228         }
229         else
230         {
231             cur=next;
232         }
233     }  
234     
235     if (nsum > 0)
236     {
237         inc_nrnb(nrnb,eNR_FSUM,nsum/3);
238     }
239 #undef next
240 #undef prev
241 }
242
243 void move_x(FILE *log,const t_commrec *cr,
244             int left,int right,rvec x[],
245             t_nrnb *nrnb)
246 {
247   move_rvecs(cr,FALSE,FALSE,left,right,x,NULL,pd_shift(cr),nrnb);
248   move_rvecs(cr,TRUE, FALSE,left,right,x,NULL,pd_bshift(cr),nrnb);
249
250   where();
251 }
252
253 void move_rborn(FILE *log,const t_commrec *cr,
254                 int left,int right,real rborn[],
255                 t_nrnb *nrnb)
256 {
257     move_reals(cr,FALSE,FALSE,left,right,rborn,NULL,pd_shift(cr),nrnb);
258     move_reals(cr,TRUE, FALSE,left,right,rborn,NULL,pd_bshift(cr),nrnb);
259     
260     where();
261 }
262
263 void move_f(FILE *log,const t_commrec *cr,
264             int left,int right,rvec f[],rvec fadd[],
265             t_nrnb *nrnb)
266 {
267   move_rvecs(cr,TRUE, TRUE,left,right,f,fadd,pd_shift(cr),nrnb);
268   move_rvecs(cr,FALSE,TRUE,left,right,f,fadd,pd_bshift(cr),nrnb);
269
270   where();
271 }
272
273 void move_gpol(FILE *log,const t_commrec *cr,
274                int left,int right,real gpol[],real gpol_add[],
275                t_nrnb *nrnb)
276 {
277     move_reals(cr,TRUE, TRUE,left,right,gpol,gpol_add,pd_shift(cr),nrnb);
278     move_reals(cr,FALSE,TRUE,left,right,gpol,gpol_add,pd_bshift(cr),nrnb);
279     
280     where();
281 }
282
283
284 void move_cgcm(FILE *log,const t_commrec *cr,rvec cg_cm[])
285 {
286   int i,start,nr;
287   int cur=cr->nodeid;
288   int *cgindex;
289
290 #define next ((cur+1) % cr->nnodes)
291   
292   cgindex = pd_cgindex(cr);
293
294   for(i=0; (i<cr->nnodes-1); i++) {
295     start = cgindex[cur];
296     nr    = cgindex[cur+1] - start;
297
298     gmx_tx(cr,GMX_LEFT, cg_cm[start], nr*sizeof(cg_cm[0]));
299 #ifdef DEBUG
300     fprintf(log,"move_cgcm: TX start=%d, nr=%d\n",start,nr);
301 #endif    
302     start = cgindex[next];
303     nr    = cgindex[next+1] - start;
304
305     gmx_rx(cr,GMX_RIGHT,cg_cm[start], nr*sizeof(cg_cm[0]));
306 #ifdef DEBUG
307     fprintf(log,"move_cgcm: RX start=%d, nr=%d\n",start,nr);
308 #endif    
309     gmx_tx_wait(cr,GMX_LEFT);
310     gmx_rx_wait(cr,GMX_RIGHT);
311
312     if (debug)
313       fprintf(debug,"cgcm[0][XX] %f\n",cg_cm[0][XX]);
314
315     cur=next;    
316   }
317 #undef next
318 }
319
320