2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2014,2015,2017, 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
47 #include "gromacs/domdec/domdec.h"
48 #include "gromacs/domdec/domdec_network.h"
49 #include "gromacs/domdec/domdec_struct.h"
50 #include "gromacs/gmxlib/network.h"
51 #include "gromacs/math/functions.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/mdlib/nsgrid.h"
54 #include "gromacs/mdtypes/commrec.h"
55 #include "gromacs/mdtypes/inputrec.h"
56 #include "gromacs/pbcutil/pbc.h"
57 #include "gromacs/utility/fatalerror.h"
59 /*! \brief Calculates the average and standard deviation in 3D of n charge groups */
60 static void calc_cgcm_av_stddev(const t_block
*cgs
, int n
, const rvec
*x
,
67 int cg
, d
, k0
, k1
, k
, nrcg
;
75 for (cg
= 0; cg
< n
; cg
++)
82 copy_rvec(x
[k0
], cg_cm
);
89 for (k
= k0
; (k
< k1
); k
++)
91 rvec_inc(cg_cm
, x
[k
]);
93 for (d
= 0; (d
< DIM
); d
++)
98 for (d
= 0; d
< DIM
; d
++)
101 s2
[d
] += cg_cm
[d
]*cg_cm
[d
];
105 if (cr_sum
!= nullptr)
107 for (d
= 0; d
< DIM
; d
++)
113 gmx_sumd(7, buf
, cr_sum
);
114 for (d
= 0; d
< DIM
; d
++)
119 n
= (int)(buf
[6] + 0.5);
122 dsvmul(1.0/n
, s1
, s1
);
123 dsvmul(1.0/n
, s2
, s2
);
125 for (d
= 0; d
< DIM
; d
++)
128 stddev
[d
] = std::sqrt(s2
[d
] - s1
[d
]*s1
[d
]);
132 /*! \brief Determines if dimensions require triclinic treatment and stores this info in ddbox */
133 static void set_tric_dir(const ivec
*dd_nc
, gmx_ddbox_t
*ddbox
, const matrix box
)
135 int npbcdim
, d
, i
, j
;
137 real dep
, inv_skew_fac2
;
139 npbcdim
= ddbox
->npbcdim
;
140 normal
= ddbox
->normal
;
141 for (d
= 0; d
< DIM
; d
++)
143 ddbox
->tric_dir
[d
] = 0;
144 for (j
= d
+1; j
< npbcdim
; j
++)
148 ddbox
->tric_dir
[d
] = 1;
149 if (dd_nc
!= nullptr && (*dd_nc
)[j
] > 1 && (*dd_nc
)[d
] == 1)
151 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",
152 (*dd_nc
)[XX
], (*dd_nc
)[YY
], (*dd_nc
)[ZZ
],
153 j
+1, box
[j
][XX
], box
[j
][YY
], box
[j
][ZZ
]);
158 /* Convert box vectors to orthogonal vectors for this dimension,
159 * for use in distance calculations.
160 * Set the trilinic skewing factor that translates
161 * the thickness of a slab perpendicular to this dimension
162 * into the real thickness of the slab.
164 if (ddbox
->tric_dir
[d
])
168 if (d
== XX
|| d
== YY
)
170 /* Normalize such that the "diagonal" is 1 */
171 svmul(1/box
[d
+1][d
+1], box
[d
+1], v
[d
+1]);
172 for (i
= 0; i
< d
; i
++)
176 inv_skew_fac2
+= gmx::square(v
[d
+1][d
]);
179 /* Normalize such that the "diagonal" is 1 */
180 svmul(1/box
[d
+2][d
+2], box
[d
+2], v
[d
+2]);
181 for (i
= 0; i
< d
; i
++)
185 /* Make vector [d+2] perpendicular to vector [d+1],
186 * this does not affect the normalization.
188 dep
= iprod(v
[d
+1], v
[d
+2])/norm2(v
[d
+1]);
189 for (i
= 0; i
< DIM
; i
++)
191 v
[d
+2][i
] -= dep
*v
[d
+1][i
];
193 inv_skew_fac2
+= gmx::square(v
[d
+2][d
]);
195 cprod(v
[d
+1], v
[d
+2], normal
[d
]);
199 /* cross product with (1,0,0) */
201 normal
[d
][YY
] = v
[d
+1][ZZ
];
202 normal
[d
][ZZ
] = -v
[d
+1][YY
];
206 fprintf(debug
, "box[%d] %.3f %.3f %.3f\n",
207 d
, box
[d
][XX
], box
[d
][YY
], box
[d
][ZZ
]);
208 for (i
= d
+1; i
< DIM
; i
++)
210 fprintf(debug
, " v[%d] %.3f %.3f %.3f\n",
211 i
, v
[i
][XX
], v
[i
][YY
], v
[i
][ZZ
]);
215 ddbox
->skew_fac
[d
] = 1.0/std::sqrt(inv_skew_fac2
);
216 /* Set the normal vector length to skew_fac */
217 dep
= ddbox
->skew_fac
[d
]/norm(normal
[d
]);
218 svmul(dep
, normal
[d
], normal
[d
]);
222 fprintf(debug
, "skew_fac[%d] = %f\n", d
, ddbox
->skew_fac
[d
]);
223 fprintf(debug
, "normal[%d] %.3f %.3f %.3f\n",
224 d
, normal
[d
][XX
], normal
[d
][YY
], normal
[d
][ZZ
]);
229 ddbox
->skew_fac
[d
] = 1;
231 for (i
= 0; i
< DIM
; i
++)
233 clear_rvec(ddbox
->v
[d
][i
]);
234 ddbox
->v
[d
][i
][i
] = 1;
236 clear_rvec(normal
[d
]);
242 /*! \brief This function calculates bounding box and pbc info and populates ddbox */
243 static void low_set_ddbox(const t_inputrec
*ir
, const ivec
*dd_nc
, const matrix box
,
244 gmx_bool bCalcUnboundedSize
, int ncg
, const t_block
*cgs
, const rvec
*x
,
252 ddbox
->npbcdim
= ePBC2npbcdim(ir
->ePBC
);
253 ddbox
->nboundeddim
= inputrec2nboundeddim(ir
);
255 for (d
= 0; d
< ddbox
->nboundeddim
; d
++)
258 ddbox
->box_size
[d
] = box
[d
][d
];
261 if (ddbox
->nboundeddim
< DIM
&& bCalcUnboundedSize
)
263 calc_cgcm_av_stddev(cgs
, ncg
, x
, av
, stddev
, cr_sum
);
265 /* GRID_STDDEV_FAC * stddev
266 * gives a uniform load for a rectangular block of cg's.
267 * For a sphere it is not a bad approximation for 4x1x1 up to 4x2x2.
269 for (d
= ddbox
->nboundeddim
; d
< DIM
; d
++)
271 b0
= av
[d
] - GRID_STDDEV_FAC
*stddev
[d
];
272 b1
= av
[d
] + GRID_STDDEV_FAC
*stddev
[d
];
275 fprintf(debug
, "Setting global DD grid boundaries to %f - %f\n",
279 ddbox
->box_size
[d
] = b1
- b0
;
283 set_tric_dir(dd_nc
, ddbox
, box
);
286 void set_ddbox(gmx_domdec_t
*dd
, gmx_bool bMasterState
, t_commrec
*cr_sum
,
287 const t_inputrec
*ir
, const matrix box
,
288 gmx_bool bCalcUnboundedSize
, const t_block
*cgs
, const rvec
*x
,
291 if (!bMasterState
|| DDMASTER(dd
))
293 low_set_ddbox(ir
, &dd
->nc
, box
, bCalcUnboundedSize
,
294 bMasterState
? cgs
->nr
: dd
->ncg_home
, cgs
, x
,
295 bMasterState
? nullptr : cr_sum
,
301 dd_bcast(dd
, sizeof(gmx_ddbox_t
), ddbox
);
305 void set_ddbox_cr(t_commrec
*cr
, const ivec
*dd_nc
,
306 const t_inputrec
*ir
, const matrix box
,
307 const t_block
*cgs
, const rvec
*x
,
312 low_set_ddbox(ir
, dd_nc
, box
, TRUE
, cgs
->nr
, cgs
, x
, nullptr, ddbox
);
315 gmx_bcast(sizeof(gmx_ddbox_t
), ddbox
, cr
);