Remove orphaned code
[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 int 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     const DDBondedChecking ddBondedChecking   = rtOptions.ddBondedChecking;
203
204     int nint = 0;
205
206     for (int ftype = 0; ftype < F_NRE; ftype++)
207     {
208         if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE))
209             || (includeConstraints && (ftype == F_CONSTR || ftype == F_CONSTRNC))
210             || (includeSettles && ftype == F_SETTLE))
211         {
212             const bool  isVSite = ((interaction_function[ftype].flags & IF_VSITE) != 0U);
213             const int   nral    = NRAL(ftype);
214             const auto& il      = il_mt[ftype];
215             for (int i = 0; i < il.size(); i += 1 + nral)
216             {
217                 const int* ia = il.iatoms.data() + i;
218                 // Virtual sites should not be linked for bonded interactions
219                 const int nlink = (atomLinkRule == AtomLinkRule::FirstAtom) ? 1 : (isVSite ? 0 : nral);
220                 for (int link = 0; link < nlink; link++)
221                 {
222                     const int a = ia[1 + link];
223                     if (assignReverseIlist)
224                     {
225                         GMX_ASSERT(!r_il.empty(), "with bAssign not allowed to be empty");
226                         GMX_ASSERT(!r_index.empty(), "with bAssign not allowed to be empty");
227                         r_il[r_index[a] + count[a]]     = (ftype == F_CONSTRNC ? F_CONSTR : ftype);
228                         r_il[r_index[a] + count[a] + 1] = ia[0];
229                         for (int j = 1; j < 1 + nral; j++)
230                         {
231                             /* Store the molecular atom number */
232                             r_il[r_index[a] + count[a] + 1 + j] = ia[j];
233                         }
234                     }
235                     if (interaction_function[ftype].flags & IF_VSITE)
236                     {
237                         if (assignReverseIlist)
238                         {
239                             /* Add an entry to iatoms for storing
240                              * which of the constructing atoms are
241                              * vsites again.
242                              */
243                             r_il[r_index[a] + count[a] + 2 + nral] = 0;
244                             for (int j = 2; j < 1 + nral; j++)
245                             {
246                                 if (atom[ia[j]].ptype == ParticleType::VSite)
247                                 {
248                                     r_il[r_index[a] + count[a] + 2 + nral] |= (2 << j);
249                                 }
250                             }
251                         }
252                     }
253                     else
254                     {
255                         /* We do not count vsites since they are always
256                          * uniquely assigned and can be assigned
257                          * to multiple nodes with recursive vsites.
258                          */
259                         if (ddBondedChecking == DDBondedChecking::All
260                             || !(interaction_function[ftype].flags & IF_LIMZERO))
261                         {
262                             nint++;
263                         }
264                     }
265                     count[a] += 2 + nral_rt(ftype);
266                 }
267             }
268         }
269     }
270
271     return nint;
272 }
273
274 int make_reverse_ilist(const InteractionLists&  ilist,
275                        const t_atoms*           atoms,
276                        const ReverseTopOptions& rtOptions,
277                        const AtomLinkRule       atomLinkRule,
278                        reverse_ilist_t*         ril_mt)
279 {
280     /* Count the interactions */
281     const int        nat_mt = atoms->nr;
282     std::vector<int> count(nat_mt);
283     low_make_reverse_ilist(ilist, atoms->atom, count.data(), rtOptions, {}, {}, atomLinkRule, FALSE);
284
285     ril_mt->index.push_back(0);
286     for (int i = 0; i < nat_mt; i++)
287     {
288         ril_mt->index.push_back(ril_mt->index[i] + count[i]);
289         count[i] = 0;
290     }
291     ril_mt->il.resize(ril_mt->index[nat_mt]);
292
293     /* Store the interactions */
294     int nint_mt = low_make_reverse_ilist(
295             ilist, atoms->atom, count.data(), rtOptions, ril_mt->index, ril_mt->il, atomLinkRule, TRUE);
296
297     ril_mt->numAtomsInMolecule = atoms->nr;
298
299     return nint_mt;
300 }
301
302 gmx_reverse_top_t::gmx_reverse_top_t(const gmx_mtop_t&        mtop,
303                                      bool                     useFreeEnergy,
304                                      const ReverseTopOptions& reverseTopOptions) :
305     impl_(std::make_unique<Impl>(mtop, useFreeEnergy, reverseTopOptions))
306 {
307 }
308
309 gmx_reverse_top_t::~gmx_reverse_top_t() {}
310
311 const ReverseTopOptions& gmx_reverse_top_t::options() const
312 {
313     return impl_->options;
314 }
315
316 const reverse_ilist_t& gmx_reverse_top_t::interactionListForMoleculeType(int moleculeType) const
317 {
318     return impl_->ril_mt[moleculeType];
319 }
320
321 int gmx_reverse_top_t::expectedNumGlobalBondedInteractions() const
322 {
323     return impl_->expectedNumGlobalBondedInteractions;
324 }
325
326 ArrayRef<const MolblockIndices> gmx_reverse_top_t::molblockIndices() const
327 {
328     return impl_->mbi;
329 }
330
331 bool gmx_reverse_top_t::hasIntermolecularInteractions() const
332 {
333     return impl_->bIntermolecularInteractions;
334 }
335
336 const reverse_ilist_t& gmx_reverse_top_t::interactionListForIntermolecularInteractions() const
337 {
338     return impl_->ril_intermol;
339 }
340
341 bool gmx_reverse_top_t::hasInterAtomicInteractions() const
342 {
343     return impl_->bInterAtomicInteractions;
344 }
345
346 bool gmx_reverse_top_t::hasPositionRestraints() const
347 {
348     return impl_->hasPositionRestraints;
349 }
350
351 ArrayRef<thread_work_t> gmx_reverse_top_t::threadWorkObjects() const
352 {
353     return impl_->th_work;
354 }
355
356 bool gmx_reverse_top_t::doSorting() const
357 {
358     return impl_->ilsort != ilsortNO_FE;
359 }
360
361 /*! \brief Generate the reverse topology */
362 gmx_reverse_top_t::Impl::Impl(const gmx_mtop_t&        mtop,
363                               const bool               useFreeEnergy,
364                               const ReverseTopOptions& reverseTopOptions) :
365     options(reverseTopOptions),
366     hasPositionRestraints(gmx_mtop_ftype_count(mtop, F_POSRES) + gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0),
367     bInterAtomicInteractions(mtop.bIntermolecularInteractions)
368 {
369     bInterAtomicInteractions = mtop.bIntermolecularInteractions;
370     ril_mt.resize(mtop.moltype.size());
371     ril_mt_tot_size = 0;
372     std::vector<int> nint_mt;
373     for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
374     {
375         const gmx_moltype_t& molt = mtop.moltype[mt];
376         if (molt.atoms.nr > 1)
377         {
378             bInterAtomicInteractions = true;
379         }
380
381         /* Make the atom to interaction list for this molecule type */
382         int numberOfInteractions = make_reverse_ilist(
383                 molt.ilist, &molt.atoms, options, AtomLinkRule::FirstAtom, &ril_mt[mt]);
384         nint_mt.push_back(numberOfInteractions);
385
386         ril_mt_tot_size += ril_mt[mt].index[molt.atoms.nr];
387     }
388     if (debug)
389     {
390         fprintf(debug, "The total size of the atom to interaction index is %d integers\n", ril_mt_tot_size);
391     }
392
393     expectedNumGlobalBondedInteractions = 0;
394     for (const gmx_molblock_t& molblock : mtop.molblock)
395     {
396         expectedNumGlobalBondedInteractions += molblock.nmol * nint_mt[molblock.type];
397     }
398
399     /* Make an intermolecular reverse top, if necessary */
400     bIntermolecularInteractions = mtop.bIntermolecularInteractions;
401     if (bIntermolecularInteractions)
402     {
403         t_atoms atoms_global;
404
405         atoms_global.nr   = mtop.natoms;
406         atoms_global.atom = nullptr; /* Only used with virtual sites */
407
408         GMX_RELEASE_ASSERT(mtop.intermolecular_ilist,
409                            "We should have an ilist when intermolecular interactions are on");
410
411         expectedNumGlobalBondedInteractions += make_reverse_ilist(
412                 *mtop.intermolecular_ilist, &atoms_global, options, AtomLinkRule::FirstAtom, &ril_intermol);
413     }
414
415     if (useFreeEnergy && gmx_mtop_bondeds_free_energy(&mtop))
416     {
417         ilsort = ilsortFE_UNSORTED;
418     }
419     else
420     {
421         ilsort = ilsortNO_FE;
422     }
423
424     /* Make a molblock index for fast searching */
425     int i = 0;
426     for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
427     {
428         const gmx_molblock_t& molb           = mtop.molblock[mb];
429         const int             numAtomsPerMol = mtop.moltype[molb.type].atoms.nr;
430         MolblockIndices       mbiMolblock;
431         mbiMolblock.a_start = i;
432         i += molb.nmol * numAtomsPerMol;
433         mbiMolblock.a_end      = i;
434         mbiMolblock.natoms_mol = numAtomsPerMol;
435         mbiMolblock.type       = molb.type;
436         mbi.push_back(mbiMolblock);
437     }
438
439     for (int th = 0; th < gmx_omp_nthreads_get(ModuleMultiThread::Domdec); th++)
440     {
441         th_work.emplace_back(mtop.ffparams);
442     }
443 }
444
445 void dd_make_reverse_top(FILE*                           fplog,
446                          gmx_domdec_t*                   dd,
447                          const gmx_mtop_t&               mtop,
448                          const gmx::VirtualSitesHandler* vsite,
449                          const t_inputrec&               inputrec,
450                          const DDBondedChecking          ddBondedChecking)
451 {
452     if (fplog)
453     {
454         fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
455     }
456
457     /* If normal and/or settle constraints act only within charge groups,
458      * we can store them in the reverse top and simply assign them to domains.
459      * Otherwise we need to assign them to multiple domains and set up
460      * the parallel version constraint algorithm(s).
461      */
462     GMX_RELEASE_ASSERT(ddBondedChecking == DDBondedChecking::ExcludeZeroLimit
463                                || ddBondedChecking == DDBondedChecking::All,
464                        "Invalid enum value for mdrun -ddcheck");
465     const ReverseTopOptions rtOptions(ddBondedChecking,
466                                       !dd->comm->systemInfo.mayHaveSplitConstraints,
467                                       !dd->comm->systemInfo.mayHaveSplitSettles);
468
469     dd->reverse_top = std::make_unique<gmx_reverse_top_t>(
470             mtop, inputrec.efep != FreeEnergyPerturbationType::No, rtOptions);
471
472     dd->haveExclusions = false;
473     for (const gmx_molblock_t& molb : mtop.molblock)
474     {
475         const int maxNumExclusionsPerAtom = getMaxNumExclusionsPerAtom(mtop.moltype[molb.type].excls);
476         // We checked above that max 1 exclusion means only self exclusions
477         if (maxNumExclusionsPerAtom > 1)
478         {
479             dd->haveExclusions = true;
480         }
481     }
482
483     const int numInterUpdategroupVirtualSites =
484             (vsite == nullptr ? 0 : vsite->numInterUpdategroupVirtualSites());
485     if (numInterUpdategroupVirtualSites > 0)
486     {
487         if (fplog)
488         {
489             fprintf(fplog,
490                     "There are %d inter update-group virtual sites,\n"
491                     "will perform an extra communication step for selected coordinates and "
492                     "forces\n",
493                     numInterUpdategroupVirtualSites);
494         }
495         init_domdec_vsites(dd, numInterUpdategroupVirtualSites);
496     }
497
498     if (dd->comm->systemInfo.mayHaveSplitConstraints || dd->comm->systemInfo.mayHaveSplitSettles)
499     {
500         init_domdec_constraints(dd, mtop);
501     }
502     if (fplog)
503     {
504         fprintf(fplog, "\n");
505     }
506 }