1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2006 - 2014, The GROMACS development team.
5  * Copyright (c) 2015,2016,2017,2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 
37 /*! \internal \file
38  *
39  * \brief This file defines functions used by the domdec module
40  * while managing the construction, use and error checking for
41  * topologies local to a DD rank.
42  *
43  * \author Berk Hess <hess@kth.se>
44  * \ingroup module_domdec
45  */
46 
47 #include "gmxpre.h"
48 
49 #include <cassert>
50 #include <cstdlib>
51 #include <cstring>
52 
53 #include <algorithm>
54 #include <memory>
55 #include <string>
56 
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_network.h"
59 #include "gromacs/domdec/ga2la.h"
60 #include "gromacs/gmxlib/network.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/mdlib/forcerec.h"
63 #include "gromacs/mdlib/gmx_omp_nthreads.h"
64 #include "gromacs/mdlib/vsite.h"
65 #include "gromacs/mdtypes/commrec.h"
66 #include "gromacs/mdtypes/forcerec.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/md_enums.h"
69 #include "gromacs/mdtypes/mdatom.h"
70 #include "gromacs/mdtypes/state.h"
71 #include "gromacs/pbcutil/mshift.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/topology/ifunc.h"
74 #include "gromacs/topology/mtop_util.h"
75 #include "gromacs/topology/topsort.h"
76 #include "gromacs/utility/cstringutil.h"
77 #include "gromacs/utility/exceptions.h"
78 #include "gromacs/utility/fatalerror.h"
79 #include "gromacs/utility/gmxassert.h"
80 #include "gromacs/utility/listoflists.h"
81 #include "gromacs/utility/logger.h"
82 #include "gromacs/utility/smalloc.h"
83 #include "gromacs/utility/strconvert.h"
84 #include "gromacs/utility/stringstream.h"
85 #include "gromacs/utility/stringutil.h"
86 #include "gromacs/utility/textwriter.h"
87 
88 #include "domdec_constraints.h"
89 #include "domdec_internal.h"
90 #include "domdec_vsite.h"
91 #include "dump.h"
92 
93 using gmx::ArrayRef;
94 using gmx::ListOfLists;
95 using gmx::RVec;
96 
97 /*! \brief The number of integer item in the local state, used for broadcasting of the state */
98 #define NITEM_DD_INIT_LOCAL_STATE 5
99 
100 struct reverse_ilist_t
101 {
102     std::vector<int> index;              /* Index for each atom into il          */
103     std::vector<int> il;                 /* ftype|type|a0|...|an|ftype|...       */
104     int              numAtomsInMolecule; /* The number of atoms in this molecule */
105 };
106 
107 struct MolblockIndices
108 {
109     int a_start;
110     int a_end;
111     int natoms_mol;
112     int type;
113 };
114 
115 /*! \brief Struct for thread local work data for local topology generation */
116 struct thread_work_t
117 {
118     /*! \brief Constructor
119      *
120      * \param[in] ffparams  The interaction parameters, the lifetime of the created object should not exceed the lifetime of the passed parameters
121      */
thread_work_tthread_work_t122     thread_work_t(const gmx_ffparams_t& ffparams) : idef(ffparams) {}
123 
124     InteractionDefinitions         idef;               /**< Partial local topology */
125     std::unique_ptr<gmx::VsitePbc> vsitePbc = nullptr; /**< vsite PBC structure */
126     int                            nbonded  = 0;       /**< The number of bondeds in this struct */
127     ListOfLists<int>               excl;               /**< List of exclusions */
128     int                            excl_count = 0;     /**< The total exclusion count for \p excl */
129 };
130 
131 /*! \brief Struct for the reverse topology: links bonded interactions to atomsx */
132 struct gmx_reverse_top_t
133 {
134     //! @cond Doxygen_Suppress
135     //! \brief Are there constraints in this revserse top?
136     bool bConstr = false;
137     //! \brief Are there settles in this revserse top?
138     bool bSettle = false;
139     //! \brief All bonded interactions have to be assigned?
140     bool bBCheck = false;
141     //! \brief Are there bondeds/exclusions between atoms?
142     bool bInterAtomicInteractions = false;
143     //! \brief Reverse ilist for all moltypes
144     std::vector<reverse_ilist_t> ril_mt;
145     //! \brief The size of ril_mt[?].index summed over all entries
146     int ril_mt_tot_size = 0;
147     //! \brief The sorting state of bondeds for free energy
148     int ilsort = ilsortUNKNOWN;
149     //! \brief molblock to global atom index for quick lookup of molblocks on atom index
150     std::vector<MolblockIndices> mbi;
151 
152     //! \brief Do we have intermolecular interactions?
153     bool bIntermolecularInteractions = false;
154     //! \brief Intermolecular reverse ilist
155     reverse_ilist_t ril_intermol;
156 
157     /* Work data structures for multi-threading */
158     //! \brief Thread work array for local topology generation
159     std::vector<thread_work_t> th_work;
160     //! @endcond
161 };
162 
163 /*! \brief Returns the number of atom entries for il in gmx_reverse_top_t */
nral_rt(int ftype)164 static int nral_rt(int ftype)
165 {
166     int nral;
167 
168     nral = NRAL(ftype);
169     if (interaction_function[ftype].flags & IF_VSITE)
170     {
171         /* With vsites the reverse topology contains an extra entry
172          * for storing if constructing atoms are vsites.
173          */
174         nral += 1;
175     }
176 
177     return nral;
178 }
179 
180 /*! \brief Return whether interactions of type \p ftype need to be assigned exactly once */
dd_check_ftype(int ftype,gmx_bool bBCheck,gmx_bool bConstr,gmx_bool bSettle)181 static gmx_bool dd_check_ftype(int ftype, gmx_bool bBCheck, gmx_bool bConstr, gmx_bool bSettle)
182 {
183     return ((((interaction_function[ftype].flags & IF_BOND) != 0U)
184              && ((interaction_function[ftype].flags & IF_VSITE) == 0U)
185              && (bBCheck || ((interaction_function[ftype].flags & IF_LIMZERO) == 0U)))
186             || (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) || (bSettle && ftype == F_SETTLE));
187 }
188 
189 /*! \brief Checks whether interactions have been assigned for one function type
190  *
191  * Loops over a list of interactions in the local topology of one function type
192  * and flags each of the interactions as assigned in the global \p isAssigned list.
193  * Exits with an inconsistency error when an interaction is assigned more than once.
194  */
flagInteractionsForType(const int ftype,const InteractionList & il,const reverse_ilist_t & ril,const gmx::Range<int> & atomRange,const int numAtomsPerMolecule,gmx::ArrayRef<const int> globalAtomIndices,gmx::ArrayRef<int> isAssigned)195 static void flagInteractionsForType(const int                ftype,
196                                     const InteractionList&   il,
197                                     const reverse_ilist_t&   ril,
198                                     const gmx::Range<int>&   atomRange,
199                                     const int                numAtomsPerMolecule,
200                                     gmx::ArrayRef<const int> globalAtomIndices,
201                                     gmx::ArrayRef<int>       isAssigned)
202 {
203     const int nril_mol = ril.index[numAtomsPerMolecule];
204     const int nral     = NRAL(ftype);
205 
206     for (int i = 0; i < il.size(); i += 1 + nral)
207     {
208         // ia[0] is the interaction type, ia[1, ...] the atom indices
209         const int* ia = il.iatoms.data() + i;
210         // Extract the global atom index of the first atom in this interaction
211         const int a0 = globalAtomIndices[ia[1]];
212         /* Check if this interaction is in
213          * the currently checked molblock.
214          */
215         if (atomRange.isInRange(a0))
216         {
217             // The molecule index in the list of this molecule type
218             const int moleculeIndex = (a0 - atomRange.begin()) / numAtomsPerMolecule;
219             const int atomOffset = (a0 - atomRange.begin()) - moleculeIndex * numAtomsPerMolecule;
220             const int globalAtomStartInMolecule = atomRange.begin() + moleculeIndex * numAtomsPerMolecule;
221             int       j_mol                     = ril.index[atomOffset];
222             bool found                          = false;
223             while (j_mol < ril.index[atomOffset + 1] && !found)
224             {
225                 const int j       = moleculeIndex * nril_mol + j_mol;
226                 const int ftype_j = ril.il[j_mol];
227                 /* Here we need to check if this interaction has
228                  * not already been assigned, since we could have
229                  * multiply defined interactions.
230                  */
231                 if (ftype == ftype_j && ia[0] == ril.il[j_mol + 1] && isAssigned[j] == 0)
232                 {
233                     /* Check the atoms */
234                     found = true;
235                     for (int a = 0; a < nral; a++)
236                     {
237                         if (globalAtomIndices[ia[1 + a]]
238                             != globalAtomStartInMolecule + ril.il[j_mol + 2 + a])
239                         {
240                             found = false;
241                         }
242                     }
243                     if (found)
244                     {
245                         isAssigned[j] = 1;
246                     }
247                 }
248                 j_mol += 2 + nral_rt(ftype_j);
249             }
250             if (!found)
251             {
252                 gmx_incons("Some interactions seem to be assigned multiple times");
253             }
254         }
255     }
256 }
257 
258 /*! \brief Help print error output when interactions are missing in a molblock
259  *
260  * \note This function needs to be called on all ranks (contains a global summation)
261  */
printMissingInteractionsMolblock(t_commrec * cr,const gmx_reverse_top_t & rt,const char * moltypename,const reverse_ilist_t & ril,const gmx::Range<int> & atomRange,const int numAtomsPerMolecule,const int numMolecules,const InteractionDefinitions & idef)262 static std::string printMissingInteractionsMolblock(t_commrec*               cr,
263                                                     const gmx_reverse_top_t& rt,
264                                                     const char*              moltypename,
265                                                     const reverse_ilist_t&   ril,
266                                                     const gmx::Range<int>&   atomRange,
267                                                     const int                numAtomsPerMolecule,
268                                                     const int                numMolecules,
269                                                     const InteractionDefinitions& idef)
270 {
271     const int               nril_mol = ril.index[numAtomsPerMolecule];
272     std::vector<int>        isAssigned(numMolecules * nril_mol, 0);
273     gmx::StringOutputStream stream;
274     gmx::TextWriter         log(&stream);
275 
276     for (int ftype = 0; ftype < F_NRE; ftype++)
277     {
278         if (dd_check_ftype(ftype, rt.bBCheck, rt.bConstr, rt.bSettle))
279         {
280             flagInteractionsForType(ftype, idef.il[ftype], ril, atomRange, numAtomsPerMolecule,
281                                     cr->dd->globalAtomIndices, isAssigned);
282         }
283     }
284 
285     gmx_sumi(isAssigned.size(), isAssigned.data(), cr);
286 
287     const int numMissingToPrint = 10;
288     int       i                 = 0;
289     for (int mol = 0; mol < numMolecules; mol++)
290     {
291         int j_mol = 0;
292         while (j_mol < nril_mol)
293         {
294             int ftype = ril.il[j_mol];
295             int nral  = NRAL(ftype);
296             int j     = mol * nril_mol + j_mol;
297             if (isAssigned[j] == 0 && !(interaction_function[ftype].flags & IF_VSITE))
298             {
299                 if (DDMASTER(cr->dd))
300                 {
301                     if (i == 0)
302                     {
303                         log.writeLineFormatted("Molecule type '%s'", moltypename);
304                         log.writeLineFormatted(
305                                 "the first %d missing interactions, except for exclusions:",
306                                 numMissingToPrint);
307                     }
308                     log.writeStringFormatted("%20s atoms", interaction_function[ftype].longname);
309                     int a;
310                     for (a = 0; a < nral; a++)
311                     {
312                         log.writeStringFormatted("%5d", ril.il[j_mol + 2 + a] + 1);
313                     }
314                     while (a < 4)
315                     {
316                         log.writeString("     ");
317                         a++;
318                     }
319                     log.writeString(" global");
320                     for (a = 0; a < nral; a++)
321                     {
322                         log.writeStringFormatted("%6d", atomRange.begin() + mol * numAtomsPerMolecule
323                                                                 + ril.il[j_mol + 2 + a] + 1);
324                     }
325                     log.ensureLineBreak();
326                 }
327                 i++;
328                 if (i >= numMissingToPrint)
329                 {
330                     break;
331                 }
332             }
333             j_mol += 2 + nral_rt(ftype);
334         }
335     }
336 
337     return stream.toString();
338 }
339 
340 /*! \brief Help print error output when interactions are missing */
printMissingInteractionsAtoms(const gmx::MDLogger & mdlog,t_commrec * cr,const gmx_mtop_t & mtop,const InteractionDefinitions & idef)341 static void printMissingInteractionsAtoms(const gmx::MDLogger&          mdlog,
342                                           t_commrec*                    cr,
343                                           const gmx_mtop_t&             mtop,
344                                           const InteractionDefinitions& idef)
345 {
346     const gmx_reverse_top_t& rt = *cr->dd->reverse_top;
347 
348     /* Print the atoms in the missing interactions per molblock */
349     int a_end = 0;
350     for (const gmx_molblock_t& molb : mtop.molblock)
351     {
352         const gmx_moltype_t& moltype = mtop.moltype[molb.type];
353         const int            a_start = a_end;
354         a_end                        = a_start + molb.nmol * moltype.atoms.nr;
355         const gmx::Range<int> atomRange(a_start, a_end);
356 
357         auto warning = printMissingInteractionsMolblock(cr, rt, *(moltype.name), rt.ril_mt[molb.type],
358                                                         atomRange, moltype.atoms.nr, molb.nmol, idef);
359 
360         GMX_LOG(mdlog.warning).appendText(warning);
361     }
362 }
363 
dd_print_missing_interactions(const gmx::MDLogger & mdlog,t_commrec * cr,int local_count,const gmx_mtop_t * top_global,const gmx_localtop_t * top_local,gmx::ArrayRef<const gmx::RVec> x,const matrix box)364 void dd_print_missing_interactions(const gmx::MDLogger&           mdlog,
365                                    t_commrec*                     cr,
366                                    int                            local_count,
367                                    const gmx_mtop_t*              top_global,
368                                    const gmx_localtop_t*          top_local,
369                                    gmx::ArrayRef<const gmx::RVec> x,
370                                    const matrix                   box)
371 {
372     int           ndiff_tot, cl[F_NRE], n, ndiff, rest_global, rest_local;
373     int           ftype, nral;
374     gmx_domdec_t* dd;
375 
376     dd = cr->dd;
377 
378     GMX_LOG(mdlog.warning)
379             .appendText(
380                     "Not all bonded interactions have been properly assigned to the domain "
381                     "decomposition cells");
382 
383     ndiff_tot = local_count - dd->nbonded_global;
384 
385     for (ftype = 0; ftype < F_NRE; ftype++)
386     {
387         nral      = NRAL(ftype);
388         cl[ftype] = top_local->idef.il[ftype].size() / (1 + nral);
389     }
390 
391     gmx_sumi(F_NRE, cl, cr);
392 
393     if (DDMASTER(dd))
394     {
395         GMX_LOG(mdlog.warning).appendText("A list of missing interactions:");
396         rest_global = dd->nbonded_global;
397         rest_local  = local_count;
398         for (ftype = 0; ftype < F_NRE; ftype++)
399         {
400             /* In the reverse and local top all constraints are merged
401              * into F_CONSTR. So in the if statement we skip F_CONSTRNC
402              * and add these constraints when doing F_CONSTR.
403              */
404             if (((interaction_function[ftype].flags & IF_BOND)
405                  && (dd->reverse_top->bBCheck || !(interaction_function[ftype].flags & IF_LIMZERO)))
406                 || (dd->reverse_top->bConstr && ftype == F_CONSTR)
407                 || (dd->reverse_top->bSettle && ftype == F_SETTLE))
408             {
409                 n = gmx_mtop_ftype_count(top_global, ftype);
410                 if (ftype == F_CONSTR)
411                 {
412                     n += gmx_mtop_ftype_count(top_global, F_CONSTRNC);
413                 }
414                 ndiff = cl[ftype] - n;
415                 if (ndiff != 0)
416                 {
417                     GMX_LOG(mdlog.warning)
418                             .appendTextFormatted("%20s of %6d missing %6d",
419                                                  interaction_function[ftype].longname, n, -ndiff);
420                 }
421                 rest_global -= n;
422                 rest_local -= cl[ftype];
423             }
424         }
425 
426         ndiff = rest_local - rest_global;
427         if (ndiff != 0)
428         {
429             GMX_LOG(mdlog.warning).appendTextFormatted("%20s of %6d missing %6d", "exclusions", rest_global, -ndiff);
430         }
431     }
432 
433     printMissingInteractionsAtoms(mdlog, cr, *top_global, top_local->idef);
434     write_dd_pdb("dd_dump_err", 0, "dump", top_global, cr, -1, as_rvec_array(x.data()), box);
435 
436     std::string errorMessage;
437 
438     if (ndiff_tot > 0)
439     {
440         errorMessage =
441                 "One or more interactions were assigned to multiple domains of the domain "
442                 "decompostion. Please report this bug.";
443     }
444     else
445     {
446         errorMessage = gmx::formatString(
447                 "%d of the %d bonded interactions could not be calculated because some atoms "
448                 "involved moved further apart than the multi-body cut-off distance (%g nm) or the "
449                 "two-body cut-off distance (%g nm), see option -rdd, for pairs and tabulated bonds "
450                 "also see option -ddcheck",
451                 -ndiff_tot, cr->dd->nbonded_global, dd_cutoff_multibody(dd), dd_cutoff_twobody(dd));
452     }
453     gmx_fatal_collective(FARGS, cr->mpi_comm_mygroup, MASTER(cr), "%s", errorMessage.c_str());
454 }
455 
456 /*! \brief Return global topology molecule information for global atom index \p i_gl */
global_atomnr_to_moltype_ind(const gmx_reverse_top_t * rt,int i_gl,int * mb,int * mt,int * mol,int * i_mol)457 static void global_atomnr_to_moltype_ind(const gmx_reverse_top_t* rt, int i_gl, int* mb, int* mt, int* mol, int* i_mol)
458 {
459     const MolblockIndices* mbi   = rt->mbi.data();
460     int                    start = 0;
461     int                    end   = rt->mbi.size(); /* exclusive */
462     int                    mid;
463 
464     /* binary search for molblock_ind */
465     while (TRUE)
466     {
467         mid = (start + end) / 2;
468         if (i_gl >= mbi[mid].a_end)
469         {
470             start = mid + 1;
471         }
472         else if (i_gl < mbi[mid].a_start)
473         {
474             end = mid;
475         }
476         else
477         {
478             break;
479         }
480     }
481 
482     *mb = mid;
483     mbi += mid;
484 
485     *mt    = mbi->type;
486     *mol   = (i_gl - mbi->a_start) / mbi->natoms_mol;
487     *i_mol = (i_gl - mbi->a_start) - (*mol) * mbi->natoms_mol;
488 }
489 
490 /*! \brief Returns the maximum number of exclusions per atom */
getMaxNumExclusionsPerAtom(const ListOfLists<int> & excls)491 static int getMaxNumExclusionsPerAtom(const ListOfLists<int>& excls)
492 {
493     int maxNumExcls = 0;
494     for (gmx::index at = 0; at < excls.ssize(); at++)
495     {
496         const auto list     = excls[at];
497         const int  numExcls = list.ssize();
498 
499         GMX_RELEASE_ASSERT(numExcls != 1 || list[0] == at,
500                            "With 1 exclusion we expect a self-exclusion");
501 
502         maxNumExcls = std::max(maxNumExcls, numExcls);
503     }
504 
505     return maxNumExcls;
506 }
507 
508 /*! \brief Run the reverse ilist generation and store it in r_il when \p bAssign = TRUE */
low_make_reverse_ilist(const InteractionLists & il_mt,const t_atom * atom,int * count,gmx_bool bConstr,gmx_bool bSettle,gmx_bool bBCheck,gmx::ArrayRef<const int> r_index,gmx::ArrayRef<int> r_il,gmx_bool bLinkToAllAtoms,gmx_bool bAssign)509 static int low_make_reverse_ilist(const InteractionLists&  il_mt,
510                                   const t_atom*            atom,
511                                   int*                     count,
512                                   gmx_bool                 bConstr,
513                                   gmx_bool                 bSettle,
514                                   gmx_bool                 bBCheck,
515                                   gmx::ArrayRef<const int> r_index,
516                                   gmx::ArrayRef<int>       r_il,
517                                   gmx_bool                 bLinkToAllAtoms,
518                                   gmx_bool                 bAssign)
519 {
520     int ftype, j, nlink, link;
521     int a;
522     int nint;
523 
524     nint = 0;
525     for (ftype = 0; ftype < F_NRE; ftype++)
526     {
527         if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE))
528             || (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) || (bSettle && ftype == F_SETTLE))
529         {
530             const bool  bVSite = ((interaction_function[ftype].flags & IF_VSITE) != 0U);
531             const int   nral   = NRAL(ftype);
532             const auto& il     = il_mt[ftype];
533             for (int i = 0; i < il.size(); i += 1 + nral)
534             {
535                 const int* ia = il.iatoms.data() + i;
536                 if (bLinkToAllAtoms)
537                 {
538                     if (bVSite)
539                     {
540                         /* We don't need the virtual sites for the cg-links */
541                         nlink = 0;
542                     }
543                     else
544                     {
545                         nlink = nral;
546                     }
547                 }
548                 else
549                 {
550                     /* Couple to the first atom in the interaction */
551                     nlink = 1;
552                 }
553                 for (link = 0; link < nlink; link++)
554                 {
555                     a = ia[1 + link];
556                     if (bAssign)
557                     {
558                         GMX_ASSERT(!r_il.empty(), "with bAssign not allowed to be empty");
559                         GMX_ASSERT(!r_index.empty(), "with bAssign not allowed to be empty");
560                         r_il[r_index[a] + count[a]]     = (ftype == F_CONSTRNC ? F_CONSTR : ftype);
561                         r_il[r_index[a] + count[a] + 1] = ia[0];
562                         for (j = 1; j < 1 + nral; j++)
563                         {
564                             /* Store the molecular atom number */
565                             r_il[r_index[a] + count[a] + 1 + j] = ia[j];
566                         }
567                     }
568                     if (interaction_function[ftype].flags & IF_VSITE)
569                     {
570                         if (bAssign)
571                         {
572                             /* Add an entry to iatoms for storing
573                              * which of the constructing atoms are
574                              * vsites again.
575                              */
576                             r_il[r_index[a] + count[a] + 2 + nral] = 0;
577                             for (j = 2; j < 1 + nral; j++)
578                             {
579                                 if (atom[ia[j]].ptype == eptVSite)
580                                 {
581                                     r_il[r_index[a] + count[a] + 2 + nral] |= (2 << j);
582                                 }
583                             }
584                         }
585                     }
586                     else
587                     {
588                         /* We do not count vsites since they are always
589                          * uniquely assigned and can be assigned
590                          * to multiple nodes with recursive vsites.
591                          */
592                         if (bBCheck || !(interaction_function[ftype].flags & IF_LIMZERO))
593                         {
594                             nint++;
595                         }
596                     }
597                     count[a] += 2 + nral_rt(ftype);
598                 }
599             }
600         }
601     }
602 
603     return nint;
604 }
605 
606 /*! \brief Make the reverse ilist: a list of bonded interactions linked to atoms */
make_reverse_ilist(const InteractionLists & ilist,const t_atoms * atoms,gmx_bool bConstr,gmx_bool bSettle,gmx_bool bBCheck,gmx_bool bLinkToAllAtoms,reverse_ilist_t * ril_mt)607 static int make_reverse_ilist(const InteractionLists& ilist,
608                               const t_atoms*          atoms,
609                               gmx_bool                bConstr,
610                               gmx_bool                bSettle,
611                               gmx_bool                bBCheck,
612                               gmx_bool                bLinkToAllAtoms,
613                               reverse_ilist_t*        ril_mt)
614 {
615     int nat_mt, *count, i, nint_mt;
616 
617     /* Count the interactions */
618     nat_mt = atoms->nr;
619     snew(count, nat_mt);
620     low_make_reverse_ilist(ilist, atoms->atom, count, bConstr, bSettle, bBCheck, {}, {},
621                            bLinkToAllAtoms, FALSE);
622 
623     ril_mt->index.push_back(0);
624     for (i = 0; i < nat_mt; i++)
625     {
626         ril_mt->index.push_back(ril_mt->index[i] + count[i]);
627         count[i] = 0;
628     }
629     ril_mt->il.resize(ril_mt->index[nat_mt]);
630 
631     /* Store the interactions */
632     nint_mt = low_make_reverse_ilist(ilist, atoms->atom, count, bConstr, bSettle, bBCheck,
633                                      ril_mt->index, ril_mt->il, bLinkToAllAtoms, TRUE);
634 
635     sfree(count);
636 
637     ril_mt->numAtomsInMolecule = atoms->nr;
638 
639     return nint_mt;
640 }
641 
642 /*! \brief Generate the reverse topology */
make_reverse_top(const gmx_mtop_t * mtop,gmx_bool bFE,gmx_bool bConstr,gmx_bool bSettle,gmx_bool bBCheck,int * nint)643 static gmx_reverse_top_t make_reverse_top(const gmx_mtop_t* mtop,
644                                           gmx_bool          bFE,
645                                           gmx_bool          bConstr,
646                                           gmx_bool          bSettle,
647                                           gmx_bool          bBCheck,
648                                           int*              nint)
649 {
650     gmx_reverse_top_t rt;
651 
652     /* Should we include constraints (for SHAKE) in rt? */
653     rt.bConstr = bConstr;
654     rt.bSettle = bSettle;
655     rt.bBCheck = bBCheck;
656 
657     rt.bInterAtomicInteractions = mtop->bIntermolecularInteractions;
658     rt.ril_mt.resize(mtop->moltype.size());
659     rt.ril_mt_tot_size = 0;
660     std::vector<int> nint_mt;
661     for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
662     {
663         const gmx_moltype_t& molt = mtop->moltype[mt];
664         if (molt.atoms.nr > 1)
665         {
666             rt.bInterAtomicInteractions = true;
667         }
668 
669         /* Make the atom to interaction list for this molecule type */
670         int numberOfInteractions = make_reverse_ilist(
671                 molt.ilist, &molt.atoms, rt.bConstr, rt.bSettle, rt.bBCheck, FALSE, &rt.ril_mt[mt]);
672         nint_mt.push_back(numberOfInteractions);
673 
674         rt.ril_mt_tot_size += rt.ril_mt[mt].index[molt.atoms.nr];
675     }
676     if (debug)
677     {
678         fprintf(debug, "The total size of the atom to interaction index is %d integers\n",
679                 rt.ril_mt_tot_size);
680     }
681 
682     *nint = 0;
683     for (const gmx_molblock_t& molblock : mtop->molblock)
684     {
685         *nint += molblock.nmol * nint_mt[molblock.type];
686     }
687 
688     /* Make an intermolecular reverse top, if necessary */
689     rt.bIntermolecularInteractions = mtop->bIntermolecularInteractions;
690     if (rt.bIntermolecularInteractions)
691     {
692         t_atoms atoms_global;
693 
694         atoms_global.nr   = mtop->natoms;
695         atoms_global.atom = nullptr; /* Only used with virtual sites */
696 
697         GMX_RELEASE_ASSERT(mtop->intermolecular_ilist,
698                            "We should have an ilist when intermolecular interactions are on");
699 
700         *nint += make_reverse_ilist(*mtop->intermolecular_ilist, &atoms_global, rt.bConstr,
701                                     rt.bSettle, rt.bBCheck, FALSE, &rt.ril_intermol);
702     }
703 
704     if (bFE && gmx_mtop_bondeds_free_energy(mtop))
705     {
706         rt.ilsort = ilsortFE_UNSORTED;
707     }
708     else
709     {
710         rt.ilsort = ilsortNO_FE;
711     }
712 
713     /* Make a molblock index for fast searching */
714     int i = 0;
715     for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
716     {
717         const gmx_molblock_t& molb           = mtop->molblock[mb];
718         const int             numAtomsPerMol = mtop->moltype[molb.type].atoms.nr;
719         MolblockIndices       mbi;
720         mbi.a_start = i;
721         i += molb.nmol * numAtomsPerMol;
722         mbi.a_end      = i;
723         mbi.natoms_mol = numAtomsPerMol;
724         mbi.type       = molb.type;
725         rt.mbi.push_back(mbi);
726     }
727 
728     for (int th = 0; th < gmx_omp_nthreads_get(emntDomdec); th++)
729     {
730         rt.th_work.emplace_back(mtop->ffparams);
731     }
732 
733     return rt;
734 }
735 
dd_make_reverse_top(FILE * fplog,gmx_domdec_t * dd,const gmx_mtop_t * mtop,const gmx::VirtualSitesHandler * vsite,const t_inputrec * ir,gmx_bool bBCheck)736 void dd_make_reverse_top(FILE*                           fplog,
737                          gmx_domdec_t*                   dd,
738                          const gmx_mtop_t*               mtop,
739                          const gmx::VirtualSitesHandler* vsite,
740                          const t_inputrec*               ir,
741                          gmx_bool                        bBCheck)
742 {
743     if (fplog)
744     {
745         fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
746     }
747 
748     /* If normal and/or settle constraints act only within charge groups,
749      * we can store them in the reverse top and simply assign them to domains.
750      * Otherwise we need to assign them to multiple domains and set up
751      * the parallel version constraint algorithm(s).
752      */
753 
754     dd->reverse_top = new gmx_reverse_top_t;
755     *dd->reverse_top =
756             make_reverse_top(mtop, ir->efep != efepNO, !dd->comm->systemInfo.haveSplitConstraints,
757                              !dd->comm->systemInfo.haveSplitSettles, bBCheck, &dd->nbonded_global);
758 
759     dd->haveExclusions = false;
760     for (const gmx_molblock_t& molb : mtop->molblock)
761     {
762         const int maxNumExclusionsPerAtom = getMaxNumExclusionsPerAtom(mtop->moltype[molb.type].excls);
763         // We checked above that max 1 exclusion means only self exclusions
764         if (maxNumExclusionsPerAtom > 1)
765         {
766             dd->haveExclusions = true;
767         }
768     }
769 
770     if (vsite && vsite->numInterUpdategroupVirtualSites() > 0)
771     {
772         if (fplog)
773         {
774             fprintf(fplog,
775                     "There are %d inter update-group virtual sites,\n"
776                     "will an extra communication step for selected coordinates and forces\n",
777                     vsite->numInterUpdategroupVirtualSites());
778         }
779         init_domdec_vsites(dd, vsite->numInterUpdategroupVirtualSites());
780     }
781 
782     if (dd->comm->systemInfo.haveSplitConstraints || dd->comm->systemInfo.haveSplitSettles)
783     {
784         init_domdec_constraints(dd, mtop);
785     }
786     if (fplog)
787     {
788         fprintf(fplog, "\n");
789     }
790 }
791 
792 /*! \brief Store a vsite interaction at the end of \p il
793  *
794  * This routine is very similar to add_ifunc, but vsites interactions
795  * have more work to do than other kinds of interactions, and the
796  * complex way nral (and thus vector contents) depends on ftype
797  * confuses static analysis tools unless we fuse the vsite
798  * atom-indexing organization code with the ifunc-adding code, so that
799  * they can see that nral is the same value. */
add_ifunc_for_vsites(t_iatom * tiatoms,const gmx_ga2la_t & ga2la,int nral,gmx_bool bHomeA,int a,int a_gl,int a_mol,const t_iatom * iatoms,InteractionList * il)800 static inline void add_ifunc_for_vsites(t_iatom*           tiatoms,
801                                         const gmx_ga2la_t& ga2la,
802                                         int                nral,
803                                         gmx_bool           bHomeA,
804                                         int                a,
805                                         int                a_gl,
806                                         int                a_mol,
807                                         const t_iatom*     iatoms,
808                                         InteractionList*   il)
809 {
810     /* Copy the type */
811     tiatoms[0] = iatoms[0];
812 
813     if (bHomeA)
814     {
815         /* We know the local index of the first atom */
816         tiatoms[1] = a;
817     }
818     else
819     {
820         /* Convert later in make_local_vsites */
821         tiatoms[1] = -a_gl - 1;
822     }
823 
824     GMX_ASSERT(nral >= 2 && nral <= 5, "Invalid nral for vsites");
825     for (int k = 2; k < 1 + nral; k++)
826     {
827         int ak_gl = a_gl + iatoms[k] - a_mol;
828         if (const int* homeIndex = ga2la.findHome(ak_gl))
829         {
830             tiatoms[k] = *homeIndex;
831         }
832         else
833         {
834             /* Copy the global index, convert later in make_local_vsites */
835             tiatoms[k] = -(ak_gl + 1);
836         }
837         // Note that ga2la_get_home always sets the third parameter if
838         // it returns TRUE
839     }
840     il->push_back(tiatoms[0], nral, tiatoms + 1);
841 }
842 
843 /*! \brief Store a position restraint in idef and iatoms, complex because the parameters are different for each entry */
add_posres(int mol,int a_mol,int numAtomsInMolecule,const gmx_molblock_t * molb,t_iatom * iatoms,const t_iparams * ip_in,InteractionDefinitions * idef)844 static void add_posres(int                     mol,
845                        int                     a_mol,
846                        int                     numAtomsInMolecule,
847                        const gmx_molblock_t*   molb,
848                        t_iatom*                iatoms,
849                        const t_iparams*        ip_in,
850                        InteractionDefinitions* idef)
851 {
852     /* This position restraint has not been added yet,
853      * so it's index is the current number of position restraints.
854      */
855     const int n = idef->il[F_POSRES].size() / 2;
856 
857     /* Get the position restraint coordinates from the molblock */
858     const int a_molb = mol * numAtomsInMolecule + a_mol;
859     GMX_ASSERT(a_molb < ssize(molb->posres_xA),
860                "We need a sufficient number of position restraint coordinates");
861     /* Copy the force constants */
862     t_iparams ip        = ip_in[iatoms[0]];
863     ip.posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
864     ip.posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
865     ip.posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
866     if (!molb->posres_xB.empty())
867     {
868         ip.posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
869         ip.posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
870         ip.posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
871     }
872     else
873     {
874         ip.posres.pos0B[XX] = ip.posres.pos0A[XX];
875         ip.posres.pos0B[YY] = ip.posres.pos0A[YY];
876         ip.posres.pos0B[ZZ] = ip.posres.pos0A[ZZ];
877     }
878     /* Set the parameter index for idef->iparams_posres */
879     iatoms[0] = n;
880     idef->iparams_posres.push_back(ip);
881 
882     GMX_ASSERT(int(idef->iparams_posres.size()) == n + 1,
883                "The index of the parameter type should match n");
884 }
885 
886 /*! \brief Store a flat-bottomed position restraint in idef and iatoms, complex because the parameters are different for each entry */
add_fbposres(int mol,int a_mol,int numAtomsInMolecule,const gmx_molblock_t * molb,t_iatom * iatoms,const t_iparams * ip_in,InteractionDefinitions * idef)887 static void add_fbposres(int                     mol,
888                          int                     a_mol,
889                          int                     numAtomsInMolecule,
890                          const gmx_molblock_t*   molb,
891                          t_iatom*                iatoms,
892                          const t_iparams*        ip_in,
893                          InteractionDefinitions* idef)
894 {
895     /* This flat-bottom position restraint has not been added yet,
896      * so it's index is the current number of position restraints.
897      */
898     const int n = idef->il[F_FBPOSRES].size() / 2;
899 
900     /* Get the position restraint coordinats from the molblock */
901     const int a_molb = mol * numAtomsInMolecule + a_mol;
902     GMX_ASSERT(a_molb < ssize(molb->posres_xA),
903                "We need a sufficient number of position restraint coordinates");
904     /* Copy the force constants */
905     t_iparams ip = ip_in[iatoms[0]];
906     /* Take reference positions from A position of normal posres */
907     ip.fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
908     ip.fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
909     ip.fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
910 
911     /* Note: no B-type for flat-bottom posres */
912 
913     /* Set the parameter index for idef->iparams_fbposres */
914     iatoms[0] = n;
915     idef->iparams_fbposres.push_back(ip);
916 
917     GMX_ASSERT(int(idef->iparams_fbposres.size()) == n + 1,
918                "The index of the parameter type should match n");
919 }
920 
921 /*! \brief Store a virtual site interaction, complex because of PBC and recursion */
add_vsite(const gmx_ga2la_t & ga2la,gmx::ArrayRef<const int> index,gmx::ArrayRef<const int> rtil,int ftype,int nral,gmx_bool bHomeA,int a,int a_gl,int a_mol,const t_iatom * iatoms,InteractionDefinitions * idef)922 static void add_vsite(const gmx_ga2la_t&       ga2la,
923                       gmx::ArrayRef<const int> index,
924                       gmx::ArrayRef<const int> rtil,
925                       int                      ftype,
926                       int                      nral,
927                       gmx_bool                 bHomeA,
928                       int                      a,
929                       int                      a_gl,
930                       int                      a_mol,
931                       const t_iatom*           iatoms,
932                       InteractionDefinitions*  idef)
933 {
934     int     k;
935     t_iatom tiatoms[1 + MAXATOMLIST];
936     int     j, ftype_r, nral_r;
937 
938     /* Add this interaction to the local topology */
939     add_ifunc_for_vsites(tiatoms, ga2la, nral, bHomeA, a, a_gl, a_mol, iatoms, &idef->il[ftype]);
940 
941     if (iatoms[1 + nral])
942     {
943         /* Check for recursion */
944         for (k = 2; k < 1 + nral; k++)
945         {
946             if ((iatoms[1 + nral] & (2 << k)) && (tiatoms[k] < 0))
947             {
948                 /* This construction atoms is a vsite and not a home atom */
949                 if (gmx_debug_at)
950                 {
951                     fprintf(debug, "Constructing atom %d of vsite atom %d is a vsite and non-home\n",
952                             iatoms[k] + 1, a_mol + 1);
953                 }
954                 /* Find the vsite construction */
955 
956                 /* Check all interactions assigned to this atom */
957                 j = index[iatoms[k]];
958                 while (j < index[iatoms[k] + 1])
959                 {
960                     ftype_r = rtil[j++];
961                     nral_r  = NRAL(ftype_r);
962                     if (interaction_function[ftype_r].flags & IF_VSITE)
963                     {
964                         /* Add this vsite (recursion) */
965                         add_vsite(ga2la, index, rtil, ftype_r, nral_r, FALSE, -1,
966                                   a_gl + iatoms[k] - iatoms[1], iatoms[k], rtil.data() + j, idef);
967                     }
968                     j += 1 + nral_rt(ftype_r);
969                 }
970             }
971         }
972     }
973 }
974 
975 /*! \brief Returns the squared distance between atoms \p i and \p j */
dd_dist2(t_pbc * pbc_null,const rvec * x,const int i,int j)976 static real dd_dist2(t_pbc* pbc_null, const rvec* x, const int i, int j)
977 {
978     rvec dx;
979 
980     if (pbc_null)
981     {
982         pbc_dx_aiuc(pbc_null, x[i], x[j], dx);
983     }
984     else
985     {
986         rvec_sub(x[i], x[j], dx);
987     }
988 
989     return norm2(dx);
990 }
991 
992 /*! \brief Append t_idef structures 1 to nsrc in src to *dest */
combine_idef(InteractionDefinitions * dest,gmx::ArrayRef<const thread_work_t> src)993 static void combine_idef(InteractionDefinitions* dest, gmx::ArrayRef<const thread_work_t> src)
994 {
995     int ftype;
996 
997     for (ftype = 0; ftype < F_NRE; ftype++)
998     {
999         int n = 0;
1000         for (gmx::index s = 1; s < src.ssize(); s++)
1001         {
1002             n += src[s].idef.il[ftype].size();
1003         }
1004         if (n > 0)
1005         {
1006             for (gmx::index s = 1; s < src.ssize(); s++)
1007             {
1008                 dest->il[ftype].append(src[s].idef.il[ftype]);
1009             }
1010 
1011             /* Position restraints need an additional treatment */
1012             if (ftype == F_POSRES || ftype == F_FBPOSRES)
1013             {
1014                 int                     nposres = dest->il[ftype].size() / 2;
1015                 std::vector<t_iparams>& iparams_dest =
1016                         (ftype == F_POSRES ? dest->iparams_posres : dest->iparams_fbposres);
1017 
1018                 /* Set nposres to the number of original position restraints in dest */
1019                 for (gmx::index s = 1; s < src.ssize(); s++)
1020                 {
1021                     nposres -= src[s].idef.il[ftype].size() / 2;
1022                 }
1023 
1024                 for (gmx::index s = 1; s < src.ssize(); s++)
1025                 {
1026                     const std::vector<t_iparams>& iparams_src =
1027                             (ftype == F_POSRES ? src[s].idef.iparams_posres : src[s].idef.iparams_fbposres);
1028                     iparams_dest.insert(iparams_dest.end(), iparams_src.begin(), iparams_src.end());
1029 
1030                     /* Correct the indices into iparams_posres */
1031                     for (int i = 0; i < src[s].idef.il[ftype].size() / 2; i++)
1032                     {
1033                         /* Correct the index into iparams_posres */
1034                         dest->il[ftype].iatoms[nposres * 2] = nposres;
1035                         nposres++;
1036                     }
1037                 }
1038                 GMX_RELEASE_ASSERT(
1039                         int(iparams_dest.size()) == nposres,
1040                         "The number of parameters should match the number of restraints");
1041             }
1042         }
1043     }
1044 }
1045 
1046 /*! \brief Check and when available assign bonded interactions for local atom i
1047  */
check_assign_interactions_atom(int i,int i_gl,int mol,int i_mol,int numAtomsInMolecule,gmx::ArrayRef<const int> index,gmx::ArrayRef<const int> rtil,gmx_bool bInterMolInteractions,int ind_start,int ind_end,const gmx_domdec_t * dd,const gmx_domdec_zones_t * zones,const gmx_molblock_t * molb,gmx_bool bRCheckMB,const ivec rcheck,gmx_bool bRCheck2B,real rc2,t_pbc * pbc_null,rvec * cg_cm,const t_iparams * ip_in,InteractionDefinitions * idef,int iz,gmx_bool bBCheck,int * nbonded_local)1048 static inline void check_assign_interactions_atom(int                       i,
1049                                                   int                       i_gl,
1050                                                   int                       mol,
1051                                                   int                       i_mol,
1052                                                   int                       numAtomsInMolecule,
1053                                                   gmx::ArrayRef<const int>  index,
1054                                                   gmx::ArrayRef<const int>  rtil,
1055                                                   gmx_bool                  bInterMolInteractions,
1056                                                   int                       ind_start,
1057                                                   int                       ind_end,
1058                                                   const gmx_domdec_t*       dd,
1059                                                   const gmx_domdec_zones_t* zones,
1060                                                   const gmx_molblock_t*     molb,
1061                                                   gmx_bool                  bRCheckMB,
1062                                                   const ivec                rcheck,
1063                                                   gmx_bool                  bRCheck2B,
1064                                                   real                      rc2,
1065                                                   t_pbc*                    pbc_null,
1066                                                   rvec*                     cg_cm,
1067                                                   const t_iparams*          ip_in,
1068                                                   InteractionDefinitions*   idef,
1069                                                   int                       iz,
1070                                                   gmx_bool                  bBCheck,
1071                                                   int*                      nbonded_local)
1072 {
1073     gmx::ArrayRef<const DDPairInteractionRanges> iZones = zones->iZones;
1074 
1075     int j = ind_start;
1076     while (j < ind_end)
1077     {
1078         t_iatom tiatoms[1 + MAXATOMLIST];
1079 
1080         const int ftype  = rtil[j++];
1081         auto      iatoms = gmx::constArrayRefFromArray(rtil.data() + j, rtil.size() - j);
1082         const int nral   = NRAL(ftype);
1083         if (interaction_function[ftype].flags & IF_VSITE)
1084         {
1085             assert(!bInterMolInteractions);
1086             /* The vsite construction goes where the vsite itself is */
1087             if (iz == 0)
1088             {
1089                 add_vsite(*dd->ga2la, index, rtil, ftype, nral, TRUE, i, i_gl, i_mol, iatoms.data(), idef);
1090             }
1091         }
1092         else
1093         {
1094             gmx_bool bUse;
1095 
1096             /* Copy the type */
1097             tiatoms[0] = iatoms[0];
1098 
1099             if (nral == 1)
1100             {
1101                 assert(!bInterMolInteractions);
1102                 /* Assign single-body interactions to the home zone */
1103                 if (iz == 0)
1104                 {
1105                     bUse       = TRUE;
1106                     tiatoms[1] = i;
1107                     if (ftype == F_POSRES)
1108                     {
1109                         add_posres(mol, i_mol, numAtomsInMolecule, molb, tiatoms, ip_in, idef);
1110                     }
1111                     else if (ftype == F_FBPOSRES)
1112                     {
1113                         add_fbposres(mol, i_mol, numAtomsInMolecule, molb, tiatoms, ip_in, idef);
1114                     }
1115                 }
1116                 else
1117                 {
1118                     bUse = FALSE;
1119                 }
1120             }
1121             else if (nral == 2)
1122             {
1123                 /* This is a two-body interaction, we can assign
1124                  * analogous to the non-bonded assignments.
1125                  */
1126                 int k_gl;
1127 
1128                 if (!bInterMolInteractions)
1129                 {
1130                     /* Get the global index using the offset in the molecule */
1131                     k_gl = i_gl + iatoms[2] - i_mol;
1132                 }
1133                 else
1134                 {
1135                     k_gl = iatoms[2];
1136                 }
1137                 if (const auto* entry = dd->ga2la->find(k_gl))
1138                 {
1139                     int kz = entry->cell;
1140                     if (kz >= zones->n)
1141                     {
1142                         kz -= zones->n;
1143                     }
1144                     /* Check zone interaction assignments */
1145                     bUse = ((iz < iZones.ssize() && iz <= kz && iZones[iz].jZoneRange.isInRange(kz))
1146                             || (kz < iZones.ssize() && iz > kz && iZones[kz].jZoneRange.isInRange(iz)));
1147                     if (bUse)
1148                     {
1149                         GMX_ASSERT(ftype != F_CONSTR || (iz == 0 && kz == 0),
1150                                    "Constraint assigned here should only involve home atoms");
1151 
1152                         tiatoms[1] = i;
1153                         tiatoms[2] = entry->la;
1154                         /* If necessary check the cgcm distance */
1155                         if (bRCheck2B && dd_dist2(pbc_null, cg_cm, tiatoms[1], tiatoms[2]) >= rc2)
1156                         {
1157                             bUse = FALSE;
1158                         }
1159                     }
1160                 }
1161                 else
1162                 {
1163                     bUse = false;
1164                 }
1165             }
1166             else
1167             {
1168                 /* Assign this multi-body bonded interaction to
1169                  * the local node if we have all the atoms involved
1170                  * (local or communicated) and the minimum zone shift
1171                  * in each dimension is zero, for dimensions
1172                  * with 2 DD cells an extra check may be necessary.
1173                  */
1174                 ivec k_zero, k_plus;
1175                 int  k;
1176 
1177                 bUse = TRUE;
1178                 clear_ivec(k_zero);
1179                 clear_ivec(k_plus);
1180                 for (k = 1; k <= nral && bUse; k++)
1181                 {
1182                     int k_gl;
1183                     if (!bInterMolInteractions)
1184                     {
1185                         /* Get the global index using the offset in the molecule */
1186                         k_gl = i_gl + iatoms[k] - i_mol;
1187                     }
1188                     else
1189                     {
1190                         k_gl = iatoms[k];
1191                     }
1192                     const auto* entry = dd->ga2la->find(k_gl);
1193                     if (entry == nullptr || entry->cell >= zones->n)
1194                     {
1195                         /* We do not have this atom of this interaction
1196                          * locally, or it comes from more than one cell
1197                          * away.
1198                          */
1199                         bUse = FALSE;
1200                     }
1201                     else
1202                     {
1203                         int d;
1204 
1205                         tiatoms[k] = entry->la;
1206                         for (d = 0; d < DIM; d++)
1207                         {
1208                             if (zones->shift[entry->cell][d] == 0)
1209                             {
1210                                 k_zero[d] = k;
1211                             }
1212                             else
1213                             {
1214                                 k_plus[d] = k;
1215                             }
1216                         }
1217                     }
1218                 }
1219                 bUse = (bUse && (k_zero[XX] != 0) && (k_zero[YY] != 0) && (k_zero[ZZ] != 0));
1220                 if (bRCheckMB)
1221                 {
1222                     int d;
1223 
1224                     for (d = 0; (d < DIM && bUse); d++)
1225                     {
1226                         /* Check if the cg_cm distance falls within
1227                          * the cut-off to avoid possible multiple
1228                          * assignments of bonded interactions.
1229                          */
1230                         if (rcheck[d] && k_plus[d]
1231                             && dd_dist2(pbc_null, cg_cm, tiatoms[k_zero[d]], tiatoms[k_plus[d]]) >= rc2)
1232                         {
1233                             bUse = FALSE;
1234                         }
1235                     }
1236                 }
1237             }
1238             if (bUse)
1239             {
1240                 /* Add this interaction to the local topology */
1241                 idef->il[ftype].push_back(tiatoms[0], nral, tiatoms + 1);
1242                 /* Sum so we can check in global_stat
1243                  * if we have everything.
1244                  */
1245                 if (bBCheck || !(interaction_function[ftype].flags & IF_LIMZERO))
1246                 {
1247                     (*nbonded_local)++;
1248                 }
1249             }
1250         }
1251         j += 1 + nral_rt(ftype);
1252     }
1253 }
1254 
1255 /*! \brief This function looks up and assigns bonded interactions for zone iz.
1256  *
1257  * With thread parallelizing each thread acts on a different atom range:
1258  * at_start to at_end.
1259  */
make_bondeds_zone(gmx_domdec_t * dd,const gmx_domdec_zones_t * zones,const std::vector<gmx_molblock_t> & molb,gmx_bool bRCheckMB,ivec rcheck,gmx_bool bRCheck2B,real rc2,t_pbc * pbc_null,rvec * cg_cm,const t_iparams * ip_in,InteractionDefinitions * idef,int izone,const gmx::Range<int> & atomRange)1260 static int make_bondeds_zone(gmx_domdec_t*                      dd,
1261                              const gmx_domdec_zones_t*          zones,
1262                              const std::vector<gmx_molblock_t>& molb,
1263                              gmx_bool                           bRCheckMB,
1264                              ivec                               rcheck,
1265                              gmx_bool                           bRCheck2B,
1266                              real                               rc2,
1267                              t_pbc*                             pbc_null,
1268                              rvec*                              cg_cm,
1269                              const t_iparams*                   ip_in,
1270                              InteractionDefinitions*            idef,
1271                              int                                izone,
1272                              const gmx::Range<int>&             atomRange)
1273 {
1274     int                mb, mt, mol, i_mol;
1275     gmx_bool           bBCheck;
1276     gmx_reverse_top_t* rt;
1277     int                nbonded_local;
1278 
1279     rt = dd->reverse_top;
1280 
1281     bBCheck = rt->bBCheck;
1282 
1283     nbonded_local = 0;
1284 
1285     for (int i : atomRange)
1286     {
1287         /* Get the global atom number */
1288         const int i_gl = dd->globalAtomIndices[i];
1289         global_atomnr_to_moltype_ind(rt, i_gl, &mb, &mt, &mol, &i_mol);
1290         /* Check all intramolecular interactions assigned to this atom */
1291         gmx::ArrayRef<const int>     index = rt->ril_mt[mt].index;
1292         gmx::ArrayRef<const t_iatom> rtil  = rt->ril_mt[mt].il;
1293 
1294         check_assign_interactions_atom(i, i_gl, mol, i_mol, rt->ril_mt[mt].numAtomsInMolecule,
1295                                        index, rtil, FALSE, index[i_mol], index[i_mol + 1], dd,
1296                                        zones, &molb[mb], bRCheckMB, rcheck, bRCheck2B, rc2,
1297                                        pbc_null, cg_cm, ip_in, idef, izone, bBCheck, &nbonded_local);
1298 
1299 
1300         if (rt->bIntermolecularInteractions)
1301         {
1302             /* Check all intermolecular interactions assigned to this atom */
1303             index = rt->ril_intermol.index;
1304             rtil  = rt->ril_intermol.il;
1305 
1306             check_assign_interactions_atom(i, i_gl, mol, i_mol, rt->ril_mt[mt].numAtomsInMolecule,
1307                                            index, rtil, TRUE, index[i_gl], index[i_gl + 1], dd, zones,
1308                                            &molb[mb], bRCheckMB, rcheck, bRCheck2B, rc2, pbc_null,
1309                                            cg_cm, ip_in, idef, izone, bBCheck, &nbonded_local);
1310         }
1311     }
1312 
1313     return nbonded_local;
1314 }
1315 
1316 /*! \brief Set the exclusion data for i-zone \p iz */
make_exclusions_zone(gmx_domdec_t * dd,gmx_domdec_zones_t * zones,const std::vector<gmx_moltype_t> & moltype,const int * cginfo,ListOfLists<int> * lexcls,int iz,int at_start,int at_end,const gmx::ArrayRef<const int> intermolecularExclusionGroup)1317 static void make_exclusions_zone(gmx_domdec_t*                     dd,
1318                                  gmx_domdec_zones_t*               zones,
1319                                  const std::vector<gmx_moltype_t>& moltype,
1320                                  const int*                        cginfo,
1321                                  ListOfLists<int>*                 lexcls,
1322                                  int                               iz,
1323                                  int                               at_start,
1324                                  int                               at_end,
1325                                  const gmx::ArrayRef<const int>    intermolecularExclusionGroup)
1326 {
1327     const gmx_ga2la_t& ga2la = *dd->ga2la;
1328 
1329     const auto& jAtomRange = zones->iZones[iz].jAtomRange;
1330 
1331     const gmx::index oldNumLists = lexcls->ssize();
1332 
1333     std::vector<int> exclusionsForAtom;
1334     for (int at = at_start; at < at_end; at++)
1335     {
1336         exclusionsForAtom.clear();
1337 
1338         if (GET_CGINFO_EXCL_INTER(cginfo[at]))
1339         {
1340             int a_gl, mb, mt, mol, a_mol;
1341 
1342             /* Copy the exclusions from the global top */
1343             a_gl = dd->globalAtomIndices[at];
1344             global_atomnr_to_moltype_ind(dd->reverse_top, a_gl, &mb, &mt, &mol, &a_mol);
1345             const auto excls = moltype[mt].excls[a_mol];
1346             for (const int aj_mol : excls)
1347             {
1348                 if (const auto* jEntry = ga2la.find(a_gl + aj_mol - a_mol))
1349                 {
1350                     /* This check is not necessary, but it can reduce
1351                      * the number of exclusions in the list, which in turn
1352                      * can speed up the pair list construction a bit.
1353                      */
1354                     if (jAtomRange.isInRange(jEntry->la))
1355                     {
1356                         exclusionsForAtom.push_back(jEntry->la);
1357                     }
1358                 }
1359             }
1360         }
1361 
1362         bool isExcludedAtom = !intermolecularExclusionGroup.empty()
1363                               && std::find(intermolecularExclusionGroup.begin(),
1364                                            intermolecularExclusionGroup.end(), dd->globalAtomIndices[at])
1365                                          != intermolecularExclusionGroup.end();
1366 
1367         if (isExcludedAtom)
1368         {
1369             for (int qmAtomGlobalIndex : intermolecularExclusionGroup)
1370             {
1371                 if (const auto* entry = dd->ga2la->find(qmAtomGlobalIndex))
1372                 {
1373                     exclusionsForAtom.push_back(entry->la);
1374                 }
1375             }
1376         }
1377 
1378         /* Append the exclusions for this atom to the topology */
1379         lexcls->pushBack(exclusionsForAtom);
1380     }
1381 
1382     GMX_RELEASE_ASSERT(
1383             lexcls->ssize() - oldNumLists == at_end - at_start,
1384             "The number of exclusion list should match the number of atoms in the range");
1385 }
1386 
1387 /*! \brief Generate and store all required local bonded interactions in \p idef and local exclusions in \p lexcls */
make_local_bondeds_excls(gmx_domdec_t * dd,gmx_domdec_zones_t * zones,const gmx_mtop_t * mtop,const int * cginfo,gmx_bool bRCheckMB,ivec rcheck,gmx_bool bRCheck2B,real rc,t_pbc * pbc_null,rvec * cg_cm,InteractionDefinitions * idef,ListOfLists<int> * lexcls,int * excl_count)1388 static int make_local_bondeds_excls(gmx_domdec_t*           dd,
1389                                     gmx_domdec_zones_t*     zones,
1390                                     const gmx_mtop_t*       mtop,
1391                                     const int*              cginfo,
1392                                     gmx_bool                bRCheckMB,
1393                                     ivec                    rcheck,
1394                                     gmx_bool                bRCheck2B,
1395                                     real                    rc,
1396                                     t_pbc*                  pbc_null,
1397                                     rvec*                   cg_cm,
1398                                     InteractionDefinitions* idef,
1399                                     ListOfLists<int>*       lexcls,
1400                                     int*                    excl_count)
1401 {
1402     int                nzone_bondeds;
1403     int                cg0, cg1;
1404     real               rc2;
1405     int                nbonded_local;
1406     gmx_reverse_top_t* rt;
1407 
1408     if (dd->reverse_top->bInterAtomicInteractions)
1409     {
1410         nzone_bondeds = zones->n;
1411     }
1412     else
1413     {
1414         /* Only single charge group (or atom) molecules, so interactions don't
1415          * cross zone boundaries and we only need to assign in the home zone.
1416          */
1417         nzone_bondeds = 1;
1418     }
1419 
1420     /* We only use exclusions from i-zones to i- and j-zones */
1421     const int numIZonesForExclusions = (dd->haveExclusions ? zones->iZones.size() : 0);
1422 
1423     rt = dd->reverse_top;
1424 
1425     rc2 = rc * rc;
1426 
1427     /* Clear the counts */
1428     idef->clear();
1429     nbonded_local = 0;
1430 
1431     lexcls->clear();
1432     *excl_count = 0;
1433 
1434     for (int izone = 0; izone < nzone_bondeds; izone++)
1435     {
1436         cg0 = zones->cg_range[izone];
1437         cg1 = zones->cg_range[izone + 1];
1438 
1439         const int numThreads = rt->th_work.size();
1440 #pragma omp parallel for num_threads(numThreads) schedule(static)
1441         for (int thread = 0; thread < numThreads; thread++)
1442         {
1443             try
1444             {
1445                 int                     cg0t, cg1t;
1446                 InteractionDefinitions* idef_t;
1447 
1448                 cg0t = cg0 + ((cg1 - cg0) * thread) / numThreads;
1449                 cg1t = cg0 + ((cg1 - cg0) * (thread + 1)) / numThreads;
1450 
1451                 if (thread == 0)
1452                 {
1453                     idef_t = idef;
1454                 }
1455                 else
1456                 {
1457                     idef_t = &rt->th_work[thread].idef;
1458                     idef_t->clear();
1459                 }
1460 
1461                 rt->th_work[thread].nbonded = make_bondeds_zone(
1462                         dd, zones, mtop->molblock, bRCheckMB, rcheck, bRCheck2B, rc2, pbc_null,
1463                         cg_cm, idef->iparams.data(), idef_t, izone, gmx::Range<int>(cg0t, cg1t));
1464 
1465                 if (izone < numIZonesForExclusions)
1466                 {
1467                     ListOfLists<int>* excl_t;
1468                     if (thread == 0)
1469                     {
1470                         // Thread 0 stores exclusions directly in the final storage
1471                         excl_t = lexcls;
1472                     }
1473                     else
1474                     {
1475                         // Threads > 0 store in temporary storage, starting at list index 0
1476                         excl_t = &rt->th_work[thread].excl;
1477                         excl_t->clear();
1478                     }
1479 
1480                     /* No charge groups and no distance check required */
1481                     make_exclusions_zone(dd, zones, mtop->moltype, cginfo, excl_t, izone, cg0t,
1482                                          cg1t, mtop->intermolecularExclusionGroup);
1483                 }
1484             }
1485             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1486         }
1487 
1488         if (rt->th_work.size() > 1)
1489         {
1490             combine_idef(idef, rt->th_work);
1491         }
1492 
1493         for (const thread_work_t& th_work : rt->th_work)
1494         {
1495             nbonded_local += th_work.nbonded;
1496         }
1497 
1498         if (izone < numIZonesForExclusions)
1499         {
1500             for (std::size_t th = 1; th < rt->th_work.size(); th++)
1501             {
1502                 lexcls->appendListOfLists(rt->th_work[th].excl);
1503             }
1504             for (const thread_work_t& th_work : rt->th_work)
1505             {
1506                 *excl_count += th_work.excl_count;
1507             }
1508         }
1509     }
1510 
1511     if (debug)
1512     {
1513         fprintf(debug, "We have %d exclusions, check count %d\n", lexcls->numElements(), *excl_count);
1514     }
1515 
1516     return nbonded_local;
1517 }
1518 
dd_make_local_top(gmx_domdec_t * dd,gmx_domdec_zones_t * zones,int npbcdim,matrix box,rvec cellsize_min,const ivec npulse,t_forcerec * fr,rvec * cgcm_or_x,const gmx_mtop_t & mtop,gmx_localtop_t * ltop)1519 void dd_make_local_top(gmx_domdec_t*       dd,
1520                        gmx_domdec_zones_t* zones,
1521                        int                 npbcdim,
1522                        matrix              box,
1523                        rvec                cellsize_min,
1524                        const ivec          npulse,
1525                        t_forcerec*         fr,
1526                        rvec*               cgcm_or_x,
1527                        const gmx_mtop_t&   mtop,
1528                        gmx_localtop_t*     ltop)
1529 {
1530     gmx_bool bRCheckMB, bRCheck2B;
1531     real     rc = -1;
1532     ivec     rcheck;
1533     int      d, nexcl;
1534     t_pbc    pbc, *pbc_null = nullptr;
1535 
1536     if (debug)
1537     {
1538         fprintf(debug, "Making local topology\n");
1539     }
1540 
1541     bRCheckMB = FALSE;
1542     bRCheck2B = FALSE;
1543 
1544     if (dd->reverse_top->bInterAtomicInteractions)
1545     {
1546         /* We need to check to which cell bondeds should be assigned */
1547         rc = dd_cutoff_twobody(dd);
1548         if (debug)
1549         {
1550             fprintf(debug, "Two-body bonded cut-off distance is %g\n", rc);
1551         }
1552 
1553         /* Should we check cg_cm distances when assigning bonded interactions? */
1554         for (d = 0; d < DIM; d++)
1555         {
1556             rcheck[d] = FALSE;
1557             /* Only need to check for dimensions where the part of the box
1558              * that is not communicated is smaller than the cut-off.
1559              */
1560             if (d < npbcdim && dd->numCells[d] > 1
1561                 && (dd->numCells[d] - npulse[d]) * cellsize_min[d] < 2 * rc)
1562             {
1563                 if (dd->numCells[d] == 2)
1564                 {
1565                     rcheck[d] = TRUE;
1566                     bRCheckMB = TRUE;
1567                 }
1568                 /* Check for interactions between two atoms,
1569                  * where we can allow interactions up to the cut-off,
1570                  * instead of up to the smallest cell dimension.
1571                  */
1572                 bRCheck2B = TRUE;
1573             }
1574             if (debug)
1575             {
1576                 fprintf(debug, "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %s\n", d,
1577                         cellsize_min[d], d, rcheck[d], gmx::boolToString(bRCheck2B));
1578             }
1579         }
1580         if (bRCheckMB || bRCheck2B)
1581         {
1582             if (fr->bMolPBC)
1583             {
1584                 pbc_null = set_pbc_dd(&pbc, fr->pbcType, dd->numCells, TRUE, box);
1585             }
1586             else
1587             {
1588                 pbc_null = nullptr;
1589             }
1590         }
1591     }
1592 
1593     dd->nbonded_local = make_local_bondeds_excls(dd, zones, &mtop, fr->cginfo.data(), bRCheckMB,
1594                                                  rcheck, bRCheck2B, rc, pbc_null, cgcm_or_x,
1595                                                  &ltop->idef, &ltop->excls, &nexcl);
1596 
1597     /* The ilist is not sorted yet,
1598      * we can only do this when we have the charge arrays.
1599      */
1600     ltop->idef.ilsort = ilsortUNKNOWN;
1601 }
1602 
dd_sort_local_top(gmx_domdec_t * dd,const t_mdatoms * mdatoms,gmx_localtop_t * ltop)1603 void dd_sort_local_top(gmx_domdec_t* dd, const t_mdatoms* mdatoms, gmx_localtop_t* ltop)
1604 {
1605     if (dd->reverse_top->ilsort == ilsortNO_FE)
1606     {
1607         ltop->idef.ilsort = ilsortNO_FE;
1608     }
1609     else
1610     {
1611         gmx_sort_ilist_fe(&ltop->idef, mdatoms->chargeA, mdatoms->chargeB);
1612     }
1613 }
1614 
dd_init_local_state(gmx_domdec_t * dd,const t_state * state_global,t_state * state_local)1615 void dd_init_local_state(gmx_domdec_t* dd, const t_state* state_global, t_state* state_local)
1616 {
1617     int buf[NITEM_DD_INIT_LOCAL_STATE];
1618 
1619     if (DDMASTER(dd))
1620     {
1621         buf[0] = state_global->flags;
1622         buf[1] = state_global->ngtc;
1623         buf[2] = state_global->nnhpres;
1624         buf[3] = state_global->nhchainlength;
1625         buf[4] = state_global->dfhist ? state_global->dfhist->nlambda : 0;
1626     }
1627     dd_bcast(dd, NITEM_DD_INIT_LOCAL_STATE * sizeof(int), buf);
1628 
1629     init_gtc_state(state_local, buf[1], buf[2], buf[3]);
1630     init_dfhist_state(state_local, buf[4]);
1631     state_local->flags = buf[0];
1632 }
1633 
1634 /*! \brief Check if a link is stored in \p link between charge groups \p cg_gl and \p cg_gl_j and if not so, store a link */
check_link(t_blocka * link,int cg_gl,int cg_gl_j)1635 static void check_link(t_blocka* link, int cg_gl, int cg_gl_j)
1636 {
1637     int      k;
1638     gmx_bool bFound;
1639 
1640     bFound = FALSE;
1641     for (k = link->index[cg_gl]; k < link->index[cg_gl + 1]; k++)
1642     {
1643         GMX_RELEASE_ASSERT(link->a, "Inconsistent NULL pointer while making charge-group links");
1644         if (link->a[k] == cg_gl_j)
1645         {
1646             bFound = TRUE;
1647         }
1648     }
1649     if (!bFound)
1650     {
1651         GMX_RELEASE_ASSERT(link->a || link->index[cg_gl + 1] + 1 > link->nalloc_a,
1652                            "Inconsistent allocation of link");
1653         /* Add this charge group link */
1654         if (link->index[cg_gl + 1] + 1 > link->nalloc_a)
1655         {
1656             link->nalloc_a = over_alloc_large(link->index[cg_gl + 1] + 1);
1657             srenew(link->a, link->nalloc_a);
1658         }
1659         link->a[link->index[cg_gl + 1]] = cg_gl_j;
1660         link->index[cg_gl + 1]++;
1661     }
1662 }
1663 
makeBondedLinks(const gmx_mtop_t & mtop,gmx::ArrayRef<cginfo_mb_t> cginfo_mb)1664 t_blocka* makeBondedLinks(const gmx_mtop_t& mtop, gmx::ArrayRef<cginfo_mb_t> cginfo_mb)
1665 {
1666     t_blocka*    link;
1667     cginfo_mb_t* cgi_mb;
1668 
1669     /* For each atom make a list of other atoms in the system
1670      * that a linked to it via bonded interactions
1671      * which are also stored in reverse_top.
1672      */
1673 
1674     reverse_ilist_t ril_intermol;
1675     if (mtop.bIntermolecularInteractions)
1676     {
1677         t_atoms atoms;
1678 
1679         atoms.nr   = mtop.natoms;
1680         atoms.atom = nullptr;
1681 
1682         GMX_RELEASE_ASSERT(mtop.intermolecular_ilist,
1683                            "We should have an ilist when intermolecular interactions are on");
1684 
1685         make_reverse_ilist(*mtop.intermolecular_ilist, &atoms, FALSE, FALSE, FALSE, TRUE, &ril_intermol);
1686     }
1687 
1688     snew(link, 1);
1689     snew(link->index, mtop.natoms + 1);
1690     link->nalloc_a = 0;
1691     link->a        = nullptr;
1692 
1693     link->index[0] = 0;
1694     int cg_offset  = 0;
1695     int ncgi       = 0;
1696     for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
1697     {
1698         const gmx_molblock_t& molb = mtop.molblock[mb];
1699         if (molb.nmol == 0)
1700         {
1701             continue;
1702         }
1703         const gmx_moltype_t& molt = mtop.moltype[molb.type];
1704         /* Make a reverse ilist in which the interactions are linked
1705          * to all atoms, not only the first atom as in gmx_reverse_top.
1706          * The constraints are discarded here.
1707          */
1708         reverse_ilist_t ril;
1709         make_reverse_ilist(molt.ilist, &molt.atoms, FALSE, FALSE, FALSE, TRUE, &ril);
1710 
1711         cgi_mb = &cginfo_mb[mb];
1712 
1713         int mol;
1714         for (mol = 0; mol < (mtop.bIntermolecularInteractions ? molb.nmol : 1); mol++)
1715         {
1716             for (int a = 0; a < molt.atoms.nr; a++)
1717             {
1718                 int cg_gl              = cg_offset + a;
1719                 link->index[cg_gl + 1] = link->index[cg_gl];
1720                 int i                  = ril.index[a];
1721                 while (i < ril.index[a + 1])
1722                 {
1723                     int ftype = ril.il[i++];
1724                     int nral  = NRAL(ftype);
1725                     /* Skip the ifunc index */
1726                     i++;
1727                     for (int j = 0; j < nral; j++)
1728                     {
1729                         int aj = ril.il[i + j];
1730                         if (aj != a)
1731                         {
1732                             check_link(link, cg_gl, cg_offset + aj);
1733                         }
1734                     }
1735                     i += nral_rt(ftype);
1736                 }
1737 
1738                 if (mtop.bIntermolecularInteractions)
1739                 {
1740                     int i = ril_intermol.index[cg_gl];
1741                     while (i < ril_intermol.index[cg_gl + 1])
1742                     {
1743                         int ftype = ril_intermol.il[i++];
1744                         int nral  = NRAL(ftype);
1745                         /* Skip the ifunc index */
1746                         i++;
1747                         for (int j = 0; j < nral; j++)
1748                         {
1749                             /* Here we assume we have no charge groups;
1750                              * this has been checked above.
1751                              */
1752                             int aj = ril_intermol.il[i + j];
1753                             check_link(link, cg_gl, aj);
1754                         }
1755                         i += nral_rt(ftype);
1756                     }
1757                 }
1758                 if (link->index[cg_gl + 1] - link->index[cg_gl] > 0)
1759                 {
1760                     SET_CGINFO_BOND_INTER(cgi_mb->cginfo[a]);
1761                     ncgi++;
1762                 }
1763             }
1764 
1765             cg_offset += molt.atoms.nr;
1766         }
1767         int nlink_mol = link->index[cg_offset] - link->index[cg_offset - molt.atoms.nr];
1768 
1769         if (debug)
1770         {
1771             fprintf(debug, "molecule type '%s' %d atoms has %d atom links through bonded interac.\n",
1772                     *molt.name, molt.atoms.nr, nlink_mol);
1773         }
1774 
1775         if (molb.nmol > mol)
1776         {
1777             /* Copy the data for the rest of the molecules in this block */
1778             link->nalloc_a += (molb.nmol - mol) * nlink_mol;
1779             srenew(link->a, link->nalloc_a);
1780             for (; mol < molb.nmol; mol++)
1781             {
1782                 for (int a = 0; a < molt.atoms.nr; a++)
1783                 {
1784                     int cg_gl              = cg_offset + a;
1785                     link->index[cg_gl + 1] = link->index[cg_gl + 1 - molt.atoms.nr] + nlink_mol;
1786                     for (int j = link->index[cg_gl]; j < link->index[cg_gl + 1]; j++)
1787                     {
1788                         link->a[j] = link->a[j - nlink_mol] + molt.atoms.nr;
1789                     }
1790                     if (link->index[cg_gl + 1] - link->index[cg_gl] > 0
1791                         && cg_gl - cgi_mb->cg_start < cgi_mb->cg_mod)
1792                     {
1793                         SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg_gl - cgi_mb->cg_start]);
1794                         ncgi++;
1795                     }
1796                 }
1797                 cg_offset += molt.atoms.nr;
1798             }
1799         }
1800     }
1801 
1802     if (debug)
1803     {
1804         fprintf(debug, "Of the %d atoms %d are linked via bonded interactions\n", mtop.natoms, ncgi);
1805     }
1806 
1807     return link;
1808 }
1809 
1810 typedef struct
1811 {
1812     real r2;
1813     int  ftype;
1814     int  a1;
1815     int  a2;
1816 } bonded_distance_t;
1817 
1818 /*! \brief Compare distance^2 \p r2 against the distance in \p bd and if larger store it along with \p ftype and atom indices \p a1 and \p a2 */
update_max_bonded_distance(real r2,int ftype,int a1,int a2,bonded_distance_t * bd)1819 static void update_max_bonded_distance(real r2, int ftype, int a1, int a2, bonded_distance_t* bd)
1820 {
1821     if (r2 > bd->r2)
1822     {
1823         bd->r2    = r2;
1824         bd->ftype = ftype;
1825         bd->a1    = a1;
1826         bd->a2    = a2;
1827     }
1828 }
1829 
1830 /*! \brief Set the distance, function type and atom indices for the longest distance between atoms of molecule type \p molt for two-body and multi-body bonded interactions */
bonded_cg_distance_mol(const gmx_moltype_t * molt,gmx_bool bBCheck,gmx_bool bExcl,ArrayRef<const RVec> x,bonded_distance_t * bd_2b,bonded_distance_t * bd_mb)1831 static void bonded_cg_distance_mol(const gmx_moltype_t* molt,
1832                                    gmx_bool             bBCheck,
1833                                    gmx_bool             bExcl,
1834                                    ArrayRef<const RVec> x,
1835                                    bonded_distance_t*   bd_2b,
1836                                    bonded_distance_t*   bd_mb)
1837 {
1838     for (int ftype = 0; ftype < F_NRE; ftype++)
1839     {
1840         if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
1841         {
1842             const auto& il   = molt->ilist[ftype];
1843             int         nral = NRAL(ftype);
1844             if (nral > 1)
1845             {
1846                 for (int i = 0; i < il.size(); i += 1 + nral)
1847                 {
1848                     for (int ai = 0; ai < nral; ai++)
1849                     {
1850                         int atomI = il.iatoms[i + 1 + ai];
1851                         for (int aj = ai + 1; aj < nral; aj++)
1852                         {
1853                             int atomJ = il.iatoms[i + 1 + aj];
1854                             if (atomI != atomJ)
1855                             {
1856                                 real rij2 = distance2(x[atomI], x[atomJ]);
1857 
1858                                 update_max_bonded_distance(rij2, ftype, atomI, atomJ,
1859                                                            (nral == 2) ? bd_2b : bd_mb);
1860                             }
1861                         }
1862                     }
1863                 }
1864             }
1865         }
1866     }
1867     if (bExcl)
1868     {
1869         const auto& excls = molt->excls;
1870         for (gmx::index ai = 0; ai < excls.ssize(); ai++)
1871         {
1872             for (const int aj : excls[ai])
1873             {
1874                 if (ai != aj)
1875                 {
1876                     real rij2 = distance2(x[ai], x[aj]);
1877 
1878                     /* There is no function type for exclusions, use -1 */
1879                     update_max_bonded_distance(rij2, -1, ai, aj, bd_2b);
1880                 }
1881             }
1882         }
1883     }
1884 }
1885 
1886 /*! \brief Set the distance, function type and atom indices for the longest atom distance involved in intermolecular interactions for two-body and multi-body bonded interactions */
bonded_distance_intermol(const InteractionLists & ilists_intermol,gmx_bool bBCheck,ArrayRef<const RVec> x,PbcType pbcType,const matrix box,bonded_distance_t * bd_2b,bonded_distance_t * bd_mb)1887 static void bonded_distance_intermol(const InteractionLists& ilists_intermol,
1888                                      gmx_bool                bBCheck,
1889                                      ArrayRef<const RVec>    x,
1890                                      PbcType                 pbcType,
1891                                      const matrix            box,
1892                                      bonded_distance_t*      bd_2b,
1893                                      bonded_distance_t*      bd_mb)
1894 {
1895     t_pbc pbc;
1896 
1897     set_pbc(&pbc, pbcType, box);
1898 
1899     for (int ftype = 0; ftype < F_NRE; ftype++)
1900     {
1901         if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
1902         {
1903             const auto& il   = ilists_intermol[ftype];
1904             int         nral = NRAL(ftype);
1905 
1906             /* No nral>1 check here, since intermol interactions always
1907              * have nral>=2 (and the code is also correct for nral=1).
1908              */
1909             for (int i = 0; i < il.size(); i += 1 + nral)
1910             {
1911                 for (int ai = 0; ai < nral; ai++)
1912                 {
1913                     int atom_i = il.iatoms[i + 1 + ai];
1914 
1915                     for (int aj = ai + 1; aj < nral; aj++)
1916                     {
1917                         rvec dx;
1918                         real rij2;
1919 
1920                         int atom_j = il.iatoms[i + 1 + aj];
1921 
1922                         pbc_dx(&pbc, x[atom_i], x[atom_j], dx);
1923 
1924                         rij2 = norm2(dx);
1925 
1926                         update_max_bonded_distance(rij2, ftype, atom_i, atom_j, (nral == 2) ? bd_2b : bd_mb);
1927                     }
1928                 }
1929             }
1930         }
1931     }
1932 }
1933 
1934 //! Returns whether \p molt has at least one virtual site
moltypeHasVsite(const gmx_moltype_t & molt)1935 static bool moltypeHasVsite(const gmx_moltype_t& molt)
1936 {
1937     bool hasVsite = false;
1938     for (int i = 0; i < F_NRE; i++)
1939     {
1940         if ((interaction_function[i].flags & IF_VSITE) && !molt.ilist[i].empty())
1941         {
1942             hasVsite = true;
1943         }
1944     }
1945 
1946     return hasVsite;
1947 }
1948 
1949 //! Returns coordinates not broken over PBC for a molecule
getWholeMoleculeCoordinates(const gmx_moltype_t * molt,const gmx_ffparams_t * ffparams,PbcType pbcType,t_graph * graph,const matrix box,ArrayRef<const RVec> x,ArrayRef<RVec> xs)1950 static void getWholeMoleculeCoordinates(const gmx_moltype_t*  molt,
1951                                         const gmx_ffparams_t* ffparams,
1952                                         PbcType               pbcType,
1953                                         t_graph*              graph,
1954                                         const matrix          box,
1955                                         ArrayRef<const RVec>  x,
1956                                         ArrayRef<RVec>        xs)
1957 {
1958     int n, i;
1959 
1960     if (pbcType != PbcType::No)
1961     {
1962         mk_mshift(nullptr, graph, pbcType, box, as_rvec_array(x.data()));
1963 
1964         shift_x(graph, box, as_rvec_array(x.data()), as_rvec_array(xs.data()));
1965         /* By doing an extra mk_mshift the molecules that are broken
1966          * because they were e.g. imported from another software
1967          * will be made whole again. Such are the healing powers
1968          * of GROMACS.
1969          */
1970         mk_mshift(nullptr, graph, pbcType, box, as_rvec_array(xs.data()));
1971     }
1972     else
1973     {
1974         /* We copy the coordinates so the original coordinates remain
1975          * unchanged, just to be 100% sure that we do not affect
1976          * binary reproducibility of simulations.
1977          */
1978         n = molt->atoms.nr;
1979         for (i = 0; i < n; i++)
1980         {
1981             copy_rvec(x[i], xs[i]);
1982         }
1983     }
1984 
1985     if (moltypeHasVsite(*molt))
1986     {
1987         gmx::constructVirtualSites(xs, ffparams->iparams, molt->ilist);
1988     }
1989 }
1990 
dd_bonded_cg_distance(const gmx::MDLogger & mdlog,const gmx_mtop_t * mtop,const t_inputrec * ir,ArrayRef<const RVec> x,const matrix box,gmx_bool bBCheck,real * r_2b,real * r_mb)1991 void dd_bonded_cg_distance(const gmx::MDLogger& mdlog,
1992                            const gmx_mtop_t*    mtop,
1993                            const t_inputrec*    ir,
1994                            ArrayRef<const RVec> x,
1995                            const matrix         box,
1996                            gmx_bool             bBCheck,
1997                            real*                r_2b,
1998                            real*                r_mb)
1999 {
2000     gmx_bool          bExclRequired;
2001     int               at_offset;
2002     bonded_distance_t bd_2b = { 0, -1, -1, -1 };
2003     bonded_distance_t bd_mb = { 0, -1, -1, -1 };
2004 
2005     bExclRequired = inputrecExclForces(ir);
2006 
2007     *r_2b     = 0;
2008     *r_mb     = 0;
2009     at_offset = 0;
2010     for (const gmx_molblock_t& molb : mtop->molblock)
2011     {
2012         const gmx_moltype_t& molt = mtop->moltype[molb.type];
2013         if (molt.atoms.nr == 1 || molb.nmol == 0)
2014         {
2015             at_offset += molb.nmol * molt.atoms.nr;
2016         }
2017         else
2018         {
2019             t_graph graph;
2020             if (ir->pbcType != PbcType::No)
2021             {
2022                 graph = mk_graph_moltype(molt);
2023             }
2024 
2025             std::vector<RVec> xs(molt.atoms.nr);
2026             for (int mol = 0; mol < molb.nmol; mol++)
2027             {
2028                 getWholeMoleculeCoordinates(&molt, &mtop->ffparams, ir->pbcType, &graph, box,
2029                                             x.subArray(at_offset, molt.atoms.nr), xs);
2030 
2031                 bonded_distance_t bd_mol_2b = { 0, -1, -1, -1 };
2032                 bonded_distance_t bd_mol_mb = { 0, -1, -1, -1 };
2033 
2034                 bonded_cg_distance_mol(&molt, bBCheck, bExclRequired, xs, &bd_mol_2b, &bd_mol_mb);
2035 
2036                 /* Process the mol data adding the atom index offset */
2037                 update_max_bonded_distance(bd_mol_2b.r2, bd_mol_2b.ftype, at_offset + bd_mol_2b.a1,
2038                                            at_offset + bd_mol_2b.a2, &bd_2b);
2039                 update_max_bonded_distance(bd_mol_mb.r2, bd_mol_mb.ftype, at_offset + bd_mol_mb.a1,
2040                                            at_offset + bd_mol_mb.a2, &bd_mb);
2041 
2042                 at_offset += molt.atoms.nr;
2043             }
2044         }
2045     }
2046 
2047     if (mtop->bIntermolecularInteractions)
2048     {
2049         GMX_RELEASE_ASSERT(mtop->intermolecular_ilist,
2050                            "We should have an ilist when intermolecular interactions are on");
2051 
2052         bonded_distance_intermol(*mtop->intermolecular_ilist, bBCheck, x, ir->pbcType, box, &bd_2b, &bd_mb);
2053     }
2054 
2055     *r_2b = sqrt(bd_2b.r2);
2056     *r_mb = sqrt(bd_mb.r2);
2057 
2058     if (*r_2b > 0 || *r_mb > 0)
2059     {
2060         GMX_LOG(mdlog.info).appendText("Initial maximum distances in bonded interactions:");
2061         if (*r_2b > 0)
2062         {
2063             GMX_LOG(mdlog.info)
2064                     .appendTextFormatted(
2065                             "    two-body bonded interactions: %5.3f nm, %s, atoms %d %d", *r_2b,
2066                             (bd_2b.ftype >= 0) ? interaction_function[bd_2b.ftype].longname : "Exclusion",
2067                             bd_2b.a1 + 1, bd_2b.a2 + 1);
2068         }
2069         if (*r_mb > 0)
2070         {
2071             GMX_LOG(mdlog.info)
2072                     .appendTextFormatted(
2073                             "  multi-body bonded interactions: %5.3f nm, %s, atoms %d %d", *r_mb,
2074                             interaction_function[bd_mb.ftype].longname, bd_mb.a1 + 1, bd_mb.a2 + 1);
2075         }
2076     }
2077 }
2078