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.
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
49 #include "gromacs/domdec/domdec.h"
50 #include "gromacs/domdec/domdec_network.h"
51 #include "gromacs/domdec/domdec_struct.h"
52 #include "gromacs/gmxlib/network.h"
53 #include "gromacs/math/functions.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/mdlib/nsgrid.h"
56 #include "gromacs/mdtypes/commrec.h"
57 #include "gromacs/mdtypes/inputrec.h"
58 #include "gromacs/pbcutil/pbc.h"
59 #include "gromacs/topology/block.h"
60 #include "gromacs/utility/basedefinitions.h"
61 #include "gromacs/utility/fatalerror.h"
63 #include "domdec_internal.h"
65 /*! \brief Calculates the average and standard deviation in 3D of n atoms */
66 static void calc_pos_av_stddev(int n, const rvec *x,
68 const MPI_Comm *mpiCommunicator)
75 for (int i = 0; i < n; i++)
77 for (int d = 0; d < DIM; d++)
80 s2[d] += x[i][d]*x[i][d];
87 constexpr int c_bufSize = 7;
88 double sendBuffer[c_bufSize];
89 double receiveBuffer[c_bufSize];
91 for (int d = 0; d < DIM; d++)
93 sendBuffer[d] = s1[d];
94 sendBuffer[DIM + d] = s2[d];
98 MPI_Allreduce(sendBuffer, receiveBuffer, c_bufSize, MPI_DOUBLE,
99 MPI_SUM, *mpiCommunicator);
101 for (int d = 0; d < DIM; d++)
103 s1[d] = receiveBuffer[d];
104 s2[d] = receiveBuffer[DIM + d];
106 n = (int)(receiveBuffer[6] + 0.5);
109 GMX_UNUSED_VALUE(mpiCommunicator);
112 dsvmul(1.0/n, s1, s1);
113 dsvmul(1.0/n, s2, s2);
115 for (int d = 0; d < DIM; d++)
118 stddev[d] = std::sqrt(s2[d] - s1[d]*s1[d]);
122 /*! \brief Determines if dimensions require triclinic treatment and stores this info in ddbox */
123 static void set_tric_dir(const ivec *dd_nc, gmx_ddbox_t *ddbox, const matrix box)
125 int npbcdim, d, i, j;
127 real dep, inv_skew_fac2;
129 npbcdim = ddbox->npbcdim;
130 normal = ddbox->normal;
131 for (d = 0; d < DIM; d++)
133 ddbox->tric_dir[d] = 0;
134 for (j = d+1; j < npbcdim; j++)
138 ddbox->tric_dir[d] = 1;
139 if (dd_nc != nullptr && (*dd_nc)[j] > 1 && (*dd_nc)[d] == 1)
141 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",
142 (*dd_nc)[XX], (*dd_nc)[YY], (*dd_nc)[ZZ],
143 j+1, box[j][XX], box[j][YY], box[j][ZZ]);
148 /* Construct vectors v for dimension d that are perpendicular
149 * to the triclinic plane of dimension d. Each vector v[i] has
150 * v[i][i]=1 and v[i][d]!=0 for triclinic dimensions, while the third
151 * component is zero. These are used for computing the distance
152 * to a triclinic plane given the distance along dimension d.
153 * Set the trilinic skewing factor that translates
154 * the thickness of a slab perpendicular to this dimension
155 * into the real thickness of the slab.
157 if (ddbox->tric_dir[d])
161 if (d == XX || d == YY)
163 /* Normalize such that the "diagonal" is 1 */
164 svmul(1/box[d+1][d+1], box[d+1], v[d+1]);
165 for (i = 0; i < d; i++)
169 inv_skew_fac2 += gmx::square(v[d+1][d]);
172 /* Normalize such that the "diagonal" is 1 */
173 svmul(1/box[d+2][d+2], box[d+2], v[d+2]);
174 /* Set v[d+2][d+1] to zero by shifting along v[d+1] */
175 dep = v[d+2][d+1]/v[d+1][d+1];
176 for (i = 0; i < DIM; i++)
178 v[d+2][i] -= dep*v[d+1][i];
180 inv_skew_fac2 += gmx::square(v[d+2][d]);
182 cprod(v[d+1], v[d+2], normal[d]);
186 /* cross product with (1,0,0) */
188 normal[d][YY] = v[d+1][ZZ];
189 normal[d][ZZ] = -v[d+1][YY];
193 fprintf(debug, "box[%d] %.3f %.3f %.3f\n",
194 d, box[d][XX], box[d][YY], box[d][ZZ]);
195 for (i = d+1; i < DIM; i++)
197 fprintf(debug, " v[%d] %.3f %.3f %.3f\n",
198 i, v[i][XX], v[i][YY], v[i][ZZ]);
202 ddbox->skew_fac[d] = 1.0/std::sqrt(inv_skew_fac2);
203 /* Set the normal vector length to skew_fac */
204 dep = ddbox->skew_fac[d]/norm(normal[d]);
205 svmul(dep, normal[d], normal[d]);
209 fprintf(debug, "skew_fac[%d] = %f\n", d, ddbox->skew_fac[d]);
210 fprintf(debug, "normal[%d] %.3f %.3f %.3f\n",
211 d, normal[d][XX], normal[d][YY], normal[d][ZZ]);
216 ddbox->skew_fac[d] = 1;
218 for (i = 0; i < DIM; i++)
220 clear_rvec(ddbox->v[d][i]);
221 ddbox->v[d][i][i] = 1;
223 clear_rvec(normal[d]);
229 /*! \brief This function calculates bounding box and pbc info and populates ddbox */
230 static void low_set_ddbox(const t_inputrec *ir, const ivec *dd_nc, const matrix box,
231 bool calculateUnboundedSize,
232 int numAtoms, const rvec *x,
233 const MPI_Comm *mpiCommunicator,
240 ddbox->npbcdim = ePBC2npbcdim(ir->ePBC);
241 ddbox->nboundeddim = inputrec2nboundeddim(ir);
243 for (d = 0; d < ddbox->nboundeddim; d++)
246 ddbox->box_size[d] = box[d][d];
249 if (ddbox->nboundeddim < DIM && calculateUnboundedSize)
251 calc_pos_av_stddev(numAtoms, x, av, stddev, mpiCommunicator);
253 /* GRID_STDDEV_FAC * stddev
254 * gives a uniform load for a rectangular block of cg's.
255 * For a sphere it is not a bad approximation for 4x1x1 up to 4x2x2.
257 for (d = ddbox->nboundeddim; d < DIM; d++)
259 b0 = av[d] - GRID_STDDEV_FAC*stddev[d];
260 b1 = av[d] + GRID_STDDEV_FAC*stddev[d];
263 fprintf(debug, "Setting global DD grid boundaries to %f - %f\n",
267 ddbox->box_size[d] = b1 - b0;
271 set_tric_dir(dd_nc, ddbox, box);
274 void set_ddbox(gmx_domdec_t *dd, bool masterRankHasTheSystemState,
275 const t_inputrec *ir, const matrix box,
276 bool calculateUnboundedSize,
277 gmx::ArrayRef<const gmx::RVec> x,
280 if (!masterRankHasTheSystemState || DDMASTER(dd))
282 bool needToReduceCoordinateData = (!masterRankHasTheSystemState && dd->nnodes > 1);
284 low_set_ddbox(ir, &dd->nc, box, calculateUnboundedSize,
285 masterRankHasTheSystemState ? x.size() : dd->comm->atomRanges.numHomeAtoms(), as_rvec_array(x.data()),
286 needToReduceCoordinateData ? &dd->mpi_comm_all : nullptr,
290 if (masterRankHasTheSystemState)
292 dd_bcast(dd, sizeof(gmx_ddbox_t), ddbox);
296 void set_ddbox_cr(const t_commrec *cr, const ivec *dd_nc,
297 const t_inputrec *ir, const matrix box,
298 gmx::ArrayRef<const gmx::RVec> x,
303 low_set_ddbox(ir, dd_nc, box, true, x.size(), as_rvec_array(x.data()), nullptr, ddbox);
306 gmx_bcast(sizeof(gmx_ddbox_t), ddbox, cr);