Apply re-formatting to C++ in src/ tree.
[alexxy/gromacs.git] / src / gromacs / domdec / domdec_vsite.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37
38 /*! \internal \file
39  *
40  * \brief This file implements functions for domdec to use
41  * while managing inter-atomic constraints.
42  *
43  * \author Berk Hess <hess@kth.se>
44  * \ingroup module_domdec
45  */
46
47 #include "gmxpre.h"
48
49 #include "domdec_vsite.h"
50
51 #include <cassert>
52
53 #include <algorithm>
54
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"
66
67 #include "domdec_specatomcomm.h"
68
69 void dd_move_f_vsites(const gmx_domdec_t& dd, gmx::ArrayRef<gmx::RVec> f, gmx::ArrayRef<gmx::RVec> fshift)
70 {
71     if (dd.vsite_comm)
72     {
73         dd_move_f_specat(&dd, dd.vsite_comm, as_rvec_array(f.data()), as_rvec_array(fshift.data()));
74     }
75 }
76
77 void dd_clear_f_vsites(const gmx_domdec_t& dd, gmx::ArrayRef<gmx::RVec> f)
78 {
79     if (dd.vsite_comm)
80     {
81         for (int i = dd.vsite_comm->at_start; i < dd.vsite_comm->at_end; i++)
82         {
83             clear_rvec(f[i]);
84         }
85     }
86 }
87
88 void dd_move_x_vsites(const gmx_domdec_t& dd, const matrix box, rvec* x)
89 {
90     if (dd.vsite_comm)
91     {
92         dd_move_x_specat(&dd, dd.vsite_comm, box, x, nullptr, FALSE);
93     }
94 }
95
96 void dd_clear_local_vsite_indices(gmx_domdec_t* dd)
97 {
98     if (dd->vsite_comm)
99     {
100         dd->ga2la_vsite->clear();
101     }
102 }
103
104 int dd_make_local_vsites(gmx_domdec_t* dd, int at_start, gmx::ArrayRef<InteractionList> lil)
105 {
106     std::vector<int>&    ireq         = dd->vsite_requestedGlobalAtomIndices;
107     gmx::HashedMap<int>* ga2la_specat = dd->ga2la_vsite;
108
109     ireq.clear();
110     /* Loop over all the home vsites */
111     for (int ftype = 0; ftype < F_NRE; ftype++)
112     {
113         if (interaction_function[ftype].flags & IF_VSITE)
114         {
115             const int              nral = NRAL(ftype);
116             const InteractionList& lilf = lil[ftype];
117             for (int i = 0; i < lilf.size(); i += 1 + nral)
118             {
119                 const int* iatoms = lilf.iatoms.data() + i;
120                 /* Check if we have the other atoms */
121                 for (int j = 1; j < 1 + nral; j++)
122                 {
123                     if (iatoms[j] < 0)
124                     {
125                         /* This is not a home atom,
126                          * we need to ask our neighbors.
127                          */
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))
131                         {
132                             /* Add this non-home atom to the list */
133                             ireq.push_back(a);
134                             /* Temporarily mark with -2,
135                              * we get the index later.
136                              */
137                             ga2la_specat->insert(a, -2);
138                         }
139                     }
140                 }
141             }
142         }
143     }
144
145     int at_end = setup_specat_communication(
146             dd, &ireq, dd->vsite_comm, ga2la_specat, at_start, 1, "vsite", "");
147
148     /* Fill in the missing indices */
149     for (int ftype = 0; ftype < F_NRE; ftype++)
150     {
151         if (interaction_function[ftype].flags & IF_VSITE)
152         {
153             const int        nral = NRAL(ftype);
154             InteractionList& lilf = lil[ftype];
155             for (int i = 0; i < lilf.size(); i += 1 + nral)
156             {
157                 t_iatom* iatoms = lilf.iatoms.data() + i;
158                 for (int j = 1; j < 1 + nral; j++)
159                 {
160                     if (iatoms[j] < 0)
161                     {
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");
164                         iatoms[j] = *a;
165                     }
166                 }
167             }
168         }
169     }
170
171     return at_end;
172 }
173
174 void init_domdec_vsites(gmx_domdec_t* dd, int n_intercg_vsite)
175 {
176     if (debug)
177     {
178         fprintf(debug, "Begin init_domdec_vsites\n");
179     }
180
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.
183      */
184     int numKeysEstimate = std::min(n_intercg_vsite / 20, n_intercg_vsite / (2 * dd->nnodes));
185     dd->ga2la_vsite     = new gmx::HashedMap<int>(numKeysEstimate);
186
187     dd->vsite_comm = new gmx_domdec_specat_comm_t;
188 }