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