Fix remaining copyright headers
[alexxy/gromacs.git] / src / gromacs / 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  * Copyright (c) 2013, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <sysstuff.h>
43 #include <string.h>
44 #include "typedefs.h"
45 #include "main.h"
46 #include "mvdata.h"
47 #include "network.h"
48 #include "smalloc.h"
49 #include "gmx_fatal.h"
50 #include "symtab.h"
51 #include "main.h"
52 #include "typedefs.h"
53 #include "vec.h"
54 #include "tgroup.h"
55 #include "nrnb.h"
56 #include "partdec.h"
57
58 void move_rvecs(const t_commrec *cr, gmx_bool bForward, gmx_bool bSum,
59                 rvec vecs[], rvec buf[],
60                 int shift, t_nrnb *nrnb)
61 {
62     int     i, j, j0 = 137, j1 = 391;
63     int     cur, nsum;
64     int    *index;
65 #define next ((cur + 1) % cr->nnodes)
66 #define prev ((cur - 1 + cr->nnodes) % cr->nnodes)
67
68 #define HOMENRI(ind, i) ((ind)[(i)+1] - (ind)[(i)])
69
70     index = pd_index(cr);
71
72     if (bSum)
73     {
74         cur = (cr->nodeid + pd_shift(cr)) % cr->nnodes;
75     }
76     else
77     {
78         cur = cr->nodeid;
79     }
80
81     nsum = 0;
82     for (i = 0; (i < shift); i++)
83     {
84         if (bSum)
85         {
86             if (bForward)
87             {
88                 j0 = index[prev];
89                 j1 = index[prev+1];
90             }
91             else
92             {
93                 j0 = index[next];
94                 j1 = index[next+1];
95             }
96             for (j = j0; (j < j1); j++)
97             {
98                 clear_rvec(buf[j]);
99             }
100         }
101         /* Forward pulse around the ring, to increasing NODE number */
102         if (bForward)
103         {
104             if (bSum)
105             {
106                 gmx_tx_rx_real(cr,
107                                GMX_RIGHT, vecs[index[cur ]], HOMENRI(index, cur )*DIM,
108                                GMX_LEFT, buf [index[prev]], HOMENRI(index, prev)*DIM);
109             }
110             else
111             {
112                 gmx_tx_rx_real(cr,
113                                GMX_RIGHT, vecs[index[cur ]], HOMENRI(index, cur )*DIM,
114                                GMX_LEFT, vecs[index[prev]], HOMENRI(index, prev)*DIM);
115             }
116             /* Wait for communication to end */
117             gmx_wait(cr);
118         }
119
120         /* Backward pulse around the ring, to decreasing NODE number */
121         else
122         {
123             if (bSum)
124             {
125                 gmx_tx_rx_real(cr,
126                                GMX_LEFT, vecs[index[cur ]], HOMENRI(index, cur )*DIM,
127                                GMX_RIGHT, buf [index[next]], HOMENRI(index, next)*DIM);
128             }
129             else
130             {
131                 gmx_tx_rx_real(cr,
132                                GMX_LEFT, vecs[index[cur ]], HOMENRI(index, cur )*DIM,
133                                GMX_RIGHT, vecs[index[next]], HOMENRI(index, next)*DIM);
134             }
135             /* Wait for communication to end */
136             gmx_wait(cr);
137         }
138
139         /* Actual summation */
140         if (bSum)
141         {
142             for (j = j0; (j < j1); j++)
143             {
144                 rvec_inc(vecs[j], buf[j]);
145             }
146             nsum += (j1-j0);
147         }
148         if (bForward)
149         {
150             cur = prev;
151         }
152         else
153         {
154             cur = next;
155         }
156     }
157     if (nsum > 0)
158     {
159         inc_nrnb(nrnb, eNR_FSUM, nsum);
160     }
161 #undef next
162 #undef prev
163 }
164
165
166 void move_x(const t_commrec *cr,
167             rvec x[], t_nrnb *nrnb)
168 {
169     move_rvecs(cr, FALSE, FALSE, x, NULL, pd_shift(cr), nrnb);
170     move_rvecs(cr, TRUE, FALSE, x, NULL, pd_bshift(cr), nrnb);
171
172     where();
173 }
174
175
176 void move_f(const t_commrec *cr,
177             rvec f[], rvec fadd[],
178             t_nrnb *nrnb)
179 {
180     move_rvecs(cr, TRUE, TRUE, f, fadd, pd_shift(cr), nrnb);
181     move_rvecs(cr, FALSE, TRUE, f, fadd, pd_bshift(cr), nrnb);
182
183     where();
184 }
185
186 void move_cgcm(FILE gmx_unused *log, const t_commrec *cr, rvec cg_cm[])
187 {
188     int  i, start, nr;
189     int  cur = cr->nodeid;
190     int *cgindex;
191
192 #define next ((cur+1) % cr->nnodes)
193
194     cgindex = pd_cgindex(cr);
195
196     for (i = 0; (i < cr->nnodes-1); i++)
197     {
198         start = cgindex[cur];
199         nr    = cgindex[cur+1] - start;
200
201         gmx_tx(cr, GMX_LEFT, cg_cm[start], nr*sizeof(cg_cm[0]));
202 #ifdef DEBUG
203         fprintf(log, "move_cgcm: TX start=%d, nr=%d\n", start, nr);
204 #endif
205         start = cgindex[next];
206         nr    = cgindex[next+1] - start;
207
208         gmx_rx(cr, GMX_RIGHT, cg_cm[start], nr*sizeof(cg_cm[0]));
209 #ifdef DEBUG
210         fprintf(log, "move_cgcm: RX start=%d, nr=%d\n", start, nr);
211 #endif
212         gmx_tx_wait(cr);
213         gmx_rx_wait(cr);
214
215         if (debug)
216         {
217             fprintf(debug, "cgcm[0][XX] %f\n", cg_cm[0][XX]);
218         }
219
220         cur = next;
221     }
222 #undef next
223 }