0601c441f7b2d59f39c489a5ee0f0e5c5f3cd8c8
[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 The sorting state of bondeds for free energy
91     int ilsort = ilsortUNKNOWN;
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     //! The number of bonded interactions computed from the full topology
101     int expectedNumGlobalBondedInteractions = 0;
102
103     /* Work data structures for multi-threading */
104     //! \brief Thread work array for local topology generation
105     std::vector<thread_work_t> th_work;
106     //! @endcond
107 };
108
109
110 int nral_rt(int ftype)
111 {
112     int nral = NRAL(ftype);
113     if (interaction_function[ftype].flags & IF_VSITE)
114     {
115         /* With vsites the reverse topology contains an extra entry
116          * for storing if constructing atoms are vsites.
117          */
118         nral += 1;
119     }
120
121     return nral;
122 }
123
124 bool dd_check_ftype(const int ftype, const ReverseTopOptions& rtOptions)
125 {
126     return ((((interaction_function[ftype].flags & IF_BOND) != 0U)
127              && ((interaction_function[ftype].flags & IF_VSITE) == 0U)
128              && ((rtOptions.ddBondedChecking == DDBondedChecking::All)
129                  || ((interaction_function[ftype].flags & IF_LIMZERO) == 0U)))
130             || (rtOptions.includeConstraints && (ftype == F_CONSTR || ftype == F_CONSTRNC))
131             || (rtOptions.includeSettles && ftype == F_SETTLE));
132 }
133
134 void global_atomnr_to_moltype_ind(ArrayRef<const MolblockIndices> molblockIndices,
135                                   int                             i_gl,
136                                   int*                            mb,
137                                   int*                            mt,
138                                   int*                            mol,
139                                   int*                            i_mol)
140 {
141     const MolblockIndices* mbi   = molblockIndices.data();
142     int                    start = 0;
143     int                    end   = molblockIndices.size(); /* exclusive */
144     int                    mid   = 0;
145
146     /* binary search for molblock_ind */
147     while (TRUE)
148     {
149         mid = (start + end) / 2;
150         if (i_gl >= mbi[mid].a_end)
151         {
152             start = mid + 1;
153         }
154         else if (i_gl < mbi[mid].a_start)
155         {
156             end = mid;
157         }
158         else
159         {
160             break;
161         }
162     }
163
164     *mb = mid;
165     mbi += mid;
166
167     *mt    = mbi->type;
168     *mol   = (i_gl - mbi->a_start) / mbi->natoms_mol;
169     *i_mol = (i_gl - mbi->a_start) - (*mol) * mbi->natoms_mol;
170 }
171
172 /*! \brief Returns the maximum number of exclusions per atom */
173 static int getMaxNumExclusionsPerAtom(const ListOfLists<int>& excls)
174 {
175     int maxNumExcls = 0;
176     for (gmx::index at = 0; at < excls.ssize(); at++)
177     {
178         const auto list     = excls[at];
179         const int  numExcls = list.ssize();
180
181         GMX_RELEASE_ASSERT(numExcls != 1 || list[0] == at,
182                            "With 1 exclusion we expect a self-exclusion");
183
184         maxNumExcls = std::max(maxNumExcls, numExcls);
185     }
186
187     return maxNumExcls;
188 }
189
190 /*! \brief Run the reverse ilist generation and store it in r_il when \p bAssign = TRUE */
191 static void low_make_reverse_ilist(const InteractionLists&  il_mt,
192                                    const t_atom*            atom,
193                                    int*                     count,
194                                    const ReverseTopOptions& rtOptions,
195                                    gmx::ArrayRef<const int> r_index,
196                                    gmx::ArrayRef<int>       r_il,
197                                    const AtomLinkRule       atomLinkRule,
198                                    const bool               assignReverseIlist)
199 {
200     const bool includeConstraints = rtOptions.includeConstraints;
201     const bool includeSettles     = rtOptions.includeSettles;
202
203     for (int ftype = 0; ftype < F_NRE; ftype++)
204     {
205         if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE))
206             || (includeConstraints && (ftype == F_CONSTR || ftype == F_CONSTRNC))
207             || (includeSettles && ftype == F_SETTLE))
208         {
209             const bool  isVSite = ((interaction_function[ftype].flags & IF_VSITE) != 0U);
210             const int   nral    = NRAL(ftype);
211             const auto& il      = il_mt[ftype];
212             for (int i = 0; i < il.size(); i += 1 + nral)
213             {
214                 const int* ia = il.iatoms.data() + i;
215                 // Virtual sites should not be linked for bonded interactions
216                 const int nlink = (atomLinkRule == AtomLinkRule::FirstAtom) ? 1 : (isVSite ? 0 : nral);
217                 for (int link = 0; link < nlink; link++)
218                 {
219                     const int a = ia[1 + link];
220                     if (assignReverseIlist)
221                     {
222                         GMX_ASSERT(!r_il.empty(), "with bAssign not allowed to be empty");
223                         GMX_ASSERT(!r_index.empty(), "with bAssign not allowed to be empty");
224                         r_il[r_index[a] + count[a]]     = (ftype == F_CONSTRNC ? F_CONSTR : ftype);
225                         r_il[r_index[a] + count[a] + 1] = ia[0];
226                         for (int j = 1; j < 1 + nral; j++)
227                         {
228                             /* Store the molecular atom number */
229                             r_il[r_index[a] + count[a] + 1 + j] = ia[j];
230                         }
231                     }
232                     if (interaction_function[ftype].flags & IF_VSITE)
233                     {
234                         if (assignReverseIlist)
235                         {
236                             /* Add an entry to iatoms for storing
237                              * which of the constructing atoms are
238                              * vsites again.
239                              */
240                             r_il[r_index[a] + count[a] + 2 + nral] = 0;
241                             for (int j = 2; j < 1 + nral; j++)
242                             {
243                                 if (atom[ia[j]].ptype == ParticleType::VSite)
244                                 {
245                                     r_il[r_index[a] + count[a] + 2 + nral] |= (2 << j);
246                                 }
247                             }
248                         }
249                     }
250                     count[a] += 2 + nral_rt(ftype);
251                 }
252             }
253         }
254     }
255 }
256
257 void make_reverse_ilist(const InteractionLists&  ilist,
258                         const t_atoms*           atoms,
259                         const ReverseTopOptions& rtOptions,
260                         const AtomLinkRule       atomLinkRule,
261                         reverse_ilist_t*         ril_mt)
262 {
263     /* Count the interactions */
264     const int        nat_mt = atoms->nr;
265     std::vector<int> count(nat_mt);
266     low_make_reverse_ilist(ilist, atoms->atom, count.data(), rtOptions, {}, {}, atomLinkRule, FALSE);
267
268     ril_mt->index.push_back(0);
269     for (int i = 0; i < nat_mt; i++)
270     {
271         ril_mt->index.push_back(ril_mt->index[i] + count[i]);
272         count[i] = 0;
273     }
274     ril_mt->il.resize(ril_mt->index[nat_mt]);
275
276     /* Store the interactions */
277     low_make_reverse_ilist(
278             ilist, atoms->atom, count.data(), rtOptions, ril_mt->index, ril_mt->il, atomLinkRule, TRUE);
279
280     ril_mt->numAtomsInMolecule = atoms->nr;
281 }
282
283 gmx_reverse_top_t::gmx_reverse_top_t(const gmx_mtop_t&        mtop,
284                                      bool                     useFreeEnergy,
285                                      const ReverseTopOptions& reverseTopOptions) :
286     impl_(std::make_unique<Impl>(mtop, useFreeEnergy, reverseTopOptions))
287 {
288 }
289
290 gmx_reverse_top_t::~gmx_reverse_top_t() {}
291
292 const ReverseTopOptions& gmx_reverse_top_t::options() const
293 {
294     return impl_->options;
295 }
296
297 const reverse_ilist_t& gmx_reverse_top_t::interactionListForMoleculeType(int moleculeType) const
298 {
299     return impl_->ril_mt[moleculeType];
300 }
301
302 int gmx_reverse_top_t::expectedNumGlobalBondedInteractions() const
303 {
304     return impl_->expectedNumGlobalBondedInteractions;
305 }
306
307 ArrayRef<const MolblockIndices> gmx_reverse_top_t::molblockIndices() const
308 {
309     return impl_->mbi;
310 }
311
312 bool gmx_reverse_top_t::hasIntermolecularInteractions() const
313 {
314     return impl_->bIntermolecularInteractions;
315 }
316
317 const reverse_ilist_t& gmx_reverse_top_t::interactionListForIntermolecularInteractions() const
318 {
319     return impl_->ril_intermol;
320 }
321
322 bool gmx_reverse_top_t::hasInterAtomicInteractions() const
323 {
324     return impl_->bInterAtomicInteractions;
325 }
326
327 bool gmx_reverse_top_t::hasPositionRestraints() const
328 {
329     return impl_->hasPositionRestraints;
330 }
331
332 ArrayRef<thread_work_t> gmx_reverse_top_t::threadWorkObjects() const
333 {
334     return impl_->th_work;
335 }
336
337 bool gmx_reverse_top_t::doSorting() const
338 {
339     return impl_->ilsort != ilsortNO_FE;
340 }
341
342 /*! \brief Compute the total bonded interaction count
343  *
344  * \param[in] mtop                The global system topology
345  * \param[in] includeConstraints  Whether constraint interactions are included in the count
346  * \param[in] includeSettles      Whether settle interactions are included in the count
347  *
348  * When using domain decomposition without update groups,
349  * constraint-type interactions can be split across domains, and so we
350  * do not consider them in this correctness check. Otherwise, we
351  * include them.
352  */
353 static int computeExpectedNumGlobalBondedInteractions(const gmx_mtop_t& mtop,
354                                                       const bool        includeConstraints,
355                                                       const bool        includeSettles)
356 {
357     int expectedNumGlobalBondedInteractions = gmx_mtop_interaction_count(mtop, IF_BOND);
358     if (includeConstraints)
359     {
360         expectedNumGlobalBondedInteractions +=
361                 gmx_mtop_ftype_count(mtop, F_CONSTR) + gmx_mtop_ftype_count(mtop, F_CONSTRNC);
362     }
363     if (includeSettles)
364     {
365         expectedNumGlobalBondedInteractions += gmx_mtop_ftype_count(mtop, F_SETTLE);
366     }
367     return expectedNumGlobalBondedInteractions;
368 }
369
370 /*! \brief Generate the reverse topology */
371 gmx_reverse_top_t::Impl::Impl(const gmx_mtop_t&        mtop,
372                               const bool               useFreeEnergy,
373                               const ReverseTopOptions& reverseTopOptions) :
374     options(reverseTopOptions),
375     hasPositionRestraints(gmx_mtop_ftype_count(mtop, F_POSRES) + gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0),
376     bInterAtomicInteractions(mtop.bIntermolecularInteractions),
377     expectedNumGlobalBondedInteractions(
378             computeExpectedNumGlobalBondedInteractions(mtop,
379                                                        reverseTopOptions.includeConstraints,
380                                                        reverseTopOptions.includeSettles))
381 {
382     bInterAtomicInteractions = mtop.bIntermolecularInteractions;
383     ril_mt.resize(mtop.moltype.size());
384     ril_mt_tot_size = 0;
385     for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
386     {
387         const gmx_moltype_t& molt = mtop.moltype[mt];
388         if (molt.atoms.nr > 1)
389         {
390             bInterAtomicInteractions = true;
391         }
392
393         /* Make the atom to interaction list for this molecule type */
394         make_reverse_ilist(molt.ilist, &molt.atoms, options, AtomLinkRule::FirstAtom, &ril_mt[mt]);
395
396         ril_mt_tot_size += ril_mt[mt].index[molt.atoms.nr];
397     }
398     if (debug)
399     {
400         fprintf(debug, "The total size of the atom to interaction index is %d integers\n", ril_mt_tot_size);
401     }
402
403     /* Make an intermolecular reverse top, if necessary */
404     bIntermolecularInteractions = mtop.bIntermolecularInteractions;
405     if (bIntermolecularInteractions)
406     {
407         t_atoms atoms_global;
408
409         atoms_global.nr   = mtop.natoms;
410         atoms_global.atom = nullptr; /* Only used with virtual sites */
411
412         GMX_RELEASE_ASSERT(mtop.intermolecular_ilist,
413                            "We should have an ilist when intermolecular interactions are on");
414
415         make_reverse_ilist(
416                 *mtop.intermolecular_ilist, &atoms_global, options, AtomLinkRule::FirstAtom, &ril_intermol);
417     }
418
419     if (useFreeEnergy && gmx_mtop_bondeds_free_energy(&mtop))
420     {
421         ilsort = ilsortFE_UNSORTED;
422     }
423     else
424     {
425         ilsort = ilsortNO_FE;
426     }
427
428     /* Make a molblock index for fast searching */
429     int i = 0;
430     for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
431     {
432         const gmx_molblock_t& molb           = mtop.molblock[mb];
433         const int             numAtomsPerMol = mtop.moltype[molb.type].atoms.nr;
434         MolblockIndices       mbiMolblock;
435         mbiMolblock.a_start = i;
436         i += molb.nmol * numAtomsPerMol;
437         mbiMolblock.a_end      = i;
438         mbiMolblock.natoms_mol = numAtomsPerMol;
439         mbiMolblock.type       = molb.type;
440         mbi.push_back(mbiMolblock);
441     }
442
443     for (int th = 0; th < gmx_omp_nthreads_get(ModuleMultiThread::Domdec); th++)
444     {
445         th_work.emplace_back(mtop.ffparams);
446     }
447 }
448
449 void dd_make_reverse_top(FILE*                           fplog,
450                          gmx_domdec_t*                   dd,
451                          const gmx_mtop_t&               mtop,
452                          const gmx::VirtualSitesHandler* vsite,
453                          const t_inputrec&               inputrec,
454                          const DDBondedChecking          ddBondedChecking)
455 {
456     if (fplog)
457     {
458         fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
459     }
460
461     /* If normal and/or settle constraints act only within charge groups,
462      * we can store them in the reverse top and simply assign them to domains.
463      * Otherwise we need to assign them to multiple domains and set up
464      * the parallel version constraint algorithm(s).
465      */
466     GMX_RELEASE_ASSERT(ddBondedChecking == DDBondedChecking::ExcludeZeroLimit
467                                || ddBondedChecking == DDBondedChecking::All,
468                        "Invalid enum value for mdrun -ddcheck");
469     const ReverseTopOptions rtOptions(ddBondedChecking,
470                                       !dd->comm->systemInfo.mayHaveSplitConstraints,
471                                       !dd->comm->systemInfo.mayHaveSplitSettles);
472
473     dd->reverse_top = std::make_unique<gmx_reverse_top_t>(
474             mtop, inputrec.efep != FreeEnergyPerturbationType::No, rtOptions);
475
476     dd->haveExclusions = false;
477     for (const gmx_molblock_t& molb : mtop.molblock)
478     {
479         const int maxNumExclusionsPerAtom = getMaxNumExclusionsPerAtom(mtop.moltype[molb.type].excls);
480         // We checked above that max 1 exclusion means only self exclusions
481         if (maxNumExclusionsPerAtom > 1)
482         {
483             dd->haveExclusions = true;
484         }
485     }
486
487     const int numInterUpdategroupVirtualSites =
488             (vsite == nullptr ? 0 : vsite->numInterUpdategroupVirtualSites());
489     if (numInterUpdategroupVirtualSites > 0)
490     {
491         if (fplog)
492         {
493             fprintf(fplog,
494                     "There are %d inter update-group virtual sites,\n"
495                     "will perform an extra communication step for selected coordinates and "
496                     "forces\n",
497                     numInterUpdategroupVirtualSites);
498         }
499         init_domdec_vsites(dd, numInterUpdategroupVirtualSites);
500     }
501
502     if (dd->comm->systemInfo.mayHaveSplitConstraints || dd->comm->systemInfo.mayHaveSplitSettles)
503     {
504         init_domdec_constraints(dd, mtop);
505     }
506     if (fplog)
507     {
508         fprintf(fplog, "\n");
509     }
510 }