Split lines with many copyright years
[alexxy/gromacs.git] / src / gromacs / domdec / box.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009,2010,2014,2015,2017 by the GROMACS development team.
5  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 /*! \internal \file
38  *
39  * \brief This file defines functions used by the domdec module
40  * for (bounding) box and pbc information generation.
41  *
42  * \author Berk Hess <hess@kth.se>
43  * \ingroup module_domdec
44  */
45
46 #include "gmxpre.h"
47
48 #include "box.h"
49
50 #include "config.h"
51
52 #include "gromacs/domdec/domdec.h"
53 #include "gromacs/domdec/domdec_network.h"
54 #include "gromacs/domdec/domdec_struct.h"
55 #include "gromacs/gmxlib/network.h"
56 #include "gromacs/math/functions.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/mdlib/nsgrid.h"
59 #include "gromacs/mdtypes/commrec.h"
60 #include "gromacs/mdtypes/inputrec.h"
61 #include "gromacs/pbcutil/pbc.h"
62 #include "gromacs/topology/block.h"
63 #include "gromacs/utility/basedefinitions.h"
64 #include "gromacs/utility/fatalerror.h"
65
66 #include "domdec_internal.h"
67
68 /*! \brief Calculates the average and standard deviation in 3D of atoms */
69 static void calc_pos_av_stddev(gmx::ArrayRef<const gmx::RVec> x, rvec av, rvec stddev, const MPI_Comm* mpiCommunicator)
70 {
71     dvec s1, s2;
72
73     clear_dvec(s1);
74     clear_dvec(s2);
75
76     for (const gmx::RVec& coord : x)
77     {
78         for (int d = 0; d < DIM; d++)
79         {
80             s1[d] += coord[d];
81             s2[d] += coord[d] * coord[d];
82         }
83     }
84
85     /* With mpiCommunicator != nullptr, x.size() is the home atom count */
86     int numAtoms = x.size();
87 #if GMX_MPI
88     if (mpiCommunicator)
89     {
90         constexpr int c_bufSize = 7;
91         double        sendBuffer[c_bufSize];
92         double        receiveBuffer[c_bufSize];
93
94         for (int d = 0; d < DIM; d++)
95         {
96             sendBuffer[d]       = s1[d];
97             sendBuffer[DIM + d] = s2[d];
98         }
99         sendBuffer[6] = numAtoms;
100
101         MPI_Allreduce(sendBuffer, receiveBuffer, c_bufSize, MPI_DOUBLE, MPI_SUM, *mpiCommunicator);
102
103         for (int d = 0; d < DIM; d++)
104         {
105             s1[d] = receiveBuffer[d];
106             s2[d] = receiveBuffer[DIM + d];
107         }
108         numAtoms = gmx::roundToInt(receiveBuffer[6]);
109     }
110 #else  // GMX_MPI
111     GMX_UNUSED_VALUE(mpiCommunicator);
112 #endif // GMX_MPI
113
114     dsvmul(1.0 / numAtoms, s1, s1);
115     dsvmul(1.0 / numAtoms, s2, s2);
116
117     for (int d = 0; d < DIM; d++)
118     {
119         av[d]     = s1[d];
120         stddev[d] = std::sqrt(s2[d] - s1[d] * s1[d]);
121     }
122 }
123
124 /*! \brief Determines if dimensions require triclinic treatment and stores this info in ddbox */
125 static void set_tric_dir(const ivec* dd_nc, gmx_ddbox_t* ddbox, const matrix box)
126 {
127     int   npbcdim, d, i, j;
128     rvec *v, *normal;
129     real  dep, inv_skew_fac2;
130
131     npbcdim = ddbox->npbcdim;
132     normal  = ddbox->normal;
133     for (d = 0; d < DIM; d++)
134     {
135         ddbox->tric_dir[d] = 0;
136         for (j = d + 1; j < npbcdim; j++)
137         {
138             if (box[j][d] != 0)
139             {
140                 ddbox->tric_dir[d] = 1;
141                 if (dd_nc != nullptr && (*dd_nc)[j] > 1 && (*dd_nc)[d] == 1)
142                 {
143                     gmx_fatal(FARGS,
144                               "Domain decomposition has not been implemented for box vectors that "
145                               "have non-zero components in directions that do not use domain "
146                               "decomposition: ncells = %d %d %d, box vector[%d] = %f %f %f",
147                               (*dd_nc)[XX], (*dd_nc)[YY], (*dd_nc)[ZZ], j + 1, box[j][XX],
148                               box[j][YY], box[j][ZZ]);
149                 }
150             }
151         }
152
153         /* Construct vectors v for dimension d that are perpendicular
154          * to the triclinic plane of dimension d. Each vector v[i] has
155          * v[i][i]=1 and v[i][d]!=0 for triclinic dimensions, while the third
156          * component is zero. These are used for computing the distance
157          * to a triclinic plane given the distance along dimension d.
158          * Set the trilinic skewing factor that translates
159          * the thickness of a slab perpendicular to this dimension
160          * into the real thickness of the slab.
161          */
162         if (ddbox->tric_dir[d])
163         {
164             inv_skew_fac2 = 1;
165             v             = ddbox->v[d];
166             if (d == XX || d == YY)
167             {
168                 /* Normalize such that the "diagonal" is 1 */
169                 svmul(1 / box[d + 1][d + 1], box[d + 1], v[d + 1]);
170                 for (i = 0; i < d; i++)
171                 {
172                     v[d + 1][i] = 0;
173                 }
174                 inv_skew_fac2 += gmx::square(v[d + 1][d]);
175                 if (d == XX)
176                 {
177                     /* Normalize such that the "diagonal" is 1 */
178                     svmul(1 / box[d + 2][d + 2], box[d + 2], v[d + 2]);
179                     /* Set v[d+2][d+1] to zero by shifting along v[d+1] */
180                     dep = v[d + 2][d + 1] / v[d + 1][d + 1];
181                     for (i = 0; i < DIM; i++)
182                     {
183                         v[d + 2][i] -= dep * v[d + 1][i];
184                     }
185                     inv_skew_fac2 += gmx::square(v[d + 2][d]);
186
187                     cprod(v[d + 1], v[d + 2], normal[d]);
188                 }
189                 else
190                 {
191                     /* cross product with (1,0,0) */
192                     normal[d][XX] = 0;
193                     normal[d][YY] = v[d + 1][ZZ];
194                     normal[d][ZZ] = -v[d + 1][YY];
195                 }
196                 if (debug)
197                 {
198                     fprintf(debug, "box[%d]  %.3f %.3f %.3f\n", d, box[d][XX], box[d][YY], box[d][ZZ]);
199                     for (i = d + 1; i < DIM; i++)
200                     {
201                         fprintf(debug, "  v[%d]  %.3f %.3f %.3f\n", i, v[i][XX], v[i][YY], v[i][ZZ]);
202                     }
203                 }
204             }
205             ddbox->skew_fac[d] = 1.0 / std::sqrt(inv_skew_fac2);
206             /* Set the normal vector length to skew_fac */
207             dep = ddbox->skew_fac[d] / norm(normal[d]);
208             svmul(dep, normal[d], normal[d]);
209
210             if (debug)
211             {
212                 fprintf(debug, "skew_fac[%d] = %f\n", d, ddbox->skew_fac[d]);
213                 fprintf(debug, "normal[%d]  %.3f %.3f %.3f\n", d, normal[d][XX], normal[d][YY],
214                         normal[d][ZZ]);
215             }
216         }
217         else
218         {
219             ddbox->skew_fac[d] = 1;
220
221             for (i = 0; i < DIM; i++)
222             {
223                 clear_rvec(ddbox->v[d][i]);
224                 ddbox->v[d][i][i] = 1;
225             }
226             clear_rvec(normal[d]);
227             normal[d][d] = 1;
228         }
229     }
230 }
231
232 /*! \brief This function calculates bounding box and pbc info and populates ddbox */
233 static void low_set_ddbox(int                            numPbcDimensions,
234                           int                            numBoundedDimensions,
235                           const ivec*                    dd_nc,
236                           const matrix                   box,
237                           bool                           calculateUnboundedSize,
238                           gmx::ArrayRef<const gmx::RVec> x,
239                           const MPI_Comm*                mpiCommunicator,
240                           gmx_ddbox_t*                   ddbox)
241 {
242     rvec av, stddev;
243     real b0, b1;
244     int  d;
245
246     ddbox->npbcdim     = numPbcDimensions;
247     ddbox->nboundeddim = numBoundedDimensions;
248
249     for (d = 0; d < numBoundedDimensions; d++)
250     {
251         ddbox->box0[d]     = 0;
252         ddbox->box_size[d] = box[d][d];
253     }
254
255     if (ddbox->nboundeddim < DIM && calculateUnboundedSize)
256     {
257         calc_pos_av_stddev(x, av, stddev, mpiCommunicator);
258
259         /* GRID_STDDEV_FAC * stddev
260          * gives a uniform load for a rectangular block of cg's.
261          * For a sphere it is not a bad approximation for 4x1x1 up to 4x2x2.
262          */
263         for (d = ddbox->nboundeddim; d < DIM; d++)
264         {
265             b0 = av[d] - GRID_STDDEV_FAC * stddev[d];
266             b1 = av[d] + GRID_STDDEV_FAC * stddev[d];
267             if (debug)
268             {
269                 fprintf(debug, "Setting global DD grid boundaries to %f - %f\n", b0, b1);
270             }
271             ddbox->box0[d]     = b0;
272             ddbox->box_size[d] = b1 - b0;
273         }
274     }
275
276     set_tric_dir(dd_nc, ddbox, box);
277 }
278
279 void set_ddbox(const gmx_domdec_t&            dd,
280                bool                           masterRankHasTheSystemState,
281                const matrix                   box,
282                bool                           calculateUnboundedSize,
283                gmx::ArrayRef<const gmx::RVec> x,
284                gmx_ddbox_t*                   ddbox)
285 {
286     if (!masterRankHasTheSystemState || DDMASTER(dd))
287     {
288         bool needToReduceCoordinateData     = (!masterRankHasTheSystemState && dd.nnodes > 1);
289         gmx::ArrayRef<const gmx::RVec> xRef = constArrayRefFromArray(
290                 x.data(), masterRankHasTheSystemState ? x.size() : dd.comm->atomRanges.numHomeAtoms());
291
292         low_set_ddbox(dd.unitCellInfo.npbcdim, dd.unitCellInfo.numBoundedDimensions, &dd.numCells,
293                       box, calculateUnboundedSize, xRef,
294                       needToReduceCoordinateData ? &dd.mpi_comm_all : nullptr, ddbox);
295     }
296
297     if (masterRankHasTheSystemState)
298     {
299         dd_bcast(&dd, sizeof(gmx_ddbox_t), ddbox);
300     }
301 }
302
303 void set_ddbox_cr(const t_commrec&               cr,
304                   const ivec*                    dd_nc,
305                   const t_inputrec&              ir,
306                   const matrix                   box,
307                   gmx::ArrayRef<const gmx::RVec> x,
308                   gmx_ddbox_t*                   ddbox)
309 {
310     if (MASTER(&cr))
311     {
312         low_set_ddbox(ePBC2npbcdim(ir.ePBC), inputrec2nboundeddim(&ir), dd_nc, box, true, x, nullptr, ddbox);
313     }
314
315     gmx_bcast(sizeof(gmx_ddbox_t), ddbox, &cr);
316 }