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