2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2014,2015,2017,2018,2019, 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.
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
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, rvec av, rvec stddev, const MPI_Comm* mpiCommunicator)
75 for (const gmx::RVec& coord : x)
77 for (int d = 0; d < DIM; d++)
80 s2[d] += coord[d] * coord[d];
84 /* With mpiCommunicator != nullptr, x.size() is the home atom count */
85 int numAtoms = x.size();
89 constexpr int c_bufSize = 7;
90 double sendBuffer[c_bufSize];
91 double receiveBuffer[c_bufSize];
93 for (int d = 0; d < DIM; d++)
95 sendBuffer[d] = s1[d];
96 sendBuffer[DIM + d] = s2[d];
98 sendBuffer[6] = numAtoms;
100 MPI_Allreduce(sendBuffer, receiveBuffer, c_bufSize, MPI_DOUBLE, MPI_SUM, *mpiCommunicator);
102 for (int d = 0; d < DIM; d++)
104 s1[d] = receiveBuffer[d];
105 s2[d] = receiveBuffer[DIM + d];
107 numAtoms = gmx::roundToInt(receiveBuffer[6]);
110 GMX_UNUSED_VALUE(mpiCommunicator);
113 dsvmul(1.0 / numAtoms, s1, s1);
114 dsvmul(1.0 / numAtoms, s2, s2);
116 for (int d = 0; d < DIM; d++)
119 stddev[d] = std::sqrt(s2[d] - s1[d] * s1[d]);
123 /*! \brief Determines if dimensions require triclinic treatment and stores this info in ddbox */
124 static void set_tric_dir(const ivec* dd_nc, gmx_ddbox_t* ddbox, const matrix box)
126 int npbcdim, d, i, j;
128 real dep, inv_skew_fac2;
130 npbcdim = ddbox->npbcdim;
131 normal = ddbox->normal;
132 for (d = 0; d < DIM; d++)
134 ddbox->tric_dir[d] = 0;
135 for (j = d + 1; j < npbcdim; j++)
139 ddbox->tric_dir[d] = 1;
140 if (dd_nc != nullptr && (*dd_nc)[j] > 1 && (*dd_nc)[d] == 1)
143 "Domain decomposition has not been implemented for box vectors that "
144 "have non-zero components in directions that do not use domain "
145 "decomposition: ncells = %d %d %d, box vector[%d] = %f %f %f",
146 (*dd_nc)[XX], (*dd_nc)[YY], (*dd_nc)[ZZ], j + 1, box[j][XX],
147 box[j][YY], box[j][ZZ]);
152 /* Construct vectors v for dimension d that are perpendicular
153 * to the triclinic plane of dimension d. Each vector v[i] has
154 * v[i][i]=1 and v[i][d]!=0 for triclinic dimensions, while the third
155 * component is zero. These are used for computing the distance
156 * to a triclinic plane given the distance along dimension d.
157 * Set the trilinic skewing factor that translates
158 * the thickness of a slab perpendicular to this dimension
159 * into the real thickness of the slab.
161 if (ddbox->tric_dir[d])
165 if (d == XX || d == YY)
167 /* Normalize such that the "diagonal" is 1 */
168 svmul(1 / box[d + 1][d + 1], box[d + 1], v[d + 1]);
169 for (i = 0; i < d; i++)
173 inv_skew_fac2 += gmx::square(v[d + 1][d]);
176 /* Normalize such that the "diagonal" is 1 */
177 svmul(1 / box[d + 2][d + 2], box[d + 2], v[d + 2]);
178 /* Set v[d+2][d+1] to zero by shifting along v[d+1] */
179 dep = v[d + 2][d + 1] / v[d + 1][d + 1];
180 for (i = 0; i < DIM; i++)
182 v[d + 2][i] -= dep * v[d + 1][i];
184 inv_skew_fac2 += gmx::square(v[d + 2][d]);
186 cprod(v[d + 1], v[d + 2], normal[d]);
190 /* cross product with (1,0,0) */
192 normal[d][YY] = v[d + 1][ZZ];
193 normal[d][ZZ] = -v[d + 1][YY];
197 fprintf(debug, "box[%d] %.3f %.3f %.3f\n", d, box[d][XX], box[d][YY], box[d][ZZ]);
198 for (i = d + 1; i < DIM; i++)
200 fprintf(debug, " v[%d] %.3f %.3f %.3f\n", i, v[i][XX], v[i][YY], v[i][ZZ]);
204 ddbox->skew_fac[d] = 1.0 / std::sqrt(inv_skew_fac2);
205 /* Set the normal vector length to skew_fac */
206 dep = ddbox->skew_fac[d] / norm(normal[d]);
207 svmul(dep, normal[d], normal[d]);
211 fprintf(debug, "skew_fac[%d] = %f\n", d, ddbox->skew_fac[d]);
212 fprintf(debug, "normal[%d] %.3f %.3f %.3f\n", d, normal[d][XX], normal[d][YY],
218 ddbox->skew_fac[d] = 1;
220 for (i = 0; i < DIM; i++)
222 clear_rvec(ddbox->v[d][i]);
223 ddbox->v[d][i][i] = 1;
225 clear_rvec(normal[d]);
231 /*! \brief This function calculates bounding box and pbc info and populates ddbox */
232 static void low_set_ddbox(int numPbcDimensions,
233 int numBoundedDimensions,
236 bool calculateUnboundedSize,
237 gmx::ArrayRef<const gmx::RVec> x,
238 const MPI_Comm* mpiCommunicator,
245 ddbox->npbcdim = numPbcDimensions;
246 ddbox->nboundeddim = numBoundedDimensions;
248 for (d = 0; d < numBoundedDimensions; d++)
251 ddbox->box_size[d] = box[d][d];
254 if (ddbox->nboundeddim < DIM && calculateUnboundedSize)
256 calc_pos_av_stddev(x, av, stddev, mpiCommunicator);
258 /* GRID_STDDEV_FAC * stddev
259 * gives a uniform load for a rectangular block of cg's.
260 * For a sphere it is not a bad approximation for 4x1x1 up to 4x2x2.
262 for (d = ddbox->nboundeddim; d < DIM; d++)
264 b0 = av[d] - GRID_STDDEV_FAC * stddev[d];
265 b1 = av[d] + GRID_STDDEV_FAC * stddev[d];
268 fprintf(debug, "Setting global DD grid boundaries to %f - %f\n", b0, b1);
271 ddbox->box_size[d] = b1 - b0;
275 set_tric_dir(dd_nc, ddbox, box);
278 void set_ddbox(const gmx_domdec_t& dd,
279 bool masterRankHasTheSystemState,
281 bool calculateUnboundedSize,
282 gmx::ArrayRef<const gmx::RVec> x,
285 if (!masterRankHasTheSystemState || DDMASTER(dd))
287 bool needToReduceCoordinateData = (!masterRankHasTheSystemState && dd.nnodes > 1);
288 gmx::ArrayRef<const gmx::RVec> xRef = constArrayRefFromArray(
289 x.data(), masterRankHasTheSystemState ? x.size() : dd.comm->atomRanges.numHomeAtoms());
291 low_set_ddbox(dd.unitCellInfo.npbcdim, dd.unitCellInfo.numBoundedDimensions, &dd.numCells,
292 box, calculateUnboundedSize, xRef,
293 needToReduceCoordinateData ? &dd.mpi_comm_all : nullptr, ddbox);
296 if (masterRankHasTheSystemState)
298 dd_bcast(&dd, sizeof(gmx_ddbox_t), ddbox);
302 void set_ddbox_cr(const t_commrec& cr,
304 const t_inputrec& ir,
306 gmx::ArrayRef<const gmx::RVec> x,
311 low_set_ddbox(ePBC2npbcdim(ir.ePBC), inputrec2nboundeddim(&ir), dd_nc, box, true, x, nullptr, ddbox);
314 gmx_bcast(sizeof(gmx_ddbox_t), ddbox, &cr);