2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006,2007,2008,2009,2010 by the GROMACS development team.
5 * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
6 * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
40 * \brief This file implements functions for domdec to use
41 * while managing inter-atomic constraints.
43 * \author Berk Hess <hess@kth.se>
44 * \ingroup module_domdec
49 #include "domdec_vsite.h"
55 #include "gromacs/domdec/domdec.h"
56 #include "gromacs/domdec/domdec_struct.h"
57 #include "gromacs/domdec/ga2la.h"
58 #include "gromacs/domdec/hashedmap.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/pbcutil/ishift.h"
62 #include "gromacs/topology/ifunc.h"
63 #include "gromacs/topology/mtop_util.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/gmxassert.h"
67 #include "domdec_specatomcomm.h"
69 void dd_move_f_vsites(gmx_domdec_t* dd, rvec* f, rvec* fshift)
73 dd_move_f_specat(dd, dd->vsite_comm, f, fshift);
77 void dd_clear_f_vsites(gmx_domdec_t* dd, rvec* f)
83 for (i = dd->vsite_comm->at_start; i < dd->vsite_comm->at_end; i++)
90 void dd_move_x_vsites(gmx_domdec_t* dd, const matrix box, rvec* x)
94 dd_move_x_specat(dd, dd->vsite_comm, box, x, nullptr, FALSE);
98 void dd_clear_local_vsite_indices(gmx_domdec_t* dd)
102 dd->ga2la_vsite->clear();
106 int dd_make_local_vsites(gmx_domdec_t* dd, int at_start, t_ilist* lil)
108 std::vector<int>& ireq = dd->vsite_requestedGlobalAtomIndices;
109 gmx::HashedMap<int>* ga2la_specat = dd->ga2la_vsite;
112 /* Loop over all the home vsites */
113 for (int ftype = 0; ftype < F_NRE; ftype++)
115 if (interaction_function[ftype].flags & IF_VSITE)
117 const int nral = NRAL(ftype);
118 const t_ilist& lilf = lil[ftype];
119 for (int i = 0; i < lilf.nr; i += 1 + nral)
121 const t_iatom* iatoms = lilf.iatoms + i;
122 /* Check if we have the other atoms */
123 for (int j = 1; j < 1 + nral; j++)
127 /* This is not a home atom,
128 * we need to ask our neighbors.
130 int a = -iatoms[j] - 1;
131 /* Check to not ask for the same atom more than once */
132 if (!dd->ga2la_vsite->find(a))
134 /* Add this non-home atom to the list */
136 /* Temporarily mark with -2,
137 * we get the index later.
139 ga2la_specat->insert(a, -2);
147 int at_end = setup_specat_communication(dd, &ireq, dd->vsite_comm, ga2la_specat, at_start, 1,
150 /* Fill in the missing indices */
151 for (int ftype = 0; ftype < F_NRE; ftype++)
153 if (interaction_function[ftype].flags & IF_VSITE)
155 const int nral = NRAL(ftype);
156 t_ilist& lilf = lil[ftype];
157 for (int i = 0; i < lilf.nr; i += 1 + nral)
159 t_iatom* iatoms = lilf.iatoms + i;
160 for (int j = 1; j < 1 + nral; j++)
164 const int* a = ga2la_specat->find(-iatoms[j] - 1);
165 GMX_ASSERT(a, "We have checked before that this atom index has been set");
176 void init_domdec_vsites(gmx_domdec_t* dd, int n_intercg_vsite)
180 fprintf(debug, "Begin init_domdec_vsites\n");
183 /* Use a hash table for the global to local index.
184 * The number of keys is a rough estimate, it will be optimized later.
186 int numKeysEstimate = std::min(n_intercg_vsite / 20, n_intercg_vsite / (2 * dd->nnodes));
187 dd->ga2la_vsite = new gmx::HashedMap<int>(numKeysEstimate);
189 dd->vsite_comm = new gmx_domdec_specat_comm_t;