2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006,2007,2008,2009,2010,2012,2013,2014,2015,2016,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 implements functions for domdec to use
39 * while managing inter-atomic constraints.
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_domdec
47 #include "domdec_vsite.h"
53 #include "gromacs/domdec/domdec.h"
54 #include "gromacs/domdec/domdec_struct.h"
55 #include "gromacs/domdec/ga2la.h"
56 #include "gromacs/domdec/hashedmap.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/mdtypes/commrec.h"
59 #include "gromacs/pbcutil/ishift.h"
60 #include "gromacs/topology/ifunc.h"
61 #include "gromacs/topology/mtop_util.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/gmxassert.h"
65 #include "domdec_specatomcomm.h"
67 void dd_move_f_vsites(gmx_domdec_t* dd, rvec* f, rvec* fshift)
71 dd_move_f_specat(dd, dd->vsite_comm, f, fshift);
75 void dd_clear_f_vsites(gmx_domdec_t* dd, rvec* f)
81 for (i = dd->vsite_comm->at_start; i < dd->vsite_comm->at_end; i++)
88 void dd_move_x_vsites(gmx_domdec_t* dd, const matrix box, rvec* x)
92 dd_move_x_specat(dd, dd->vsite_comm, box, x, nullptr, FALSE);
96 void dd_clear_local_vsite_indices(gmx_domdec_t* dd)
100 dd->ga2la_vsite->clear();
104 int dd_make_local_vsites(gmx_domdec_t* dd, int at_start, t_ilist* lil)
106 std::vector<int>& ireq = dd->vsite_requestedGlobalAtomIndices;
107 gmx::HashedMap<int>* ga2la_specat = dd->ga2la_vsite;
110 /* Loop over all the home vsites */
111 for (int ftype = 0; ftype < F_NRE; ftype++)
113 if (interaction_function[ftype].flags & IF_VSITE)
115 const int nral = NRAL(ftype);
116 const t_ilist& lilf = lil[ftype];
117 for (int i = 0; i < lilf.nr; i += 1 + nral)
119 const t_iatom* iatoms = lilf.iatoms + i;
120 /* Check if we have the other atoms */
121 for (int j = 1; j < 1 + nral; j++)
125 /* This is not a home atom,
126 * we need to ask our neighbors.
128 int a = -iatoms[j] - 1;
129 /* Check to not ask for the same atom more than once */
130 if (!dd->ga2la_vsite->find(a))
132 /* Add this non-home atom to the list */
134 /* Temporarily mark with -2,
135 * we get the index later.
137 ga2la_specat->insert(a, -2);
145 int at_end = setup_specat_communication(dd, &ireq, dd->vsite_comm, ga2la_specat, at_start, 1,
148 /* Fill in the missing indices */
149 for (int ftype = 0; ftype < F_NRE; ftype++)
151 if (interaction_function[ftype].flags & IF_VSITE)
153 const int nral = NRAL(ftype);
154 t_ilist& lilf = lil[ftype];
155 for (int i = 0; i < lilf.nr; i += 1 + nral)
157 t_iatom* iatoms = lilf.iatoms + i;
158 for (int j = 1; j < 1 + nral; j++)
162 const int* a = ga2la_specat->find(-iatoms[j] - 1);
163 GMX_ASSERT(a, "We have checked before that this atom index has been set");
174 void init_domdec_vsites(gmx_domdec_t* dd, int n_intercg_vsite)
178 fprintf(debug, "Begin init_domdec_vsites\n");
181 /* Use a hash table for the global to local index.
182 * The number of keys is a rough estimate, it will be optimized later.
184 int numKeysEstimate = std::min(n_intercg_vsite / 20, n_intercg_vsite / (2 * dd->nnodes));
185 dd->ga2la_vsite = new gmx::HashedMap<int>(numKeysEstimate);
187 dd->vsite_comm = new gmx_domdec_specat_comm_t;