a9d9f6b28dcb7bcefa42c0619546e3e47002dcb7
[alexxy/gromacs.git] / src / gromacs / mdlib / calcvir.cpp
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,2015,2016,2018,2019, 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 "calcvir.h"
41
42 #include "config.h" /* for GMX_MAX_OPENMP_THREADS */
43
44 #include <algorithm>
45
46 #include "gromacs/math/vec.h"
47 #include "gromacs/math/vectypes.h"
48 #include "gromacs/mdlib/gmx_omp_nthreads.h"
49 #include "gromacs/pbcutil/ishift.h"
50 #include "gromacs/pbcutil/mshift.h"
51 #include "gromacs/pbcutil/pbc.h"
52 #include "gromacs/utility/gmxassert.h"
53
54 #define XXXX 0
55 #define XXYY 1
56 #define XXZZ 2
57 #define YYXX 3
58 #define YYYY 4
59 #define YYZZ 5
60 #define ZZXX 6
61 #define ZZYY 7
62 #define ZZZZ 8
63
64 static void upd_vir(rvec vir, real dvx, real dvy, real dvz)
65 {
66     vir[XX] -= 0.5 * dvx;
67     vir[YY] -= 0.5 * dvy;
68     vir[ZZ] -= 0.5 * dvz;
69 }
70
71 static void calc_x_times_f(int nxf, const rvec x[], const rvec f[], gmx_bool bScrewPBC, const matrix box, matrix x_times_f)
72 {
73     clear_mat(x_times_f);
74
75     for (int i = 0; i < nxf; i++)
76     {
77         for (int d = 0; d < DIM; d++)
78         {
79             for (int n = 0; n < DIM; n++)
80             {
81                 x_times_f[d][n] += x[i][d] * f[i][n];
82             }
83         }
84
85         if (bScrewPBC)
86         {
87             int isx = IS2X(i);
88             /* We should correct all odd x-shifts, but the range of isx is -2 to 2 */
89             if (isx == 1 || isx == -1)
90             {
91                 for (int d = 0; d < DIM; d++)
92                 {
93                     for (int n = 0; n < DIM; n++)
94                     {
95                         x_times_f[d][n] += box[d][d] * f[i][n];
96                     }
97                 }
98             }
99         }
100     }
101 }
102
103 void calc_vir(int nxf, const rvec x[], const rvec f[], tensor vir, bool bScrewPBC, const matrix box)
104 {
105     matrix x_times_f;
106
107     int nthreads = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, nxf * 9);
108
109     GMX_ASSERT(nthreads >= 1, "Avoids uninitialized x_times_f (warning)");
110
111     if (nthreads == 1)
112     {
113         calc_x_times_f(nxf, x, f, bScrewPBC, box, x_times_f);
114     }
115     else
116     {
117         /* Use a buffer on the stack for storing thread-local results.
118          * We use 2 extra elements (=18 reals) per thread to separate thread
119          * local data by at least a cache line. Element 0 is not used.
120          */
121         matrix xf_buf[GMX_OPENMP_MAX_THREADS * 3];
122
123 #pragma omp parallel for num_threads(nthreads) schedule(static)
124         for (int thread = 0; thread < nthreads; thread++)
125         {
126             int start = (nxf * thread) / nthreads;
127             int end   = std::min(nxf * (thread + 1) / nthreads, nxf);
128
129             calc_x_times_f(end - start, x + start, f + start, bScrewPBC, box,
130                            thread == 0 ? x_times_f : xf_buf[thread * 3]);
131         }
132
133         for (int thread = 1; thread < nthreads; thread++)
134         {
135             m_add(x_times_f, xf_buf[thread * 3], x_times_f);
136         }
137     }
138
139     for (int d = 0; d < DIM; d++)
140     {
141         upd_vir(vir[d], x_times_f[d][XX], x_times_f[d][YY], x_times_f[d][ZZ]);
142     }
143 }
144
145
146 static void
147 lo_fcv(int i0, int i1, const real x[], const real f[], tensor vir, const int is[], const real box[], gmx_bool bTriclinic)
148 {
149     int  i, i3, tx, ty, tz;
150     real xx, yy, zz;
151     real dvxx = 0, dvxy = 0, dvxz = 0, dvyx = 0, dvyy = 0, dvyz = 0, dvzx = 0, dvzy = 0, dvzz = 0;
152
153     if (bTriclinic)
154     {
155         for (i = i0; (i < i1); i++)
156         {
157             i3 = DIM * i;
158             tx = is[i3 + XX];
159             ty = is[i3 + YY];
160             tz = is[i3 + ZZ];
161
162             xx = x[i3 + XX] - tx * box[XXXX] - ty * box[YYXX] - tz * box[ZZXX];
163             dvxx += xx * f[i3 + XX];
164             dvxy += xx * f[i3 + YY];
165             dvxz += xx * f[i3 + ZZ];
166
167             yy = x[i3 + YY] - ty * box[YYYY] - tz * box[ZZYY];
168             dvyx += yy * f[i3 + XX];
169             dvyy += yy * f[i3 + YY];
170             dvyz += yy * f[i3 + ZZ];
171
172             zz = x[i3 + ZZ] - tz * box[ZZZZ];
173             dvzx += zz * f[i3 + XX];
174             dvzy += zz * f[i3 + YY];
175             dvzz += zz * f[i3 + ZZ];
176         }
177     }
178     else
179     {
180         for (i = i0; (i < i1); i++)
181         {
182             i3 = DIM * i;
183             tx = is[i3 + XX];
184             ty = is[i3 + YY];
185             tz = is[i3 + ZZ];
186
187             xx = x[i3 + XX] - tx * box[XXXX];
188             dvxx += xx * f[i3 + XX];
189             dvxy += xx * f[i3 + YY];
190             dvxz += xx * f[i3 + ZZ];
191
192             yy = x[i3 + YY] - ty * box[YYYY];
193             dvyx += yy * f[i3 + XX];
194             dvyy += yy * f[i3 + YY];
195             dvyz += yy * f[i3 + ZZ];
196
197             zz = x[i3 + ZZ] - tz * box[ZZZZ];
198             dvzx += zz * f[i3 + XX];
199             dvzy += zz * f[i3 + YY];
200             dvzz += zz * f[i3 + ZZ];
201         }
202     }
203
204     upd_vir(vir[XX], dvxx, dvxy, dvxz);
205     upd_vir(vir[YY], dvyx, dvyy, dvyz);
206     upd_vir(vir[ZZ], dvzx, dvzy, dvzz);
207 }
208
209 void f_calc_vir(int i0, int i1, const rvec x[], const rvec f[], tensor vir, const t_graph* g, const matrix box)
210 {
211     int start, end;
212
213     if (g && (g->nnodes > 0))
214     {
215         /* Calculate virial for bonded forces only when they belong to
216          * this node.
217          */
218         start = std::max(i0, g->at_start);
219         end   = std::min(i1, g->at_end);
220         lo_fcv(start, end, x[0], f[0], vir, g->ishift[0], box[0], TRICLINIC(box));
221
222         /* If not all atoms are bonded, calculate their virial contribution
223          * anyway, without shifting back their coordinates.
224          * Note the nifty pointer arithmetic...
225          */
226         if (start > i0)
227         {
228             calc_vir(start - i0, x + i0, f + i0, vir, FALSE, box);
229         }
230         if (end < i1)
231         {
232             calc_vir(i1 - end, x + end, f + end, vir, FALSE, box);
233         }
234     }
235     else
236     {
237         calc_vir(i1 - i0, x + i0, f + i0, vir, FALSE, box);
238     }
239 }