2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2014, 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.
43 #include "types/commrec.h"
45 #include "domdec_network.h"
49 #include "gmx_fatal.h"
51 static void calc_cgcm_av_stddev(t_block *cgs, int n, rvec *x, rvec av, rvec stddev,
57 int cg, d, k0, k1, k, nrcg;
65 for (cg = 0; cg < n; cg++)
72 copy_rvec(x[k0], cg_cm);
79 for (k = k0; (k < k1); k++)
81 rvec_inc(cg_cm, x[k]);
83 for (d = 0; (d < DIM); d++)
88 for (d = 0; d < DIM; d++)
91 s2[d] += cg_cm[d]*cg_cm[d];
97 for (d = 0; d < DIM; d++)
103 gmx_sumd(7, buf, cr_sum);
104 for (d = 0; d < DIM; d++)
109 n = (int)(buf[6] + 0.5);
112 dsvmul(1.0/n, s1, s1);
113 dsvmul(1.0/n, s2, s2);
115 for (d = 0; d < DIM; d++)
118 stddev[d] = sqrt(s2[d] - s1[d]*s1[d]);
122 static void set_tric_dir(ivec *dd_nc, gmx_ddbox_t *ddbox, matrix box)
124 int npbcdim, d, i, j;
126 real dep, inv_skew_fac2;
128 npbcdim = ddbox->npbcdim;
129 normal = ddbox->normal;
130 for (d = 0; d < DIM; d++)
132 ddbox->tric_dir[d] = 0;
133 for (j = d+1; j < npbcdim; j++)
137 ddbox->tric_dir[d] = 1;
138 if (dd_nc != NULL && (*dd_nc)[j] > 1 && (*dd_nc)[d] == 1)
140 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",
141 dd_nc[XX], dd_nc[YY], dd_nc[ZZ],
142 j+1, box[j][XX], box[j][YY], box[j][ZZ]);
147 /* Convert box vectors to orthogonal vectors for this dimension,
148 * for use in distance calculations.
149 * Set the trilinic skewing factor that translates
150 * the thickness of a slab perpendicular to this dimension
151 * into the real thickness of the slab.
153 if (ddbox->tric_dir[d])
157 if (d == XX || d == YY)
159 /* Normalize such that the "diagonal" is 1 */
160 svmul(1/box[d+1][d+1], box[d+1], v[d+1]);
161 for (i = 0; i < d; i++)
165 inv_skew_fac2 += sqr(v[d+1][d]);
168 /* Normalize such that the "diagonal" is 1 */
169 svmul(1/box[d+2][d+2], box[d+2], v[d+2]);
170 for (i = 0; i < d; i++)
174 /* Make vector [d+2] perpendicular to vector [d+1],
175 * this does not affect the normalization.
177 dep = iprod(v[d+1], v[d+2])/norm2(v[d+1]);
178 for (i = 0; i < DIM; i++)
180 v[d+2][i] -= dep*v[d+1][i];
182 inv_skew_fac2 += sqr(v[d+2][d]);
184 cprod(v[d+1], v[d+2], normal[d]);
188 /* cross product with (1,0,0) */
190 normal[d][YY] = v[d+1][ZZ];
191 normal[d][ZZ] = -v[d+1][YY];
195 fprintf(debug, "box[%d] %.3f %.3f %.3f\n",
196 d, box[d][XX], box[d][YY], box[d][ZZ]);
197 for (i = d+1; i < DIM; i++)
199 fprintf(debug, " v[%d] %.3f %.3f %.3f\n",
200 i, v[i][XX], v[i][YY], v[i][ZZ]);
204 ddbox->skew_fac[d] = 1.0/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",
213 d, normal[d][XX], normal[d][YY], normal[d][ZZ]);
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 static void low_set_ddbox(t_inputrec *ir, ivec *dd_nc, matrix box,
232 gmx_bool bCalcUnboundedSize, int ncg, t_block *cgs, rvec *x,
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 && bCalcUnboundedSize)
251 calc_cgcm_av_stddev(cgs, ncg, x, av, stddev, cr_sum);
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, gmx_bool bMasterState, t_commrec *cr_sum,
275 t_inputrec *ir, matrix box,
276 gmx_bool bCalcUnboundedSize, t_block *cgs, rvec *x,
279 if (!bMasterState || DDMASTER(dd))
281 low_set_ddbox(ir, &dd->nc, box, bCalcUnboundedSize,
282 bMasterState ? cgs->nr : dd->ncg_home, cgs, x,
283 bMasterState ? NULL : cr_sum,
289 dd_bcast(dd, sizeof(gmx_ddbox_t), ddbox);
293 void set_ddbox_cr(t_commrec *cr, ivec *dd_nc,
294 t_inputrec *ir, matrix box, t_block *cgs, rvec *x,
299 low_set_ddbox(ir, dd_nc, box, TRUE, cgs->nr, cgs, x, NULL, ddbox);
302 gmx_bcast(sizeof(gmx_ddbox_t), ddbox, cr);