Merge t_pbc definition into pbc.h
[alexxy/gromacs.git] / src / gromacs / mdlib / calcvir.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,2014, 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 "force.h"
43 #include "macros.h"
44
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/pbcutil/mshift.h"
47 #include "pbc.h"
48
49 #define XXXX    0
50 #define XXYY    1
51 #define XXZZ    2
52 #define YYXX    3
53 #define YYYY    4
54 #define YYZZ    5
55 #define ZZXX    6
56 #define ZZYY    7
57 #define ZZZZ    8
58
59 static void upd_vir(rvec vir, real dvx, real dvy, real dvz)
60 {
61     vir[XX] -= 0.5*dvx;
62     vir[YY] -= 0.5*dvy;
63     vir[ZZ] -= 0.5*dvz;
64 }
65
66 void calc_vir(int nxf, rvec x[], rvec f[], tensor vir,
67               gmx_bool bScrewPBC, matrix box)
68 {
69     int      i, isx;
70     double   dvxx = 0, dvxy = 0, dvxz = 0, dvyx = 0, dvyy = 0, dvyz = 0, dvzx = 0, dvzy = 0, dvzz = 0;
71
72     for (i = 0; (i < nxf); i++)
73     {
74         dvxx += x[i][XX]*f[i][XX];
75         dvxy += x[i][XX]*f[i][YY];
76         dvxz += x[i][XX]*f[i][ZZ];
77
78         dvyx += x[i][YY]*f[i][XX];
79         dvyy += x[i][YY]*f[i][YY];
80         dvyz += x[i][YY]*f[i][ZZ];
81
82         dvzx += x[i][ZZ]*f[i][XX];
83         dvzy += x[i][ZZ]*f[i][YY];
84         dvzz += x[i][ZZ]*f[i][ZZ];
85
86         if (bScrewPBC)
87         {
88             isx = IS2X(i);
89             /* We should correct all odd x-shifts, but the range of isx is -2 to 2 */
90             if (isx == 1 || isx == -1)
91             {
92                 dvyy += box[YY][YY]*f[i][YY];
93                 dvyz += box[YY][YY]*f[i][ZZ];
94
95                 dvzy += box[ZZ][ZZ]*f[i][YY];
96                 dvzz += box[ZZ][ZZ]*f[i][ZZ];
97             }
98         }
99     }
100
101     upd_vir(vir[XX], dvxx, dvxy, dvxz);
102     upd_vir(vir[YY], dvyx, dvyy, dvyz);
103     upd_vir(vir[ZZ], dvzx, dvzy, dvzz);
104 }
105
106
107 static void lo_fcv(int i0, int i1,
108                    real x[], real f[], tensor vir,
109                    int is[], real box[], gmx_bool bTriclinic)
110 {
111     int      i, i3, tx, ty, tz;
112     real     xx, yy, zz;
113     real     dvxx = 0, dvxy = 0, dvxz = 0, dvyx = 0, dvyy = 0, dvyz = 0, dvzx = 0, dvzy = 0, dvzz = 0;
114
115     if (bTriclinic)
116     {
117         for (i = i0; (i < i1); i++)
118         {
119             i3 = DIM*i;
120             tx = is[i3+XX];
121             ty = is[i3+YY];
122             tz = is[i3+ZZ];
123
124             xx    = x[i3+XX]-tx*box[XXXX]-ty*box[YYXX]-tz*box[ZZXX];
125             dvxx += xx*f[i3+XX];
126             dvxy += xx*f[i3+YY];
127             dvxz += xx*f[i3+ZZ];
128
129             yy    = x[i3+YY]-ty*box[YYYY]-tz*box[ZZYY];
130             dvyx += yy*f[i3+XX];
131             dvyy += yy*f[i3+YY];
132             dvyz += yy*f[i3+ZZ];
133
134             zz    = x[i3+ZZ]-tz*box[ZZZZ];
135             dvzx += zz*f[i3+XX];
136             dvzy += zz*f[i3+YY];
137             dvzz += zz*f[i3+ZZ];
138         }
139     }
140     else
141     {
142         for (i = i0; (i < i1); i++)
143         {
144             i3 = DIM*i;
145             tx = is[i3+XX];
146             ty = is[i3+YY];
147             tz = is[i3+ZZ];
148
149             xx    = x[i3+XX]-tx*box[XXXX];
150             dvxx += xx*f[i3+XX];
151             dvxy += xx*f[i3+YY];
152             dvxz += xx*f[i3+ZZ];
153
154             yy    = x[i3+YY]-ty*box[YYYY];
155             dvyx += yy*f[i3+XX];
156             dvyy += yy*f[i3+YY];
157             dvyz += yy*f[i3+ZZ];
158
159             zz    = x[i3+ZZ]-tz*box[ZZZZ];
160             dvzx += zz*f[i3+XX];
161             dvzy += zz*f[i3+YY];
162             dvzz += zz*f[i3+ZZ];
163         }
164     }
165
166     upd_vir(vir[XX], dvxx, dvxy, dvxz);
167     upd_vir(vir[YY], dvyx, dvyy, dvyz);
168     upd_vir(vir[ZZ], dvzx, dvzy, dvzz);
169 }
170
171 static void lo_fcv2(int i0, int i1,
172                     rvec x[], rvec f[], tensor vir,
173                     ivec is[], matrix box, gmx_bool bTriclinic)
174 {
175     int      i, gg, tx, ty, tz;
176     real     xx, yy, zz;
177     real     dvxx = 0, dvxy = 0, dvxz = 0, dvyx = 0, dvyy = 0, dvyz = 0, dvzx = 0, dvzy = 0, dvzz = 0;
178
179     if (bTriclinic)
180     {
181         for (i = i0, gg = 0; (i < i1); i++, gg++)
182         {
183             tx = is[gg][XX];
184             ty = is[gg][YY];
185             tz = is[gg][ZZ];
186
187             xx    = x[i][XX]-tx*box[XX][XX]-ty*box[YY][XX]-tz*box[ZZ][XX];
188             dvxx += xx*f[i][XX];
189             dvxy += xx*f[i][YY];
190             dvxz += xx*f[i][ZZ];
191
192             yy    = x[i][YY]-ty*box[YY][YY]-tz*box[ZZ][YY];
193             dvyx += yy*f[i][XX];
194             dvyy += yy*f[i][YY];
195             dvyz += yy*f[i][ZZ];
196
197             zz    = x[i][ZZ]-tz*box[ZZ][ZZ];
198             dvzx += zz*f[i][XX];
199             dvzy += zz*f[i][YY];
200             dvzz += zz*f[i][ZZ];
201         }
202     }
203     else
204     {
205         for (i = i0, gg = 0; (i < i1); i++, gg++)
206         {
207             tx = is[gg][XX];
208             ty = is[gg][YY];
209             tz = is[gg][ZZ];
210
211             xx    = x[i][XX]-tx*box[XX][XX];
212             dvxx += xx*f[i][XX];
213             dvxy += xx*f[i][YY];
214             dvxz += xx*f[i][ZZ];
215
216             yy    = x[i][YY]-ty*box[YY][YY];
217             dvyx += yy*f[i][XX];
218             dvyy += yy*f[i][YY];
219             dvyz += yy*f[i][ZZ];
220
221             zz    = x[i][ZZ]-tz*box[ZZ][ZZ];
222             dvzx += zz*f[i][XX];
223             dvzy += zz*f[i][YY];
224             dvzz += zz*f[i][ZZ];
225         }
226     }
227
228     upd_vir(vir[XX], dvxx, dvxy, dvxz);
229     upd_vir(vir[YY], dvyx, dvyy, dvyz);
230     upd_vir(vir[ZZ], dvzx, dvzy, dvzz);
231 }
232
233 void f_calc_vir(int i0, int i1, rvec x[], rvec f[], tensor vir,
234                 t_graph *g, matrix box)
235 {
236     int start, end;
237
238     if (g && (g->nnodes > 0))
239     {
240         /* Calculate virial for bonded forces only when they belong to
241          * this node.
242          */
243         start = max(i0, g->at_start);
244         end   = min(i1, g->at_end);
245 #ifdef SAFE
246         lo_fcv2(start, end, x, f, vir, g->ishift, box, TRICLINIC(box));
247 #else
248         lo_fcv(start, end, x[0], f[0], vir, g->ishift[0], box[0], TRICLINIC(box));
249 #endif
250
251         /* If not all atoms are bonded, calculate their virial contribution
252          * anyway, without shifting back their coordinates.
253          * Note the nifty pointer arithmetic...
254          */
255         if (start > i0)
256         {
257             calc_vir(start-i0, x + i0, f + i0, vir, FALSE, box);
258         }
259         if (end < i1)
260         {
261             calc_vir(i1-end, x + end, f + end, vir, FALSE, box);
262         }
263     }
264     else
265     {
266         calc_vir(i1-i0, x + i0, f + i0, vir, FALSE, box);
267     }
268 }