Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / domdec / box.cpp
blob050ebdccb655237183151cf5b16352dd7492f404
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2014,2015,2017,2018, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 /*! \internal \file
38 * \brief This file defines functions used by the domdec module
39 * for (bounding) box and pbc information generation.
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_domdec
45 #include "gmxpre.h"
47 #include "box.h"
49 #include "config.h"
51 #include "gromacs/domdec/domdec.h"
52 #include "gromacs/domdec/domdec_network.h"
53 #include "gromacs/domdec/domdec_struct.h"
54 #include "gromacs/gmxlib/network.h"
55 #include "gromacs/math/functions.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/mdlib/nsgrid.h"
58 #include "gromacs/mdtypes/commrec.h"
59 #include "gromacs/mdtypes/inputrec.h"
60 #include "gromacs/pbcutil/pbc.h"
61 #include "gromacs/topology/block.h"
62 #include "gromacs/utility/basedefinitions.h"
63 #include "gromacs/utility/fatalerror.h"
65 #include "domdec_internal.h"
67 /*! \brief Calculates the average and standard deviation in 3D of atoms */
68 static void calc_pos_av_stddev(gmx::ArrayRef<const gmx::RVec> x,
69 rvec av,
70 rvec stddev,
71 const MPI_Comm *mpiCommunicator)
73 dvec s1, s2;
75 clear_dvec(s1);
76 clear_dvec(s2);
78 for (const gmx::RVec &coord : x)
80 for (int d = 0; d < DIM; d++)
82 s1[d] += coord[d];
83 s2[d] += coord[d]*coord[d];
87 /* With mpiCommunicator != nullptr, x.size() is the home atom count */
88 int numAtoms = x.size();
89 #if GMX_MPI
90 if (mpiCommunicator)
92 constexpr int c_bufSize = 7;
93 double sendBuffer[c_bufSize];
94 double receiveBuffer[c_bufSize];
96 for (int d = 0; d < DIM; d++)
98 sendBuffer[d] = s1[d];
99 sendBuffer[DIM + d] = s2[d];
101 sendBuffer[6] = numAtoms;
103 MPI_Allreduce(sendBuffer, receiveBuffer, c_bufSize, MPI_DOUBLE,
104 MPI_SUM, *mpiCommunicator);
106 for (int d = 0; d < DIM; d++)
108 s1[d] = receiveBuffer[d];
109 s2[d] = receiveBuffer[DIM + d];
111 numAtoms = gmx::roundToInt(receiveBuffer[6]);
113 #else // GMX_MPI
114 GMX_UNUSED_VALUE(mpiCommunicator);
115 #endif // GMX_MPI
117 dsvmul(1.0/numAtoms, s1, s1);
118 dsvmul(1.0/numAtoms, s2, s2);
120 for (int d = 0; d < DIM; d++)
122 av[d] = s1[d];
123 stddev[d] = std::sqrt(s2[d] - s1[d]*s1[d]);
127 /*! \brief Determines if dimensions require triclinic treatment and stores this info in ddbox */
128 static void set_tric_dir(const ivec *dd_nc, gmx_ddbox_t *ddbox, const matrix box)
130 int npbcdim, d, i, j;
131 rvec *v, *normal;
132 real dep, inv_skew_fac2;
134 npbcdim = ddbox->npbcdim;
135 normal = ddbox->normal;
136 for (d = 0; d < DIM; d++)
138 ddbox->tric_dir[d] = 0;
139 for (j = d+1; j < npbcdim; j++)
141 if (box[j][d] != 0)
143 ddbox->tric_dir[d] = 1;
144 if (dd_nc != nullptr && (*dd_nc)[j] > 1 && (*dd_nc)[d] == 1)
146 gmx_fatal(FARGS, "Domain decomposition has not been implemented for box vectors that have non-zero components in directions that do not use domain decomposition: ncells = %d %d %d, box vector[%d] = %f %f %f",
147 (*dd_nc)[XX], (*dd_nc)[YY], (*dd_nc)[ZZ],
148 j+1, box[j][XX], box[j][YY], box[j][ZZ]);
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.
162 if (ddbox->tric_dir[d])
164 inv_skew_fac2 = 1;
165 v = ddbox->v[d];
166 if (d == XX || d == YY)
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++)
172 v[d+1][i] = 0;
174 inv_skew_fac2 += gmx::square(v[d+1][d]);
175 if (d == XX)
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++)
183 v[d+2][i] -= dep*v[d+1][i];
185 inv_skew_fac2 += gmx::square(v[d+2][d]);
187 cprod(v[d+1], v[d+2], normal[d]);
189 else
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];
196 if (debug)
198 fprintf(debug, "box[%d] %.3f %.3f %.3f\n",
199 d, box[d][XX], box[d][YY], box[d][ZZ]);
200 for (i = d+1; i < DIM; i++)
202 fprintf(debug, " v[%d] %.3f %.3f %.3f\n",
203 i, v[i][XX], v[i][YY], v[i][ZZ]);
207 ddbox->skew_fac[d] = 1.0/std::sqrt(inv_skew_fac2);
208 /* Set the normal vector length to skew_fac */
209 dep = ddbox->skew_fac[d]/norm(normal[d]);
210 svmul(dep, normal[d], normal[d]);
212 if (debug)
214 fprintf(debug, "skew_fac[%d] = %f\n", d, ddbox->skew_fac[d]);
215 fprintf(debug, "normal[%d] %.3f %.3f %.3f\n",
216 d, normal[d][XX], normal[d][YY], normal[d][ZZ]);
219 else
221 ddbox->skew_fac[d] = 1;
223 for (i = 0; i < DIM; i++)
225 clear_rvec(ddbox->v[d][i]);
226 ddbox->v[d][i][i] = 1;
228 clear_rvec(normal[d]);
229 normal[d][d] = 1;
234 /*! \brief This function calculates bounding box and pbc info and populates ddbox */
235 static void low_set_ddbox(int numPbcDimensions,
236 int numBoundedDimensions,
237 const ivec *dd_nc,
238 const matrix box,
239 bool calculateUnboundedSize,
240 gmx::ArrayRef<const gmx::RVec> x,
241 const MPI_Comm *mpiCommunicator,
242 gmx_ddbox_t *ddbox)
244 rvec av, stddev;
245 real b0, b1;
246 int d;
248 ddbox->npbcdim = numPbcDimensions;
249 ddbox->nboundeddim = numBoundedDimensions;
251 for (d = 0; d < numBoundedDimensions; d++)
253 ddbox->box0[d] = 0;
254 ddbox->box_size[d] = box[d][d];
257 if (ddbox->nboundeddim < DIM && calculateUnboundedSize)
259 calc_pos_av_stddev(x, av, stddev, mpiCommunicator);
261 /* GRID_STDDEV_FAC * stddev
262 * gives a uniform load for a rectangular block of cg's.
263 * For a sphere it is not a bad approximation for 4x1x1 up to 4x2x2.
265 for (d = ddbox->nboundeddim; d < DIM; d++)
267 b0 = av[d] - GRID_STDDEV_FAC*stddev[d];
268 b1 = av[d] + GRID_STDDEV_FAC*stddev[d];
269 if (debug)
271 fprintf(debug, "Setting global DD grid boundaries to %f - %f\n",
272 b0, b1);
274 ddbox->box0[d] = b0;
275 ddbox->box_size[d] = b1 - b0;
279 set_tric_dir(dd_nc, ddbox, box);
282 void set_ddbox(const gmx_domdec_t &dd,
283 bool masterRankHasTheSystemState,
284 const matrix box,
285 bool calculateUnboundedSize,
286 gmx::ArrayRef<const gmx::RVec> x,
287 gmx_ddbox_t *ddbox)
289 if (!masterRankHasTheSystemState || DDMASTER(dd))
291 bool needToReduceCoordinateData =
292 (!masterRankHasTheSystemState && dd.nnodes > 1);
293 gmx::ArrayRef<const gmx::RVec> xRef =
294 constArrayRefFromArray(x.data(), masterRankHasTheSystemState ? x.size() : dd.comm->atomRanges.numHomeAtoms());
296 low_set_ddbox(dd.npbcdim, dd.numBoundedDimensions,
297 &dd.nc, box, calculateUnboundedSize, xRef,
298 needToReduceCoordinateData ? &dd.mpi_comm_all : nullptr,
299 ddbox);
302 if (masterRankHasTheSystemState)
304 dd_bcast(&dd, sizeof(gmx_ddbox_t), ddbox);
308 void set_ddbox_cr(const t_commrec &cr,
309 const ivec *dd_nc,
310 const t_inputrec &ir,
311 const matrix box,
312 gmx::ArrayRef<const gmx::RVec> x,
313 gmx_ddbox_t *ddbox)
315 if (MASTER(&cr))
317 low_set_ddbox(ePBC2npbcdim(ir.ePBC), inputrec2nboundeddim(&ir),
318 dd_nc, box, true, x, nullptr, ddbox);
321 gmx_bcast(sizeof(gmx_ddbox_t), ddbox, &cr);
324 bool dynamic_dd_box(const gmx_domdec_t &dd)
326 return (dd.numBoundedDimensions < DIM || dd.haveDynamicBox);