30ca66aa2dbb1ab431c9e34c180c7bcfc6bf4fc3
[alexxy/gromacs.git] / src / gromacs / domdec / reversetopology.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2021, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  *
38  * \brief This file defines functions used in making the
39  * reverse topology.
40  *
41  * \author Berk Hess <hess@kth.se>
42  * \ingroup module_domdec
43  */
44
45 #include "gmxpre.h"
46
47 #include "gromacs/domdec/reversetopology.h"
48
49 #include <cstdio>
50
51 #include <memory>
52 #include <vector>
53
54 #include "gromacs/math/vec.h"
55 #include "gromacs/domdec/domdec_constraints.h"
56 #include "gromacs/domdec/domdec_internal.h"
57 #include "gromacs/domdec/domdec_vsite.h"
58 #include "gromacs/domdec/options.h"
59 #include "gromacs/mdlib/gmx_omp_nthreads.h"
60 #include "gromacs/mdlib/vsite.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/topology/atoms.h"
63 #include "gromacs/topology/mtop_util.h"
64 #include "gromacs/topology/topsort.h"
65 #include "gromacs/utility/arrayref.h"
66 #include "gromacs/utility/fatalerror.h"
67
68 using gmx::ArrayRef;
69 using gmx::DDBondedChecking;
70 using gmx::ListOfLists;
71 using gmx::RVec;
72
73 /*! \brief Struct for the reverse topology: links bonded interactions to atoms */
74 struct gmx_reverse_top_t::Impl
75 {
76     //! Constructs a reverse topology from \p mtop
77     Impl(const gmx_mtop_t& mtop, bool useFreeEnergy, const ReverseTopOptions& reverseTopOptions);
78
79     //! @cond Doxygen_Suppress
80     //! Options for the setup of this reverse topology
81     const ReverseTopOptions options;
82     //! Are there interaction of type F_POSRES and/or F_FBPOSRES
83     bool hasPositionRestraints;
84     //! \brief Are there bondeds/exclusions between atoms?
85     bool bInterAtomicInteractions = false;
86     //! \brief Reverse ilist for all moltypes
87     std::vector<reverse_ilist_t> ril_mt;
88     //! \brief The size of ril_mt[?].index summed over all entries
89     int ril_mt_tot_size = 0;
90     //! \brief Whether listed-force interaction lists should be sorted for free energy
91     bool doListedForcesSorting;
92     //! \brief molblock to global atom index for quick lookup of molblocks on atom index
93     std::vector<MolblockIndices> mbi;
94
95     //! \brief Do we have intermolecular interactions?
96     bool bIntermolecularInteractions = false;
97     //! \brief Intermolecular reverse ilist
98     reverse_ilist_t ril_intermol;
99
100     /* Work data structures for multi-threading */
101     //! \brief Thread work array for local topology generation
102     std::vector<thread_work_t> th_work;
103     //! @endcond
104 };
105
106
107 int nral_rt(int ftype)
108 {
109     int nral = NRAL(ftype);
110     if (interaction_function[ftype].flags & IF_VSITE)
111     {
112         /* With vsites the reverse topology contains an extra entry
113          * for storing if constructing atoms are vsites.
114          */
115         nral += 1;
116     }
117
118     return nral;
119 }
120
121 bool dd_check_ftype(const int ftype, const ReverseTopOptions& rtOptions)
122 {
123     return ((((interaction_function[ftype].flags & IF_BOND) != 0U)
124              && ((interaction_function[ftype].flags & IF_VSITE) == 0U)
125              && ((rtOptions.ddBondedChecking == DDBondedChecking::All)
126                  || ((interaction_function[ftype].flags & IF_LIMZERO) == 0U)))
127             || (rtOptions.includeConstraints && (ftype == F_CONSTR || ftype == F_CONSTRNC))
128             || (rtOptions.includeSettles && ftype == F_SETTLE));
129 }
130
131 void global_atomnr_to_moltype_ind(ArrayRef<const MolblockIndices> molblockIndices,
132                                   int                             i_gl,
133                                   int*                            mb,
134                                   int*                            mt,
135                                   int*                            mol,
136                                   int*                            i_mol)
137 {
138     const MolblockIndices* mbi   = molblockIndices.data();
139     int                    start = 0;
140     int                    end   = molblockIndices.size(); /* exclusive */
141     int                    mid   = 0;
142
143     /* binary search for molblock_ind */
144     while (TRUE)
145     {
146         mid = (start + end) / 2;
147         if (i_gl >= mbi[mid].a_end)
148         {
149             start = mid + 1;
150         }
151         else if (i_gl < mbi[mid].a_start)
152         {
153             end = mid;
154         }
155         else
156         {
157             break;
158         }
159     }
160
161     *mb = mid;
162     mbi += mid;
163
164     *mt    = mbi->type;
165     *mol   = (i_gl - mbi->a_start) / mbi->natoms_mol;
166     *i_mol = (i_gl - mbi->a_start) - (*mol) * mbi->natoms_mol;
167 }
168
169 /*! \brief Returns the maximum number of exclusions per atom */
170 static int getMaxNumExclusionsPerAtom(const ListOfLists<int>& excls)
171 {
172     int maxNumExcls = 0;
173     for (gmx::index at = 0; at < excls.ssize(); at++)
174     {
175         const auto list     = excls[at];
176         const int  numExcls = list.ssize();
177
178         GMX_RELEASE_ASSERT(numExcls != 1 || list[0] == at,
179                            "With 1 exclusion we expect a self-exclusion");
180
181         maxNumExcls = std::max(maxNumExcls, numExcls);
182     }
183
184     return maxNumExcls;
185 }
186
187 /*! \brief Run the reverse ilist generation and store it in r_il when \p bAssign = TRUE */
188 static void low_make_reverse_ilist(const InteractionLists&  il_mt,
189                                    const t_atom*            atom,
190                                    int*                     count,
191                                    const ReverseTopOptions& rtOptions,
192                                    gmx::ArrayRef<const int> r_index,
193                                    gmx::ArrayRef<int>       r_il,
194                                    const AtomLinkRule       atomLinkRule,
195                                    const bool               assignReverseIlist)
196 {
197     const bool includeConstraints = rtOptions.includeConstraints;
198     const bool includeSettles     = rtOptions.includeSettles;
199
200     for (int ftype = 0; ftype < F_NRE; ftype++)
201     {
202         if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE))
203             || (includeConstraints && (ftype == F_CONSTR || ftype == F_CONSTRNC))
204             || (includeSettles && ftype == F_SETTLE))
205         {
206             const bool  isVSite = ((interaction_function[ftype].flags & IF_VSITE) != 0U);
207             const int   nral    = NRAL(ftype);
208             const auto& il      = il_mt[ftype];
209             for (int i = 0; i < il.size(); i += 1 + nral)
210             {
211                 const int* ia = il.iatoms.data() + i;
212                 // Virtual sites should not be linked for bonded interactions
213                 const int nlink = (atomLinkRule == AtomLinkRule::FirstAtom) ? 1 : (isVSite ? 0 : nral);
214                 for (int link = 0; link < nlink; link++)
215                 {
216                     const int a = ia[1 + link];
217                     if (assignReverseIlist)
218                     {
219                         GMX_ASSERT(!r_il.empty(), "with bAssign not allowed to be empty");
220                         GMX_ASSERT(!r_index.empty(), "with bAssign not allowed to be empty");
221                         r_il[r_index[a] + count[a]]     = (ftype == F_CONSTRNC ? F_CONSTR : ftype);
222                         r_il[r_index[a] + count[a] + 1] = ia[0];
223                         for (int j = 1; j < 1 + nral; j++)
224                         {
225                             /* Store the molecular atom number */
226                             r_il[r_index[a] + count[a] + 1 + j] = ia[j];
227                         }
228                     }
229                     if (interaction_function[ftype].flags & IF_VSITE)
230                     {
231                         if (assignReverseIlist)
232                         {
233                             /* Add an entry to iatoms for storing
234                              * which of the constructing atoms are
235                              * vsites again.
236                              */
237                             r_il[r_index[a] + count[a] + 2 + nral] = 0;
238                             for (int j = 2; j < 1 + nral; j++)
239                             {
240                                 if (atom[ia[j]].ptype == ParticleType::VSite)
241                                 {
242                                     r_il[r_index[a] + count[a] + 2 + nral] |= (2 << j);
243                                 }
244                             }
245                         }
246                     }
247                     count[a] += 2 + nral_rt(ftype);
248                 }
249             }
250         }
251     }
252 }
253
254 void make_reverse_ilist(const InteractionLists&  ilist,
255                         const t_atoms*           atoms,
256                         const ReverseTopOptions& rtOptions,
257                         const AtomLinkRule       atomLinkRule,
258                         reverse_ilist_t*         ril_mt)
259 {
260     /* Count the interactions */
261     const int        nat_mt = atoms->nr;
262     std::vector<int> count(nat_mt);
263     low_make_reverse_ilist(ilist, atoms->atom, count.data(), rtOptions, {}, {}, atomLinkRule, FALSE);
264
265     ril_mt->index.push_back(0);
266     for (int i = 0; i < nat_mt; i++)
267     {
268         ril_mt->index.push_back(ril_mt->index[i] + count[i]);
269         count[i] = 0;
270     }
271     ril_mt->il.resize(ril_mt->index[nat_mt]);
272
273     /* Store the interactions */
274     low_make_reverse_ilist(
275             ilist, atoms->atom, count.data(), rtOptions, ril_mt->index, ril_mt->il, atomLinkRule, TRUE);
276
277     ril_mt->numAtomsInMolecule = atoms->nr;
278 }
279
280 gmx_reverse_top_t::gmx_reverse_top_t(const gmx_mtop_t&        mtop,
281                                      bool                     useFreeEnergy,
282                                      const ReverseTopOptions& reverseTopOptions) :
283     impl_(std::make_unique<Impl>(mtop, useFreeEnergy, reverseTopOptions))
284 {
285 }
286
287 gmx_reverse_top_t::~gmx_reverse_top_t() {}
288
289 const ReverseTopOptions& gmx_reverse_top_t::options() const
290 {
291     return impl_->options;
292 }
293
294 const reverse_ilist_t& gmx_reverse_top_t::interactionListForMoleculeType(int moleculeType) const
295 {
296     return impl_->ril_mt[moleculeType];
297 }
298
299 ArrayRef<const MolblockIndices> gmx_reverse_top_t::molblockIndices() const
300 {
301     return impl_->mbi;
302 }
303
304 bool gmx_reverse_top_t::hasIntermolecularInteractions() const
305 {
306     return impl_->bIntermolecularInteractions;
307 }
308
309 const reverse_ilist_t& gmx_reverse_top_t::interactionListForIntermolecularInteractions() const
310 {
311     return impl_->ril_intermol;
312 }
313
314 bool gmx_reverse_top_t::hasInterAtomicInteractions() const
315 {
316     return impl_->bInterAtomicInteractions;
317 }
318
319 bool gmx_reverse_top_t::hasPositionRestraints() const
320 {
321     return impl_->hasPositionRestraints;
322 }
323
324 ArrayRef<thread_work_t> gmx_reverse_top_t::threadWorkObjects() const
325 {
326     return impl_->th_work;
327 }
328
329 bool gmx_reverse_top_t::doListedForcesSorting() const
330 {
331     return impl_->doListedForcesSorting;
332 }
333
334 /*! \brief Generate the reverse topology */
335 gmx_reverse_top_t::Impl::Impl(const gmx_mtop_t&        mtop,
336                               const bool               useFreeEnergy,
337                               const ReverseTopOptions& reverseTopOptions) :
338     options(reverseTopOptions),
339     hasPositionRestraints(gmx_mtop_ftype_count(mtop, F_POSRES) + gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0),
340     bInterAtomicInteractions(mtop.bIntermolecularInteractions)
341 {
342     bInterAtomicInteractions = mtop.bIntermolecularInteractions;
343     ril_mt.resize(mtop.moltype.size());
344     ril_mt_tot_size = 0;
345     for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
346     {
347         const gmx_moltype_t& molt = mtop.moltype[mt];
348         if (molt.atoms.nr > 1)
349         {
350             bInterAtomicInteractions = true;
351         }
352
353         /* Make the atom to interaction list for this molecule type */
354         make_reverse_ilist(molt.ilist, &molt.atoms, options, AtomLinkRule::FirstAtom, &ril_mt[mt]);
355
356         ril_mt_tot_size += ril_mt[mt].index[molt.atoms.nr];
357     }
358     if (debug)
359     {
360         fprintf(debug, "The total size of the atom to interaction index is %d integers\n", ril_mt_tot_size);
361     }
362
363     /* Make an intermolecular reverse top, if necessary */
364     bIntermolecularInteractions = mtop.bIntermolecularInteractions;
365     if (bIntermolecularInteractions)
366     {
367         t_atoms atoms_global;
368
369         atoms_global.nr   = mtop.natoms;
370         atoms_global.atom = nullptr; /* Only used with virtual sites */
371
372         GMX_RELEASE_ASSERT(mtop.intermolecular_ilist,
373                            "We should have an ilist when intermolecular interactions are on");
374
375         make_reverse_ilist(
376                 *mtop.intermolecular_ilist, &atoms_global, options, AtomLinkRule::FirstAtom, &ril_intermol);
377     }
378
379     doListedForcesSorting = useFreeEnergy && gmx_mtop_bondeds_free_energy(&mtop);
380
381     /* Make a molblock index for fast searching */
382     int i = 0;
383     for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
384     {
385         const gmx_molblock_t& molb           = mtop.molblock[mb];
386         const int             numAtomsPerMol = mtop.moltype[molb.type].atoms.nr;
387         MolblockIndices       mbiMolblock;
388         mbiMolblock.a_start = i;
389         i += molb.nmol * numAtomsPerMol;
390         mbiMolblock.a_end      = i;
391         mbiMolblock.natoms_mol = numAtomsPerMol;
392         mbiMolblock.type       = molb.type;
393         mbi.push_back(mbiMolblock);
394     }
395
396     for (int th = 0; th < gmx_omp_nthreads_get(ModuleMultiThread::Domdec); th++)
397     {
398         th_work.emplace_back(mtop.ffparams);
399     }
400 }
401
402 void dd_make_reverse_top(FILE*                           fplog,
403                          gmx_domdec_t*                   dd,
404                          const gmx_mtop_t&               mtop,
405                          const gmx::VirtualSitesHandler* vsite,
406                          const t_inputrec&               inputrec,
407                          const DDBondedChecking          ddBondedChecking)
408 {
409     if (fplog)
410     {
411         fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
412     }
413
414     /* If normal and/or settle constraints act only within charge groups,
415      * we can store them in the reverse top and simply assign them to domains.
416      * Otherwise we need to assign them to multiple domains and set up
417      * the parallel version constraint algorithm(s).
418      */
419     GMX_RELEASE_ASSERT(ddBondedChecking == DDBondedChecking::ExcludeZeroLimit
420                                || ddBondedChecking == DDBondedChecking::All,
421                        "Invalid enum value for mdrun -ddcheck");
422     const ReverseTopOptions rtOptions(ddBondedChecking,
423                                       !dd->comm->systemInfo.mayHaveSplitConstraints,
424                                       !dd->comm->systemInfo.mayHaveSplitSettles);
425
426     dd->reverse_top = std::make_unique<gmx_reverse_top_t>(
427             mtop, inputrec.efep != FreeEnergyPerturbationType::No, rtOptions);
428
429     dd->haveExclusions = false;
430     for (const gmx_molblock_t& molb : mtop.molblock)
431     {
432         const int maxNumExclusionsPerAtom = getMaxNumExclusionsPerAtom(mtop.moltype[molb.type].excls);
433         // We checked above that max 1 exclusion means only self exclusions
434         if (maxNumExclusionsPerAtom > 1)
435         {
436             dd->haveExclusions = true;
437         }
438     }
439
440     const int numInterUpdategroupVirtualSites =
441             (vsite == nullptr ? 0 : vsite->numInterUpdategroupVirtualSites());
442     if (numInterUpdategroupVirtualSites > 0)
443     {
444         if (fplog)
445         {
446             fprintf(fplog,
447                     "There are %d inter update-group virtual sites,\n"
448                     "will perform an extra communication step for selected coordinates and "
449                     "forces\n",
450                     numInterUpdategroupVirtualSites);
451         }
452         init_domdec_vsites(dd, numInterUpdategroupVirtualSites);
453     }
454
455     if (dd->comm->systemInfo.mayHaveSplitConstraints || dd->comm->systemInfo.mayHaveSplitSettles)
456     {
457         init_domdec_constraints(dd, mtop);
458     }
459     if (fplog)
460     {
461         fprintf(fplog, "\n");
462     }
463 }