2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006 - 2014, The GROMACS development team.
5 * Copyright (c) 2015,2016,2017,2018,2019, 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.
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.
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.
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.
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.
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.
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.
43 * \author Berk Hess <hess@kth.se>
44 * \ingroup module_domdec
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_network.h"
59 #include "gromacs/domdec/ga2la.h"
60 #include "gromacs/gmxlib/chargegroup.h"
61 #include "gromacs/gmxlib/network.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdlib/forcerec.h"
64 #include "gromacs/mdlib/gmx_omp_nthreads.h"
65 #include "gromacs/mdtypes/commrec.h"
66 #include "gromacs/mdtypes/inputrec.h"
67 #include "gromacs/mdtypes/md_enums.h"
68 #include "gromacs/mdtypes/mdatom.h"
69 #include "gromacs/mdtypes/state.h"
70 #include "gromacs/pbcutil/mshift.h"
71 #include "gromacs/pbcutil/pbc.h"
72 #include "gromacs/topology/ifunc.h"
73 #include "gromacs/topology/mtop_util.h"
74 #include "gromacs/topology/topsort.h"
75 #include "gromacs/utility/cstringutil.h"
76 #include "gromacs/utility/exceptions.h"
77 #include "gromacs/utility/fatalerror.h"
78 #include "gromacs/utility/gmxassert.h"
79 #include "gromacs/utility/logger.h"
80 #include "gromacs/utility/smalloc.h"
81 #include "gromacs/utility/strconvert.h"
82 #include "gromacs/utility/stringstream.h"
83 #include "gromacs/utility/stringutil.h"
84 #include "gromacs/utility/textwriter.h"
86 #include "domdec_constraints.h"
87 #include "domdec_internal.h"
88 #include "domdec_vsite.h"
91 /*! \brief The number of integer item in the local state, used for broadcasting of the state */
92 #define NITEM_DD_INIT_LOCAL_STATE 5
94 struct reverse_ilist_t
96 std::vector
<int> index
; /* Index for each atom into il */
97 std::vector
<int> il
; /* ftype|type|a0|...|an|ftype|... */
98 int numAtomsInMolecule
; /* The number of atoms in this molecule */
101 struct MolblockIndices
109 /*! \brief Struct for thread local work data for local topology generation */
112 t_idef idef
; /**< Partial local topology */
113 std::unique_ptr
<VsitePbc
> vsitePbc
; /**< vsite PBC structure */
114 int nbonded
; /**< The number of bondeds in this struct */
115 t_blocka excl
; /**< List of exclusions */
116 int excl_count
; /**< The total exclusion count for \p excl */
119 /*! \brief Struct for the reverse topology: links bonded interactions to atomsx */
120 struct gmx_reverse_top_t
122 //! @cond Doxygen_Suppress
123 //! \brief Do we require all exclusions to be assigned?
124 bool bExclRequired
= false;
125 //! \brief The maximum number of exclusions one atom can have
126 int n_excl_at_max
= 0;
127 //! \brief Are there constraints in this revserse top?
128 bool bConstr
= false;
129 //! \brief Are there settles in this revserse top?
130 bool bSettle
= false;
131 //! \brief All bonded interactions have to be assigned?
132 bool bBCheck
= false;
133 //! \brief Are there bondeds/exclusions between charge-groups?
134 bool bInterCGInteractions
= false;
135 //! \brief Reverse ilist for all moltypes
136 std::vector
<reverse_ilist_t
> ril_mt
;
137 //! \brief The size of ril_mt[?].index summed over all entries
138 int ril_mt_tot_size
= 0;
139 //! \brief The sorting state of bondeds for free energy
140 int ilsort
= ilsortUNKNOWN
;
141 //! \brief molblock to global atom index for quick lookup of molblocks on atom index
142 std::vector
<MolblockIndices
> mbi
;
144 //! \brief Do we have intermolecular interactions?
145 bool bIntermolecularInteractions
= false;
146 //! \brief Intermolecular reverse ilist
147 reverse_ilist_t ril_intermol
;
149 /* Work data structures for multi-threading */
150 //! \brief Thread work array for local topology generation
151 std::vector
<thread_work_t
> th_work
;
155 /*! \brief Returns the number of atom entries for il in gmx_reverse_top_t */
156 static int nral_rt(int ftype
)
161 if (interaction_function
[ftype
].flags
& IF_VSITE
)
163 /* With vsites the reverse topology contains an extra entry
164 * for storing if constructing atoms are vsites.
172 /*! \brief Return whether interactions of type \p ftype need to be assigned exactly once */
173 static gmx_bool
dd_check_ftype(int ftype
, gmx_bool bBCheck
,
174 gmx_bool bConstr
, gmx_bool bSettle
)
176 return ((((interaction_function
[ftype
].flags
& IF_BOND
) != 0u) &&
177 ((interaction_function
[ftype
].flags
& IF_VSITE
) == 0u) &&
178 (bBCheck
|| ((interaction_function
[ftype
].flags
& IF_LIMZERO
) == 0u))) ||
179 (bConstr
&& (ftype
== F_CONSTR
|| ftype
== F_CONSTRNC
)) ||
180 (bSettle
&& ftype
== F_SETTLE
));
183 /*! \brief Help print error output when interactions are missing */
185 print_missing_interactions_mb(t_commrec
*cr
,
186 const gmx_reverse_top_t
*rt
,
187 const char *moltypename
,
188 const reverse_ilist_t
*ril
,
189 int a_start
, int a_end
,
190 int nat_mol
, int nmol
,
194 int nril_mol
= ril
->index
[nat_mol
];
195 snew(assigned
, nmol
*nril_mol
);
196 gmx::StringOutputStream stream
;
197 gmx::TextWriter
log(&stream
);
199 gmx::ArrayRef
<const int> gatindex
= cr
->dd
->globalAtomIndices
;
200 for (int ftype
= 0; ftype
< F_NRE
; ftype
++)
202 if (dd_check_ftype(ftype
, rt
->bBCheck
, rt
->bConstr
, rt
->bSettle
))
204 int nral
= NRAL(ftype
);
205 const t_ilist
*il
= &idef
->il
[ftype
];
206 const t_iatom
*ia
= il
->iatoms
;
207 for (int i
= 0; i
< il
->nr
; i
+= 1+nral
)
209 int a0
= gatindex
[ia
[1]];
210 /* Check if this interaction is in
211 * the currently checked molblock.
213 if (a0
>= a_start
&& a0
< a_end
)
215 int mol
= (a0
- a_start
)/nat_mol
;
216 int a0_mol
= (a0
- a_start
) - mol
*nat_mol
;
217 int j_mol
= ril
->index
[a0_mol
];
219 while (j_mol
< ril
->index
[a0_mol
+1] && !found
)
221 int j
= mol
*nril_mol
+ j_mol
;
222 int ftype_j
= ril
->il
[j_mol
];
223 /* Here we need to check if this interaction has
224 * not already been assigned, since we could have
225 * multiply defined interactions.
227 if (ftype
== ftype_j
&& ia
[0] == ril
->il
[j_mol
+1] &&
230 /* Check the atoms */
232 for (int a
= 0; a
< nral
; a
++)
234 if (gatindex
[ia
[1+a
]] !=
235 a_start
+ mol
*nat_mol
+ ril
->il
[j_mol
+2+a
])
245 j_mol
+= 2 + nral_rt(ftype_j
);
249 gmx_incons("Some interactions seem to be assigned multiple times");
257 gmx_sumi(nmol
*nril_mol
, assigned
, cr
);
261 for (int mol
= 0; mol
< nmol
; mol
++)
264 while (j_mol
< nril_mol
)
266 int ftype
= ril
->il
[j_mol
];
267 int nral
= NRAL(ftype
);
268 int j
= mol
*nril_mol
+ j_mol
;
269 if (assigned
[j
] == 0 &&
270 !(interaction_function
[ftype
].flags
& IF_VSITE
))
272 if (DDMASTER(cr
->dd
))
276 log
.writeLineFormatted("Molecule type '%s'", moltypename
);
277 log
.writeLineFormatted(
278 "the first %d missing interactions, except for exclusions:", nprint
);
280 log
.writeStringFormatted("%20s atoms",
281 interaction_function
[ftype
].longname
);
283 for (a
= 0; a
< nral
; a
++)
285 log
.writeStringFormatted("%5d", ril
->il
[j_mol
+2+a
]+1);
289 log
.writeString(" ");
292 log
.writeString(" global");
293 for (a
= 0; a
< nral
; a
++)
295 log
.writeStringFormatted("%6d",
296 a_start
+mol
*nat_mol
+ril
->il
[j_mol
+2+a
]+1);
298 log
.ensureLineBreak();
306 j_mol
+= 2 + nral_rt(ftype
);
311 return stream
.toString();
314 /*! \brief Help print error output when interactions are missing */
315 static void print_missing_interactions_atoms(const gmx::MDLogger
&mdlog
,
317 const gmx_mtop_t
*mtop
,
320 const gmx_reverse_top_t
*rt
= cr
->dd
->reverse_top
;
322 /* Print the atoms in the missing interactions per molblock */
324 for (const gmx_molblock_t
&molb
: mtop
->molblock
)
326 const gmx_moltype_t
&moltype
= mtop
->moltype
[molb
.type
];
328 a_end
= a_start
+ molb
.nmol
*moltype
.atoms
.nr
;
330 GMX_LOG(mdlog
.warning
).appendText(
331 print_missing_interactions_mb(cr
, rt
,
333 &rt
->ril_mt
[molb
.type
],
334 a_start
, a_end
, moltype
.atoms
.nr
,
340 void dd_print_missing_interactions(const gmx::MDLogger
&mdlog
,
343 const gmx_mtop_t
*top_global
,
344 const gmx_localtop_t
*top_local
,
345 const t_state
*state_local
)
347 int ndiff_tot
, cl
[F_NRE
], n
, ndiff
, rest_global
, rest_local
;
353 GMX_LOG(mdlog
.warning
).appendText(
354 "Not all bonded interactions have been properly assigned to the domain decomposition cells");
356 ndiff_tot
= local_count
- dd
->nbonded_global
;
358 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
361 cl
[ftype
] = top_local
->idef
.il
[ftype
].nr
/(1+nral
);
364 gmx_sumi(F_NRE
, cl
, cr
);
368 GMX_LOG(mdlog
.warning
).appendText("A list of missing interactions:");
369 rest_global
= dd
->nbonded_global
;
370 rest_local
= local_count
;
371 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
373 /* In the reverse and local top all constraints are merged
374 * into F_CONSTR. So in the if statement we skip F_CONSTRNC
375 * and add these constraints when doing F_CONSTR.
377 if (((interaction_function
[ftype
].flags
& IF_BOND
) &&
378 (dd
->reverse_top
->bBCheck
379 || !(interaction_function
[ftype
].flags
& IF_LIMZERO
)))
380 || (dd
->reverse_top
->bConstr
&& ftype
== F_CONSTR
)
381 || (dd
->reverse_top
->bSettle
&& ftype
== F_SETTLE
))
383 n
= gmx_mtop_ftype_count(top_global
, ftype
);
384 if (ftype
== F_CONSTR
)
386 n
+= gmx_mtop_ftype_count(top_global
, F_CONSTRNC
);
388 ndiff
= cl
[ftype
] - n
;
391 GMX_LOG(mdlog
.warning
).appendTextFormatted(
392 "%20s of %6d missing %6d",
393 interaction_function
[ftype
].longname
, n
, -ndiff
);
396 rest_local
-= cl
[ftype
];
400 ndiff
= rest_local
- rest_global
;
403 GMX_LOG(mdlog
.warning
).appendTextFormatted(
404 "%20s of %6d missing %6d", "exclusions",
405 rest_global
, -ndiff
);
409 print_missing_interactions_atoms(mdlog
, cr
, top_global
, &top_local
->idef
);
410 write_dd_pdb("dd_dump_err", 0, "dump", top_global
, cr
,
411 -1, state_local
->x
.rvec_array(), state_local
->box
);
413 std::string errorMessage
;
417 errorMessage
= "One or more interactions were assigned to multiple domains of the domain decompostion. Please report this bug.";
421 errorMessage
= gmx::formatString("%d of the %d bonded interactions could not be calculated because some atoms involved moved further apart than the multi-body cut-off distance (%g nm) or the two-body cut-off distance (%g nm), see option -rdd, for pairs and tabulated bonds also see option -ddcheck", -ndiff_tot
, cr
->dd
->nbonded_global
, dd_cutoff_multibody(dd
), dd_cutoff_twobody(dd
));
423 gmx_fatal_collective(FARGS
, cr
->mpi_comm_mygroup
, MASTER(cr
), "%s", errorMessage
.c_str());
426 /*! \brief Return global topology molecule information for global atom index \p i_gl */
427 static void global_atomnr_to_moltype_ind(const gmx_reverse_top_t
*rt
,
429 int *mb
, int *mt
, int *mol
, int *i_mol
)
431 const MolblockIndices
*mbi
= rt
->mbi
.data();
433 int end
= rt
->mbi
.size(); /* exclusive */
436 /* binary search for molblock_ind */
440 if (i_gl
>= mbi
[mid
].a_end
)
444 else if (i_gl
< mbi
[mid
].a_start
)
458 *mol
= (i_gl
- mbi
->a_start
) / mbi
->natoms_mol
;
459 *i_mol
= (i_gl
- mbi
->a_start
) - (*mol
)*mbi
->natoms_mol
;
462 /*! \brief Count the exclusions for all atoms in \p cgs */
463 static void count_excls(const t_block
*cgs
, const t_blocka
*excls
,
464 int *n_excl
, int *n_intercg_excl
, int *n_excl_at_max
)
466 int cg
, at0
, at1
, at
, excl
, atj
;
471 for (cg
= 0; cg
< cgs
->nr
; cg
++)
473 at0
= cgs
->index
[cg
];
474 at1
= cgs
->index
[cg
+1];
475 for (at
= at0
; at
< at1
; at
++)
477 for (excl
= excls
->index
[at
]; excl
< excls
->index
[at
+1]; excl
++)
479 atj
= excls
->a
[excl
];
483 if (atj
< at0
|| atj
>= at1
)
490 *n_excl_at_max
= std::max(*n_excl_at_max
,
491 excls
->index
[at
+1] - excls
->index
[at
]);
496 /*! \brief Run the reverse ilist generation and store it in r_il when \p bAssign = TRUE */
497 static int low_make_reverse_ilist(const InteractionLists
&il_mt
,
500 gmx_bool bConstr
, gmx_bool bSettle
,
502 gmx::ArrayRef
<const int> r_index
,
503 gmx::ArrayRef
<int> r_il
,
504 gmx_bool bLinkToAllAtoms
,
507 int ftype
, j
, nlink
, link
;
512 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
514 if ((interaction_function
[ftype
].flags
& (IF_BOND
| IF_VSITE
)) ||
515 (bConstr
&& (ftype
== F_CONSTR
|| ftype
== F_CONSTRNC
)) ||
516 (bSettle
&& ftype
== F_SETTLE
))
518 const bool bVSite
= ((interaction_function
[ftype
].flags
& IF_VSITE
) != 0u);
519 const int nral
= NRAL(ftype
);
520 const auto &il
= il_mt
[ftype
];
521 for (int i
= 0; i
< il
.size(); i
+= 1+nral
)
523 const int* ia
= il
.iatoms
.data() + i
;
528 /* We don't need the virtual sites for the cg-links */
538 /* Couple to the first atom in the interaction */
541 for (link
= 0; link
< nlink
; link
++)
546 GMX_ASSERT(!r_il
.empty(), "with bAssign not allowed to be empty");
547 GMX_ASSERT(!r_index
.empty(), "with bAssign not allowed to be empty");
548 r_il
[r_index
[a
]+count
[a
]] =
549 (ftype
== F_CONSTRNC
? F_CONSTR
: ftype
);
550 r_il
[r_index
[a
]+count
[a
]+1] = ia
[0];
551 for (j
= 1; j
< 1+nral
; j
++)
553 /* Store the molecular atom number */
554 r_il
[r_index
[a
]+count
[a
]+1+j
] = ia
[j
];
557 if (interaction_function
[ftype
].flags
& IF_VSITE
)
561 /* Add an entry to iatoms for storing
562 * which of the constructing atoms are
565 r_il
[r_index
[a
]+count
[a
]+2+nral
] = 0;
566 for (j
= 2; j
< 1+nral
; j
++)
568 if (atom
[ia
[j
]].ptype
== eptVSite
)
570 r_il
[r_index
[a
]+count
[a
]+2+nral
] |= (2<<j
);
577 /* We do not count vsites since they are always
578 * uniquely assigned and can be assigned
579 * to multiple nodes with recursive vsites.
582 !(interaction_function
[ftype
].flags
& IF_LIMZERO
))
587 count
[a
] += 2 + nral_rt(ftype
);
596 /*! \brief Make the reverse ilist: a list of bonded interactions linked to atoms */
597 static int make_reverse_ilist(const InteractionLists
&ilist
,
598 const t_atoms
*atoms
,
599 gmx_bool bConstr
, gmx_bool bSettle
,
601 gmx_bool bLinkToAllAtoms
,
602 reverse_ilist_t
*ril_mt
)
604 int nat_mt
, *count
, i
, nint_mt
;
606 /* Count the interactions */
609 low_make_reverse_ilist(ilist
, atoms
->atom
,
611 bConstr
, bSettle
, bBCheck
,
613 bLinkToAllAtoms
, FALSE
);
615 ril_mt
->index
.push_back(0);
616 for (i
= 0; i
< nat_mt
; i
++)
618 ril_mt
->index
.push_back(ril_mt
->index
[i
] + count
[i
]);
621 ril_mt
->il
.resize(ril_mt
->index
[nat_mt
]);
623 /* Store the interactions */
625 low_make_reverse_ilist(ilist
, atoms
->atom
,
627 bConstr
, bSettle
, bBCheck
,
628 ril_mt
->index
, ril_mt
->il
,
629 bLinkToAllAtoms
, TRUE
);
633 ril_mt
->numAtomsInMolecule
= atoms
->nr
;
638 /*! \brief Generate the reverse topology */
639 static gmx_reverse_top_t
make_reverse_top(const gmx_mtop_t
*mtop
, gmx_bool bFE
,
640 gmx_bool bConstr
, gmx_bool bSettle
,
641 gmx_bool bBCheck
, int *nint
)
643 gmx_reverse_top_t rt
;
645 /* Should we include constraints (for SHAKE) in rt? */
646 rt
.bConstr
= bConstr
;
647 rt
.bSettle
= bSettle
;
648 rt
.bBCheck
= bBCheck
;
650 rt
.bInterCGInteractions
= mtop
->bIntermolecularInteractions
;
651 rt
.ril_mt
.resize(mtop
->moltype
.size());
652 rt
.ril_mt_tot_size
= 0;
653 std::vector
<int> nint_mt
;
654 for (size_t mt
= 0; mt
< mtop
->moltype
.size(); mt
++)
656 const gmx_moltype_t
&molt
= mtop
->moltype
[mt
];
659 rt
.bInterCGInteractions
= true;
662 /* Make the atom to interaction list for this molecule type */
663 int numberOfInteractions
=
664 make_reverse_ilist(molt
.ilist
, &molt
.atoms
,
665 rt
.bConstr
, rt
.bSettle
, rt
.bBCheck
, FALSE
,
667 nint_mt
.push_back(numberOfInteractions
);
669 rt
.ril_mt_tot_size
+= rt
.ril_mt
[mt
].index
[molt
.atoms
.nr
];
673 fprintf(debug
, "The total size of the atom to interaction index is %d integers\n", rt
.ril_mt_tot_size
);
677 for (const gmx_molblock_t
&molblock
: mtop
->molblock
)
679 *nint
+= molblock
.nmol
*nint_mt
[molblock
.type
];
682 /* Make an intermolecular reverse top, if necessary */
683 rt
.bIntermolecularInteractions
= mtop
->bIntermolecularInteractions
;
684 if (rt
.bIntermolecularInteractions
)
686 t_atoms atoms_global
;
688 atoms_global
.nr
= mtop
->natoms
;
689 atoms_global
.atom
= nullptr; /* Only used with virtual sites */
691 GMX_RELEASE_ASSERT(mtop
->intermolecular_ilist
, "We should have an ilist when intermolecular interactions are on");
694 make_reverse_ilist(*mtop
->intermolecular_ilist
,
696 rt
.bConstr
, rt
.bSettle
, rt
.bBCheck
, FALSE
,
700 if (bFE
&& gmx_mtop_bondeds_free_energy(mtop
))
702 rt
.ilsort
= ilsortFE_UNSORTED
;
706 rt
.ilsort
= ilsortNO_FE
;
709 /* Make a molblock index for fast searching */
711 for (size_t mb
= 0; mb
< mtop
->molblock
.size(); mb
++)
713 const gmx_molblock_t
&molb
= mtop
->molblock
[mb
];
714 const int numAtomsPerMol
= mtop
->moltype
[molb
.type
].atoms
.nr
;
717 i
+= molb
.nmol
*numAtomsPerMol
;
719 mbi
.natoms_mol
= numAtomsPerMol
;
720 mbi
.type
= molb
.type
;
721 rt
.mbi
.push_back(mbi
);
724 rt
.th_work
.resize(gmx_omp_nthreads_get(emntDomdec
));
729 void dd_make_reverse_top(FILE *fplog
,
730 gmx_domdec_t
*dd
, const gmx_mtop_t
*mtop
,
731 const gmx_vsite_t
*vsite
,
732 const t_inputrec
*ir
, gmx_bool bBCheck
)
736 fprintf(fplog
, "\nLinking all bonded interactions to atoms\n");
739 /* If normal and/or settle constraints act only within charge groups,
740 * we can store them in the reverse top and simply assign them to domains.
741 * Otherwise we need to assign them to multiple domains and set up
742 * the parallel version constraint algorithm(s).
745 dd
->reverse_top
= new gmx_reverse_top_t
;
747 make_reverse_top(mtop
, ir
->efep
!= efepNO
,
748 !dd
->splitConstraints
, !dd
->splitSettles
,
749 bBCheck
, &dd
->nbonded_global
);
751 gmx_reverse_top_t
*rt
= dd
->reverse_top
;
753 /* With the Verlet scheme, exclusions are handled in the non-bonded
754 * kernels and only exclusions inside the cut-off lead to exclusion
755 * forces. Since each atom pair is treated at most once in the non-bonded
756 * kernels, it doesn't matter if the exclusions for the same atom pair
757 * appear multiple times in the exclusion list. In contrast, the, old,
758 * group cut-off scheme loops over a list of exclusions, so there each
759 * excluded pair should appear exactly once.
761 rt
->bExclRequired
= (ir
->cutoff_scheme
== ecutsGROUP
&&
762 inputrecExclForces(ir
));
765 dd
->n_intercg_excl
= 0;
766 rt
->n_excl_at_max
= 0;
767 for (const gmx_molblock_t
&molb
: mtop
->molblock
)
769 int n_excl_mol
, n_excl_icg
, n_excl_at_max
;
771 const gmx_moltype_t
&molt
= mtop
->moltype
[molb
.type
];
772 count_excls(&molt
.cgs
, &molt
.excls
,
773 &n_excl_mol
, &n_excl_icg
, &n_excl_at_max
);
774 nexcl
+= molb
.nmol
*n_excl_mol
;
775 dd
->n_intercg_excl
+= molb
.nmol
*n_excl_icg
;
776 rt
->n_excl_at_max
= std::max(rt
->n_excl_at_max
, n_excl_at_max
);
778 if (rt
->bExclRequired
)
780 dd
->nbonded_global
+= nexcl
;
781 if (EEL_FULL(ir
->coulombtype
) && dd
->n_intercg_excl
> 0 && fplog
)
783 fprintf(fplog
, "There are %d inter charge-group exclusions,\n"
784 "will use an extra communication step for exclusion forces for %s\n",
785 dd
->n_intercg_excl
, eel_names
[ir
->coulombtype
]);
789 if (vsite
&& vsite
->numInterUpdategroupVsites
> 0)
793 fprintf(fplog
, "There are %d inter update-group virtual sites,\n"
794 "will an extra communication step for selected coordinates and forces\n",
795 vsite
->numInterUpdategroupVsites
);
797 init_domdec_vsites(dd
, vsite
->numInterUpdategroupVsites
);
800 if (dd
->splitConstraints
|| dd
->splitSettles
)
802 init_domdec_constraints(dd
, mtop
);
806 fprintf(fplog
, "\n");
810 /*! \brief Store a vsite interaction at the end of \p il
812 * This routine is very similar to add_ifunc, but vsites interactions
813 * have more work to do than other kinds of interactions, and the
814 * complex way nral (and thus vector contents) depends on ftype
815 * confuses static analysis tools unless we fuse the vsite
816 * atom-indexing organization code with the ifunc-adding code, so that
817 * they can see that nral is the same value. */
819 add_ifunc_for_vsites(t_iatom
*tiatoms
, const gmx_ga2la_t
&ga2la
,
820 int nral
, gmx_bool bHomeA
,
821 int a
, int a_gl
, int a_mol
,
822 const t_iatom
*iatoms
,
827 if (il
->nr
+1+nral
> il
->nalloc
)
829 il
->nalloc
= over_alloc_large(il
->nr
+1+nral
);
830 srenew(il
->iatoms
, il
->nalloc
);
832 liatoms
= il
->iatoms
+ il
->nr
;
836 tiatoms
[0] = iatoms
[0];
840 /* We know the local index of the first atom */
845 /* Convert later in make_local_vsites */
846 tiatoms
[1] = -a_gl
- 1;
849 for (int k
= 2; k
< 1+nral
; k
++)
851 int ak_gl
= a_gl
+ iatoms
[k
] - a_mol
;
852 if (const int *homeIndex
= ga2la
.findHome(ak_gl
))
854 tiatoms
[k
] = *homeIndex
;
858 /* Copy the global index, convert later in make_local_vsites */
859 tiatoms
[k
] = -(ak_gl
+ 1);
861 // Note that ga2la_get_home always sets the third parameter if
864 for (int k
= 0; k
< 1+nral
; k
++)
866 liatoms
[k
] = tiatoms
[k
];
870 /*! \brief Store a bonded interaction at the end of \p il */
871 static inline void add_ifunc(int nral
, const t_iatom
*tiatoms
, t_ilist
*il
)
876 if (il
->nr
+1+nral
> il
->nalloc
)
878 il
->nalloc
= over_alloc_large(il
->nr
+1+nral
);
879 srenew(il
->iatoms
, il
->nalloc
);
881 liatoms
= il
->iatoms
+ il
->nr
;
882 for (k
= 0; k
<= nral
; k
++)
884 liatoms
[k
] = tiatoms
[k
];
889 /*! \brief Store a position restraint in idef and iatoms, complex because the parameters are different for each entry */
890 static void add_posres(int mol
, int a_mol
, int numAtomsInMolecule
,
891 const gmx_molblock_t
*molb
,
892 t_iatom
*iatoms
, const t_iparams
*ip_in
,
898 /* This position restraint has not been added yet,
899 * so it's index is the current number of position restraints.
901 n
= idef
->il
[F_POSRES
].nr
/2;
902 if (n
+1 > idef
->iparams_posres_nalloc
)
904 idef
->iparams_posres_nalloc
= over_alloc_dd(n
+1);
905 srenew(idef
->iparams_posres
, idef
->iparams_posres_nalloc
);
907 ip
= &idef
->iparams_posres
[n
];
908 /* Copy the force constants */
909 *ip
= ip_in
[iatoms
[0]];
911 /* Get the position restraint coordinates from the molblock */
912 a_molb
= mol
*numAtomsInMolecule
+ a_mol
;
913 GMX_ASSERT(a_molb
< ssize(molb
->posres_xA
), "We need a sufficient number of position restraint coordinates");
914 ip
->posres
.pos0A
[XX
] = molb
->posres_xA
[a_molb
][XX
];
915 ip
->posres
.pos0A
[YY
] = molb
->posres_xA
[a_molb
][YY
];
916 ip
->posres
.pos0A
[ZZ
] = molb
->posres_xA
[a_molb
][ZZ
];
917 if (!molb
->posres_xB
.empty())
919 ip
->posres
.pos0B
[XX
] = molb
->posres_xB
[a_molb
][XX
];
920 ip
->posres
.pos0B
[YY
] = molb
->posres_xB
[a_molb
][YY
];
921 ip
->posres
.pos0B
[ZZ
] = molb
->posres_xB
[a_molb
][ZZ
];
925 ip
->posres
.pos0B
[XX
] = ip
->posres
.pos0A
[XX
];
926 ip
->posres
.pos0B
[YY
] = ip
->posres
.pos0A
[YY
];
927 ip
->posres
.pos0B
[ZZ
] = ip
->posres
.pos0A
[ZZ
];
929 /* Set the parameter index for idef->iparams_posre */
933 /*! \brief Store a flat-bottomed position restraint in idef and iatoms, complex because the parameters are different for each entry */
934 static void add_fbposres(int mol
, int a_mol
, int numAtomsInMolecule
,
935 const gmx_molblock_t
*molb
,
936 t_iatom
*iatoms
, const t_iparams
*ip_in
,
942 /* This flat-bottom position restraint has not been added yet,
943 * so it's index is the current number of position restraints.
945 n
= idef
->il
[F_FBPOSRES
].nr
/2;
946 if (n
+1 > idef
->iparams_fbposres_nalloc
)
948 idef
->iparams_fbposres_nalloc
= over_alloc_dd(n
+1);
949 srenew(idef
->iparams_fbposres
, idef
->iparams_fbposres_nalloc
);
951 ip
= &idef
->iparams_fbposres
[n
];
952 /* Copy the force constants */
953 *ip
= ip_in
[iatoms
[0]];
955 /* Get the position restraint coordinats from the molblock */
956 a_molb
= mol
*numAtomsInMolecule
+ a_mol
;
957 GMX_ASSERT(a_molb
< ssize(molb
->posres_xA
), "We need a sufficient number of position restraint coordinates");
958 /* Take reference positions from A position of normal posres */
959 ip
->fbposres
.pos0
[XX
] = molb
->posres_xA
[a_molb
][XX
];
960 ip
->fbposres
.pos0
[YY
] = molb
->posres_xA
[a_molb
][YY
];
961 ip
->fbposres
.pos0
[ZZ
] = molb
->posres_xA
[a_molb
][ZZ
];
963 /* Note: no B-type for flat-bottom posres */
965 /* Set the parameter index for idef->iparams_posre */
969 /*! \brief Store a virtual site interaction, complex because of PBC and recursion */
970 static void add_vsite(const gmx_ga2la_t
&ga2la
,
971 gmx::ArrayRef
<const int> index
,
972 gmx::ArrayRef
<const int> rtil
,
974 gmx_bool bHomeA
, int a
, int a_gl
, int a_mol
,
975 const t_iatom
*iatoms
,
979 t_iatom tiatoms
[1+MAXATOMLIST
];
980 int j
, ftype_r
, nral_r
;
982 /* Add this interaction to the local topology */
983 add_ifunc_for_vsites(tiatoms
, ga2la
, nral
, bHomeA
, a
, a_gl
, a_mol
, iatoms
, &idef
->il
[ftype
]);
987 /* Check for recursion */
988 for (k
= 2; k
< 1+nral
; k
++)
990 if ((iatoms
[1+nral
] & (2<<k
)) && (tiatoms
[k
] < 0))
992 /* This construction atoms is a vsite and not a home atom */
995 fprintf(debug
, "Constructing atom %d of vsite atom %d is a vsite and non-home\n", iatoms
[k
]+1, a_mol
+1);
997 /* Find the vsite construction */
999 /* Check all interactions assigned to this atom */
1000 j
= index
[iatoms
[k
]];
1001 while (j
< index
[iatoms
[k
]+1])
1003 ftype_r
= rtil
[j
++];
1004 nral_r
= NRAL(ftype_r
);
1005 if (interaction_function
[ftype_r
].flags
& IF_VSITE
)
1007 /* Add this vsite (recursion) */
1008 add_vsite(ga2la
, index
, rtil
, ftype_r
, nral_r
,
1009 FALSE
, -1, a_gl
+iatoms
[k
]-iatoms
[1], iatoms
[k
],
1013 j
+= 1 + nral_rt(ftype_r
);
1020 /*! \brief Build the index that maps each local atom to its local atom group */
1021 static void makeLocalAtomGroupsFromAtoms(gmx_domdec_t
*dd
)
1023 const gmx::RangePartitioning
&atomGrouping
= dd
->atomGrouping();
1025 dd
->localAtomGroupFromAtom
.clear();
1027 for (size_t g
= 0; g
< dd
->globalAtomGroupIndices
.size(); g
++)
1029 for (int gmx_unused a
: atomGrouping
.block(g
))
1031 dd
->localAtomGroupFromAtom
.push_back(g
);
1036 /*! \brief Returns the squared distance between charge groups \p i and \p j */
1037 static real
dd_dist2(t_pbc
*pbc_null
, rvec
*cg_cm
, const int *la2lc
, int i
, int j
)
1043 pbc_dx_aiuc(pbc_null
, cg_cm
[la2lc
[i
]], cg_cm
[la2lc
[j
]], dx
);
1047 rvec_sub(cg_cm
[la2lc
[i
]], cg_cm
[la2lc
[j
]], dx
);
1053 /*! \brief Append t_blocka block structures 1 to nsrc in src to *dest */
1054 static void combine_blocka(t_blocka
*dest
,
1055 gmx::ArrayRef
<const thread_work_t
> src
)
1057 int ni
= src
.back().excl
.nr
;
1059 for (const thread_work_t
&th_work
: src
)
1061 na
+= th_work
.excl
.nra
;
1063 if (ni
+ 1 > dest
->nalloc_index
)
1065 dest
->nalloc_index
= over_alloc_large(ni
+1);
1066 srenew(dest
->index
, dest
->nalloc_index
);
1068 if (dest
->nra
+ na
> dest
->nalloc_a
)
1070 dest
->nalloc_a
= over_alloc_large(dest
->nra
+na
);
1071 srenew(dest
->a
, dest
->nalloc_a
);
1073 for (gmx::index s
= 1; s
< src
.ssize(); s
++)
1075 for (int i
= dest
->nr
+ 1; i
< src
[s
].excl
.nr
+ 1; i
++)
1077 dest
->index
[i
] = dest
->nra
+ src
[s
].excl
.index
[i
];
1079 for (int i
= 0; i
< src
[s
].excl
.nra
; i
++)
1081 dest
->a
[dest
->nra
+i
] = src
[s
].excl
.a
[i
];
1083 dest
->nr
= src
[s
].excl
.nr
;
1084 dest
->nra
+= src
[s
].excl
.nra
;
1088 /*! \brief Append t_idef structures 1 to nsrc in src to *dest */
1089 static void combine_idef(t_idef
*dest
,
1090 gmx::ArrayRef
<const thread_work_t
> src
)
1094 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
1097 for (gmx::index s
= 1; s
< src
.ssize(); s
++)
1099 n
+= src
[s
].idef
.il
[ftype
].nr
;
1105 ild
= &dest
->il
[ftype
];
1107 if (ild
->nr
+ n
> ild
->nalloc
)
1109 ild
->nalloc
= over_alloc_large(ild
->nr
+n
);
1110 srenew(ild
->iatoms
, ild
->nalloc
);
1113 for (gmx::index s
= 1; s
< src
.ssize(); s
++)
1115 const t_ilist
&ils
= src
[s
].idef
.il
[ftype
];
1117 for (int i
= 0; i
< ils
.nr
; i
++)
1119 ild
->iatoms
[ild
->nr
+ i
] = ils
.iatoms
[i
];
1125 /* Position restraints need an additional treatment */
1126 if (ftype
== F_POSRES
|| ftype
== F_FBPOSRES
)
1128 int nposres
= dest
->il
[ftype
].nr
/2;
1129 // TODO: Simplify this code using std::vector
1130 t_iparams
* &iparams_dest
= (ftype
== F_POSRES
? dest
->iparams_posres
: dest
->iparams_fbposres
);
1131 int &posres_nalloc
= (ftype
== F_POSRES
? dest
->iparams_posres_nalloc
: dest
->iparams_fbposres_nalloc
);
1132 if (nposres
> posres_nalloc
)
1134 posres_nalloc
= over_alloc_large(nposres
);
1135 srenew(iparams_dest
, posres_nalloc
);
1138 /* Set nposres to the number of original position restraints in dest */
1139 for (gmx::index s
= 1; s
< src
.ssize(); s
++)
1141 nposres
-= src
[s
].idef
.il
[ftype
].nr
/2;
1144 for (gmx::index s
= 1; s
< src
.ssize(); s
++)
1146 const t_iparams
*iparams_src
= (ftype
== F_POSRES
? src
[s
].idef
.iparams_posres
: src
[s
].idef
.iparams_fbposres
);
1148 for (int i
= 0; i
< src
[s
].idef
.il
[ftype
].nr
/2; i
++)
1150 /* Correct the index into iparams_posres */
1151 dest
->il
[ftype
].iatoms
[nposres
*2] = nposres
;
1152 /* Copy the position restraint force parameters */
1153 iparams_dest
[nposres
] = iparams_src
[i
];
1162 /*! \brief Check and when available assign bonded interactions for local atom i
1165 check_assign_interactions_atom(int i
, int i_gl
,
1167 int numAtomsInMolecule
,
1168 gmx::ArrayRef
<const int> index
,
1169 gmx::ArrayRef
<const int> rtil
,
1170 gmx_bool bInterMolInteractions
,
1171 int ind_start
, int ind_end
,
1172 const gmx_domdec_t
*dd
,
1173 const gmx_domdec_zones_t
*zones
,
1174 const gmx_molblock_t
*molb
,
1175 gmx_bool bRCheckMB
, const ivec rcheck
, gmx_bool bRCheck2B
,
1180 const t_iparams
*ip_in
,
1191 t_iatom tiatoms
[1 + MAXATOMLIST
];
1193 const int ftype
= rtil
[j
++];
1194 auto iatoms
= gmx::constArrayRefFromArray(rtil
.data() + j
, rtil
.size() - j
);
1195 const int nral
= NRAL(ftype
);
1196 if (interaction_function
[ftype
].flags
& IF_VSITE
)
1198 assert(!bInterMolInteractions
);
1199 /* The vsite construction goes where the vsite itself is */
1202 add_vsite(*dd
->ga2la
, index
, rtil
, ftype
, nral
,
1203 TRUE
, i
, i_gl
, i_mol
,
1204 iatoms
.data(), idef
);
1213 tiatoms
[0] = iatoms
[0];
1217 assert(!bInterMolInteractions
);
1218 /* Assign single-body interactions to the home zone */
1223 if (ftype
== F_POSRES
)
1225 add_posres(mol
, i_mol
, numAtomsInMolecule
,
1226 molb
, tiatoms
, ip_in
, idef
);
1228 else if (ftype
== F_FBPOSRES
)
1230 add_fbposres(mol
, i_mol
, numAtomsInMolecule
,
1231 molb
, tiatoms
, ip_in
, idef
);
1241 /* This is a two-body interaction, we can assign
1242 * analogous to the non-bonded assignments.
1246 if (!bInterMolInteractions
)
1248 /* Get the global index using the offset in the molecule */
1249 k_gl
= i_gl
+ iatoms
[2] - i_mol
;
1255 if (const auto *entry
= dd
->ga2la
->find(k_gl
))
1257 int kz
= entry
->cell
;
1262 /* Check zone interaction assignments */
1263 bUse
= ((iz
< zones
->nizone
&&
1265 kz
>= zones
->izone
[iz
].j0
&&
1266 kz
< zones
->izone
[iz
].j1
) ||
1267 (kz
< zones
->nizone
&&
1269 iz
>= zones
->izone
[kz
].j0
&&
1270 iz
< zones
->izone
[kz
].j1
));
1273 GMX_ASSERT(ftype
!= F_CONSTR
|| (iz
== 0 && kz
== 0),
1274 "Constraint assigned here should only involve home atoms");
1277 tiatoms
[2] = entry
->la
;
1278 /* If necessary check the cgcm distance */
1280 dd_dist2(pbc_null
, cg_cm
, la2lc
,
1281 tiatoms
[1], tiatoms
[2]) >= rc2
)
1294 /* Assign this multi-body bonded interaction to
1295 * the local node if we have all the atoms involved
1296 * (local or communicated) and the minimum zone shift
1297 * in each dimension is zero, for dimensions
1298 * with 2 DD cells an extra check may be necessary.
1300 ivec k_zero
, k_plus
;
1306 for (k
= 1; k
<= nral
&& bUse
; k
++)
1309 if (!bInterMolInteractions
)
1311 /* Get the global index using the offset in the molecule */
1312 k_gl
= i_gl
+ iatoms
[k
] - i_mol
;
1318 const auto *entry
= dd
->ga2la
->find(k_gl
);
1319 if (entry
== nullptr || entry
->cell
>= zones
->n
)
1321 /* We do not have this atom of this interaction
1322 * locally, or it comes from more than one cell
1331 tiatoms
[k
] = entry
->la
;
1332 for (d
= 0; d
< DIM
; d
++)
1334 if (zones
->shift
[entry
->cell
][d
] == 0)
1346 (k_zero
[XX
] != 0) && (k_zero
[YY
] != 0) && (k_zero
[ZZ
] != 0));
1351 for (d
= 0; (d
< DIM
&& bUse
); d
++)
1353 /* Check if the cg_cm distance falls within
1354 * the cut-off to avoid possible multiple
1355 * assignments of bonded interactions.
1359 dd_dist2(pbc_null
, cg_cm
, la2lc
,
1360 tiatoms
[k_zero
[d
]], tiatoms
[k_plus
[d
]]) >= rc2
)
1369 /* Add this interaction to the local topology */
1370 add_ifunc(nral
, tiatoms
, &idef
->il
[ftype
]);
1371 /* Sum so we can check in global_stat
1372 * if we have everything.
1375 !(interaction_function
[ftype
].flags
& IF_LIMZERO
))
1385 /*! \brief This function looks up and assigns bonded interactions for zone iz.
1387 * With thread parallelizing each thread acts on a different atom range:
1388 * at_start to at_end.
1390 static int make_bondeds_zone(gmx_domdec_t
*dd
,
1391 const gmx_domdec_zones_t
*zones
,
1392 const std::vector
<gmx_molblock_t
> &molb
,
1393 gmx_bool bRCheckMB
, ivec rcheck
, gmx_bool bRCheck2B
,
1395 int *la2lc
, t_pbc
*pbc_null
, rvec
*cg_cm
,
1396 const t_iparams
*ip_in
,
1399 gmx::RangePartitioning::Block atomRange
)
1401 int mb
, mt
, mol
, i_mol
;
1403 gmx_reverse_top_t
*rt
;
1406 rt
= dd
->reverse_top
;
1408 bBCheck
= rt
->bBCheck
;
1412 for (int i
: atomRange
)
1414 /* Get the global atom number */
1415 const int i_gl
= dd
->globalAtomIndices
[i
];
1416 global_atomnr_to_moltype_ind(rt
, i_gl
, &mb
, &mt
, &mol
, &i_mol
);
1417 /* Check all intramolecular interactions assigned to this atom */
1418 gmx::ArrayRef
<const int> index
= rt
->ril_mt
[mt
].index
;
1419 gmx::ArrayRef
<const t_iatom
> rtil
= rt
->ril_mt
[mt
].il
;
1421 check_assign_interactions_atom(i
, i_gl
, mol
, i_mol
,
1422 rt
->ril_mt
[mt
].numAtomsInMolecule
,
1424 index
[i_mol
], index
[i_mol
+1],
1427 bRCheckMB
, rcheck
, bRCheck2B
, rc2
,
1438 if (rt
->bIntermolecularInteractions
)
1440 /* Check all intermolecular interactions assigned to this atom */
1441 index
= rt
->ril_intermol
.index
;
1442 rtil
= rt
->ril_intermol
.il
;
1444 check_assign_interactions_atom(i
, i_gl
, mol
, i_mol
,
1445 rt
->ril_mt
[mt
].numAtomsInMolecule
,
1447 index
[i_gl
], index
[i_gl
+ 1],
1450 bRCheckMB
, rcheck
, bRCheck2B
, rc2
,
1462 return nbonded_local
;
1465 /*! \brief Set the exclusion data for i-zone \p iz for the case of no exclusions */
1466 static void set_no_exclusions_zone(const gmx_domdec_t
*dd
,
1467 const gmx_domdec_zones_t
*zones
,
1471 const auto zone
= dd
->atomGrouping().subRange(zones
->cg_range
[iz
],
1472 zones
->cg_range
[iz
+ 1]);
1476 lexcls
->index
[a
+ 1] = lexcls
->nra
;
1480 /*! \brief Set the exclusion data for i-zone \p iz
1482 * This is a legacy version for the group scheme of the same routine below.
1483 * Here charge groups and distance checks to ensure unique exclusions
1486 static int make_exclusions_zone_cg(gmx_domdec_t
*dd
, gmx_domdec_zones_t
*zones
,
1487 const std::vector
<gmx_moltype_t
> &moltype
,
1488 gmx_bool bRCheck
, real rc2
,
1489 int *la2lc
, t_pbc
*pbc_null
, rvec
*cg_cm
,
1493 int cg_start
, int cg_end
)
1497 const t_blocka
*excls
;
1499 const gmx_ga2la_t
&ga2la
= *dd
->ga2la
;
1502 dd
->atomGrouping().subRange(zones
->izone
[iz
].jcg0
,
1503 zones
->izone
[iz
].jcg1
);
1505 n_excl_at_max
= dd
->reverse_top
->n_excl_at_max
;
1507 /* We set the end index, but note that we might not start at zero here */
1508 lexcls
->nr
= dd
->atomGrouping().subRange(0, cg_end
).size();
1510 int n
= lexcls
->nra
;
1512 for (int cg
= cg_start
; cg
< cg_end
; cg
++)
1514 if (n
+ (cg_end
- cg_start
)*n_excl_at_max
> lexcls
->nalloc_a
)
1516 lexcls
->nalloc_a
= over_alloc_large(n
+ (cg_end
- cg_start
)*n_excl_at_max
);
1517 srenew(lexcls
->a
, lexcls
->nalloc_a
);
1519 const auto atomGroup
= dd
->atomGrouping().block(cg
);
1520 if (GET_CGINFO_EXCL_INTER(cginfo
[cg
]) ||
1521 !GET_CGINFO_EXCL_INTRA(cginfo
[cg
]))
1523 /* Copy the exclusions from the global top */
1524 for (int la
: atomGroup
)
1526 lexcls
->index
[la
] = n
;
1527 int a_gl
= dd
->globalAtomIndices
[la
];
1529 global_atomnr_to_moltype_ind(dd
->reverse_top
, a_gl
, &mb
, &mt
, &mol
, &a_mol
);
1530 excls
= &moltype
[mt
].excls
;
1531 for (int j
= excls
->index
[a_mol
]; j
< excls
->index
[a_mol
+1]; j
++)
1533 int aj_mol
= excls
->a
[j
];
1534 /* This computation of jla is only correct intra-cg */
1535 int jla
= la
+ aj_mol
- a_mol
;
1536 if (atomGroup
.inRange(jla
))
1538 /* This is an intra-cg exclusion. We can skip
1539 * the global indexing and distance checking.
1541 /* Intra-cg exclusions are only required
1542 * for the home zone.
1546 lexcls
->a
[n
++] = jla
;
1547 /* Check to avoid double counts */
1556 /* This is a inter-cg exclusion */
1557 /* Since exclusions are pair interactions,
1558 * just like non-bonded interactions,
1559 * they can be assigned properly up
1560 * to the DD cutoff (not cutoff_min as
1561 * for the other bonded interactions).
1563 if (const auto *jEntry
= ga2la
.find(a_gl
+ aj_mol
- a_mol
))
1565 if (iz
== 0 && jEntry
->cell
== 0)
1567 lexcls
->a
[n
++] = jEntry
->la
;
1568 /* Check to avoid double counts */
1569 if (jEntry
->la
> la
)
1574 else if (jRange
.inRange(jEntry
->la
) &&
1576 dd_dist2(pbc_null
, cg_cm
, la2lc
, la
, jEntry
->la
) < rc2
))
1578 /* jla > la, since jRange.begin() > la */
1579 lexcls
->a
[n
++] = jEntry
->la
;
1589 /* There are no inter-cg excls and this cg is self-excluded.
1590 * These exclusions are only required for zone 0,
1591 * since other zones do not see themselves.
1595 for (int la
: atomGroup
)
1597 lexcls
->index
[la
] = n
;
1598 for (int j
: atomGroup
)
1603 count
+= (atomGroup
.size()*(atomGroup
.size() - 1))/2;
1607 /* We don't need exclusions for this cg */
1608 for (int la
: atomGroup
)
1610 lexcls
->index
[la
] = n
;
1616 lexcls
->index
[lexcls
->nr
] = n
;
1622 /*! \brief Set the exclusion data for i-zone \p iz */
1623 static void make_exclusions_zone(gmx_domdec_t
*dd
, gmx_domdec_zones_t
*zones
,
1624 const std::vector
<gmx_moltype_t
> &moltype
,
1625 const int *cginfo
, t_blocka
*lexcls
, int iz
,
1626 int at_start
, int at_end
,
1627 const gmx::ArrayRef
<const int> intermolecularExclusionGroup
)
1629 int n_excl_at_max
, n
, at
;
1631 const gmx_ga2la_t
&ga2la
= *dd
->ga2la
;
1634 dd
->atomGrouping().subRange(zones
->izone
[iz
].jcg0
,
1635 zones
->izone
[iz
].jcg1
);
1637 n_excl_at_max
= dd
->reverse_top
->n_excl_at_max
;
1639 /* We set the end index, but note that we might not start at zero here */
1640 lexcls
->nr
= at_end
;
1643 for (at
= at_start
; at
< at_end
; at
++)
1645 if (n
+ 1000 > lexcls
->nalloc_a
)
1647 lexcls
->nalloc_a
= over_alloc_large(n
+ 1000);
1648 srenew(lexcls
->a
, lexcls
->nalloc_a
);
1651 if (GET_CGINFO_EXCL_INTER(cginfo
[at
]))
1653 int a_gl
, mb
, mt
, mol
, a_mol
, j
;
1654 const t_blocka
*excls
;
1656 if (n
+ n_excl_at_max
> lexcls
->nalloc_a
)
1658 lexcls
->nalloc_a
= over_alloc_large(n
+ n_excl_at_max
);
1659 srenew(lexcls
->a
, lexcls
->nalloc_a
);
1662 /* Copy the exclusions from the global top */
1663 lexcls
->index
[at
] = n
;
1664 a_gl
= dd
->globalAtomIndices
[at
];
1665 global_atomnr_to_moltype_ind(dd
->reverse_top
, a_gl
,
1666 &mb
, &mt
, &mol
, &a_mol
);
1667 excls
= &moltype
[mt
].excls
;
1668 for (j
= excls
->index
[a_mol
]; j
< excls
->index
[a_mol
+ 1]; j
++)
1670 const int aj_mol
= excls
->a
[j
];
1672 if (const auto *jEntry
= ga2la
.find(a_gl
+ aj_mol
- a_mol
))
1674 /* This check is not necessary, but it can reduce
1675 * the number of exclusions in the list, which in turn
1676 * can speed up the pair list construction a bit.
1678 if (jRange
.inRange(jEntry
->la
))
1680 lexcls
->a
[n
++] = jEntry
->la
;
1687 /* We don't need exclusions for this atom */
1688 lexcls
->index
[at
] = n
;
1691 bool isExcludedAtom
= !intermolecularExclusionGroup
.empty() &&
1692 std::find(intermolecularExclusionGroup
.begin(),
1693 intermolecularExclusionGroup
.end(),
1694 dd
->globalAtomIndices
[at
]) !=
1695 intermolecularExclusionGroup
.end();
1699 if (n
+ intermolecularExclusionGroup
.ssize() > lexcls
->nalloc_a
)
1702 over_alloc_large(n
+ intermolecularExclusionGroup
.size());
1703 srenew(lexcls
->a
, lexcls
->nalloc_a
);
1705 for (int qmAtomGlobalIndex
: intermolecularExclusionGroup
)
1707 if (const auto *entry
= dd
->ga2la
->find(qmAtomGlobalIndex
))
1709 lexcls
->a
[n
++] = entry
->la
;
1715 lexcls
->index
[lexcls
->nr
] = n
;
1720 /*! \brief Ensure we have enough space in \p ba for \p nindex_max indices */
1721 static void check_alloc_index(t_blocka
*ba
, int nindex_max
)
1723 if (nindex_max
+1 > ba
->nalloc_index
)
1725 ba
->nalloc_index
= over_alloc_dd(nindex_max
+1);
1726 srenew(ba
->index
, ba
->nalloc_index
);
1730 /*! \brief Ensure that we have enough space for exclusion storate in \p lexcls */
1731 static void check_exclusions_alloc(gmx_domdec_t
*dd
, gmx_domdec_zones_t
*zones
,
1734 int nr
= dd
->atomGrouping().subRange(0, zones
->izone
[zones
->nizone
- 1].cg1
).size();
1736 check_alloc_index(lexcls
, nr
);
1738 for (size_t thread
= 1; thread
< dd
->reverse_top
->th_work
.size(); thread
++)
1740 check_alloc_index(&dd
->reverse_top
->th_work
[thread
].excl
, nr
);
1744 /*! \brief Set the total count indexes for the local exclusions, needed by several functions */
1745 static void finish_local_exclusions(gmx_domdec_t
*dd
, gmx_domdec_zones_t
*zones
,
1748 const auto nonhomeIzonesAtomRange
=
1749 dd
->atomGrouping().subRange(zones
->izone
[0].cg1
,
1750 zones
->izone
[zones
->nizone
- 1].cg1
);
1752 if (dd
->n_intercg_excl
== 0)
1754 /* There are no exclusions involving non-home charge groups,
1755 * but we need to set the indices for neighborsearching.
1757 for (int la
: nonhomeIzonesAtomRange
)
1759 lexcls
->index
[la
] = lexcls
->nra
;
1762 /* nr is only used to loop over the exclusions for Ewald and RF,
1763 * so we can set it to the number of home atoms for efficiency.
1765 lexcls
->nr
= nonhomeIzonesAtomRange
.begin();
1769 lexcls
->nr
= nonhomeIzonesAtomRange
.end();
1773 /*! \brief Clear a t_idef struct */
1774 static void clear_idef(t_idef
*idef
)
1778 /* Clear the counts */
1779 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
1781 idef
->il
[ftype
].nr
= 0;
1785 /*! \brief Generate and store all required local bonded interactions in \p idef and local exclusions in \p lexcls */
1786 static int make_local_bondeds_excls(gmx_domdec_t
*dd
,
1787 gmx_domdec_zones_t
*zones
,
1788 const gmx_mtop_t
*mtop
,
1790 gmx_bool bRCheckMB
, ivec rcheck
, gmx_bool bRCheck2B
,
1792 int *la2lc
, t_pbc
*pbc_null
, rvec
*cg_cm
,
1794 t_blocka
*lexcls
, int *excl_count
)
1796 int nzone_bondeds
, nzone_excl
;
1797 int izone
, cg0
, cg1
;
1801 gmx_reverse_top_t
*rt
;
1803 if (dd
->reverse_top
->bInterCGInteractions
)
1805 nzone_bondeds
= zones
->n
;
1809 /* Only single charge group (or atom) molecules, so interactions don't
1810 * cross zone boundaries and we only need to assign in the home zone.
1815 if (dd
->n_intercg_excl
> 0)
1817 /* We only use exclusions from i-zones to i- and j-zones */
1818 nzone_excl
= zones
->nizone
;
1822 /* There are no inter-cg exclusions and only zone 0 sees itself */
1826 check_exclusions_alloc(dd
, zones
, lexcls
);
1828 rt
= dd
->reverse_top
;
1832 /* Clear the counts */
1840 for (izone
= 0; izone
< nzone_bondeds
; izone
++)
1842 cg0
= zones
->cg_range
[izone
];
1843 cg1
= zones
->cg_range
[izone
+ 1];
1845 const int numThreads
= rt
->th_work
.size();
1846 #pragma omp parallel for num_threads(numThreads) schedule(static)
1847 for (thread
= 0; thread
< numThreads
; thread
++)
1855 cg0t
= cg0
+ ((cg1
- cg0
)* thread
)/numThreads
;
1856 cg1t
= cg0
+ ((cg1
- cg0
)*(thread
+1))/numThreads
;
1864 idef_t
= &rt
->th_work
[thread
].idef
;
1868 rt
->th_work
[thread
].nbonded
=
1869 make_bondeds_zone(dd
, zones
,
1871 bRCheckMB
, rcheck
, bRCheck2B
, rc2
,
1872 la2lc
, pbc_null
, cg_cm
, idef
->iparams
,
1875 dd
->atomGrouping().subRange(cg0t
, cg1t
));
1877 if (izone
< nzone_excl
)
1885 excl_t
= &rt
->th_work
[thread
].excl
;
1890 if (dd
->atomGrouping().allBlocksHaveSizeOne() &&
1893 /* No charge groups and no distance check required */
1894 make_exclusions_zone(dd
, zones
, mtop
->moltype
, cginfo
,
1895 excl_t
, izone
, cg0t
,
1897 mtop
->intermolecularExclusionGroup
);
1901 rt
->th_work
[thread
].excl_count
=
1902 make_exclusions_zone_cg(dd
, zones
,
1903 mtop
->moltype
, bRCheck2B
, rc2
,
1904 la2lc
, pbc_null
, cg_cm
, cginfo
,
1911 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
1914 if (rt
->th_work
.size() > 1)
1916 combine_idef(idef
, rt
->th_work
);
1919 for (const thread_work_t
&th_work
: rt
->th_work
)
1921 nbonded_local
+= th_work
.nbonded
;
1924 if (izone
< nzone_excl
)
1926 if (rt
->th_work
.size() > 1)
1928 combine_blocka(lexcls
, rt
->th_work
);
1931 for (const thread_work_t
&th_work
: rt
->th_work
)
1933 *excl_count
+= th_work
.excl_count
;
1938 /* Some zones might not have exclusions, but some code still needs to
1939 * loop over the index, so we set the indices here.
1941 for (izone
= nzone_excl
; izone
< zones
->nizone
; izone
++)
1943 set_no_exclusions_zone(dd
, zones
, izone
, lexcls
);
1946 finish_local_exclusions(dd
, zones
, lexcls
);
1949 fprintf(debug
, "We have %d exclusions, check count %d\n",
1950 lexcls
->nra
, *excl_count
);
1953 return nbonded_local
;
1956 void dd_make_local_cgs(gmx_domdec_t
*dd
, t_block
*lcgs
)
1958 lcgs
->nr
= dd
->globalAtomGroupIndices
.size();
1959 lcgs
->index
= dd
->atomGrouping_
.rawIndex().data();
1962 void dd_make_local_top(gmx_domdec_t
*dd
, gmx_domdec_zones_t
*zones
,
1963 int npbcdim
, matrix box
,
1964 rvec cellsize_min
, const ivec npulse
,
1967 const gmx_mtop_t
&mtop
, gmx_localtop_t
*ltop
)
1969 gmx_bool bRCheckMB
, bRCheck2B
;
1973 t_pbc pbc
, *pbc_null
= nullptr;
1977 fprintf(debug
, "Making local topology\n");
1980 dd_make_local_cgs(dd
, <op
->cgs
);
1985 if (dd
->reverse_top
->bInterCGInteractions
)
1987 /* We need to check to which cell bondeds should be assigned */
1988 rc
= dd_cutoff_twobody(dd
);
1991 fprintf(debug
, "Two-body bonded cut-off distance is %g\n", rc
);
1994 /* Should we check cg_cm distances when assigning bonded interactions? */
1995 for (d
= 0; d
< DIM
; d
++)
1998 /* Only need to check for dimensions where the part of the box
1999 * that is not communicated is smaller than the cut-off.
2001 if (d
< npbcdim
&& dd
->nc
[d
] > 1 &&
2002 (dd
->nc
[d
] - npulse
[d
])*cellsize_min
[d
] < 2*rc
)
2009 /* Check for interactions between two atoms,
2010 * where we can allow interactions up to the cut-off,
2011 * instead of up to the smallest cell dimension.
2018 "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %s\n",
2019 d
, cellsize_min
[d
], d
, rcheck
[d
], gmx::boolToString(bRCheck2B
));
2022 if (bRCheckMB
|| bRCheck2B
)
2024 makeLocalAtomGroupsFromAtoms(dd
);
2027 pbc_null
= set_pbc_dd(&pbc
, fr
->ePBC
, dd
->nc
, TRUE
, box
);
2037 make_local_bondeds_excls(dd
, zones
, &mtop
, fr
->cginfo
,
2038 bRCheckMB
, rcheck
, bRCheck2B
, rc
,
2039 dd
->localAtomGroupFromAtom
.data(),
2040 pbc_null
, cgcm_or_x
,
2042 <op
->excls
, &nexcl
);
2044 /* The ilist is not sorted yet,
2045 * we can only do this when we have the charge arrays.
2047 ltop
->idef
.ilsort
= ilsortUNKNOWN
;
2049 if (dd
->reverse_top
->bExclRequired
)
2051 dd
->nbonded_local
+= nexcl
;
2054 ltop
->atomtypes
= mtop
.atomtypes
;
2057 void dd_sort_local_top(gmx_domdec_t
*dd
, const t_mdatoms
*mdatoms
,
2058 gmx_localtop_t
*ltop
)
2060 if (dd
->reverse_top
->ilsort
== ilsortNO_FE
)
2062 ltop
->idef
.ilsort
= ilsortNO_FE
;
2066 gmx_sort_ilist_fe(<op
->idef
, mdatoms
->chargeA
, mdatoms
->chargeB
);
2070 void dd_init_local_top(const gmx_mtop_t
&top_global
,
2071 gmx_localtop_t
*top
)
2073 /* TODO: Get rid of the const casts below, e.g. by using a reference */
2074 top
->idef
.ntypes
= top_global
.ffparams
.numTypes();
2075 top
->idef
.atnr
= top_global
.ffparams
.atnr
;
2076 top
->idef
.functype
= const_cast<t_functype
*>(top_global
.ffparams
.functype
.data());
2077 top
->idef
.iparams
= const_cast<t_iparams
*>(top_global
.ffparams
.iparams
.data());
2078 top
->idef
.fudgeQQ
= top_global
.ffparams
.fudgeQQ
;
2079 top
->idef
.cmap_grid
= new gmx_cmap_t
;
2080 *top
->idef
.cmap_grid
= top_global
.ffparams
.cmap_grid
;
2082 top
->idef
.ilsort
= ilsortUNKNOWN
;
2083 top
->useInDomainDecomp_
= true;
2086 void dd_init_local_state(gmx_domdec_t
*dd
,
2087 const t_state
*state_global
, t_state
*state_local
)
2089 int buf
[NITEM_DD_INIT_LOCAL_STATE
];
2093 buf
[0] = state_global
->flags
;
2094 buf
[1] = state_global
->ngtc
;
2095 buf
[2] = state_global
->nnhpres
;
2096 buf
[3] = state_global
->nhchainlength
;
2097 buf
[4] = state_global
->dfhist
? state_global
->dfhist
->nlambda
: 0;
2099 dd_bcast(dd
, NITEM_DD_INIT_LOCAL_STATE
*sizeof(int), buf
);
2101 init_gtc_state(state_local
, buf
[1], buf
[2], buf
[3]);
2102 init_dfhist_state(state_local
, buf
[4]);
2103 state_local
->flags
= buf
[0];
2106 /*! \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 */
2107 static void check_link(t_blocka
*link
, int cg_gl
, int cg_gl_j
)
2113 for (k
= link
->index
[cg_gl
]; k
< link
->index
[cg_gl
+1]; k
++)
2115 GMX_RELEASE_ASSERT(link
->a
, "Inconsistent NULL pointer while making charge-group links");
2116 if (link
->a
[k
] == cg_gl_j
)
2123 GMX_RELEASE_ASSERT(link
->a
|| link
->index
[cg_gl
+1]+1 > link
->nalloc_a
,
2124 "Inconsistent allocation of link");
2125 /* Add this charge group link */
2126 if (link
->index
[cg_gl
+1]+1 > link
->nalloc_a
)
2128 link
->nalloc_a
= over_alloc_large(link
->index
[cg_gl
+1]+1);
2129 srenew(link
->a
, link
->nalloc_a
);
2131 link
->a
[link
->index
[cg_gl
+1]] = cg_gl_j
;
2132 link
->index
[cg_gl
+1]++;
2136 /*! \brief Return a vector of the charge group index for all atoms */
2137 static std::vector
<int> make_at2cg(const t_block
&cgs
)
2139 std::vector
<int> at2cg(cgs
.index
[cgs
.nr
]);
2140 for (int cg
= 0; cg
< cgs
.nr
; cg
++)
2142 for (int a
= cgs
.index
[cg
]; a
< cgs
.index
[cg
+ 1]; a
++)
2151 t_blocka
*make_charge_group_links(const gmx_mtop_t
*mtop
, gmx_domdec_t
*dd
,
2152 cginfo_mb_t
*cginfo_mb
)
2154 gmx_bool bExclRequired
;
2156 cginfo_mb_t
*cgi_mb
;
2158 /* For each charge group make a list of other charge groups
2159 * in the system that a linked to it via bonded interactions
2160 * which are also stored in reverse_top.
2163 bExclRequired
= dd
->reverse_top
->bExclRequired
;
2165 reverse_ilist_t ril_intermol
;
2166 if (mtop
->bIntermolecularInteractions
)
2168 if (ncg_mtop(mtop
) < mtop
->natoms
)
2170 gmx_fatal(FARGS
, "The combination of intermolecular interactions, charge groups and domain decomposition is not supported. Use cutoff-scheme=Verlet (which removes the charge groups) or run without domain decomposition.");
2175 atoms
.nr
= mtop
->natoms
;
2176 atoms
.atom
= nullptr;
2178 GMX_RELEASE_ASSERT(mtop
->intermolecular_ilist
, "We should have an ilist when intermolecular interactions are on");
2180 make_reverse_ilist(*mtop
->intermolecular_ilist
,
2182 FALSE
, FALSE
, FALSE
, TRUE
, &ril_intermol
);
2186 snew(link
->index
, ncg_mtop(mtop
)+1);
2193 for (size_t mb
= 0; mb
< mtop
->molblock
.size(); mb
++)
2195 const gmx_molblock_t
&molb
= mtop
->molblock
[mb
];
2200 const gmx_moltype_t
&molt
= mtop
->moltype
[molb
.type
];
2201 const t_block
&cgs
= molt
.cgs
;
2202 const t_blocka
&excls
= molt
.excls
;
2203 std::vector
<int> a2c
= make_at2cg(cgs
);
2204 /* Make a reverse ilist in which the interactions are linked
2205 * to all atoms, not only the first atom as in gmx_reverse_top.
2206 * The constraints are discarded here.
2208 reverse_ilist_t ril
;
2209 make_reverse_ilist(molt
.ilist
, &molt
.atoms
,
2210 FALSE
, FALSE
, FALSE
, TRUE
, &ril
);
2212 cgi_mb
= &cginfo_mb
[mb
];
2215 for (mol
= 0; mol
< (mtop
->bIntermolecularInteractions
? molb
.nmol
: 1); mol
++)
2217 for (int cg
= 0; cg
< cgs
.nr
; cg
++)
2219 int cg_gl
= cg_offset
+ cg
;
2220 link
->index
[cg_gl
+1] = link
->index
[cg_gl
];
2221 for (int a
= cgs
.index
[cg
]; a
< cgs
.index
[cg
+ 1]; a
++)
2223 int i
= ril
.index
[a
];
2224 while (i
< ril
.index
[a
+1])
2226 int ftype
= ril
.il
[i
++];
2227 int nral
= NRAL(ftype
);
2228 /* Skip the ifunc index */
2230 for (int j
= 0; j
< nral
; j
++)
2232 int aj
= ril
.il
[i
+ j
];
2235 check_link(link
, cg_gl
, cg_offset
+a2c
[aj
]);
2238 i
+= nral_rt(ftype
);
2242 /* Exclusions always go both ways */
2243 for (int j
= excls
.index
[a
]; j
< excls
.index
[a
+ 1]; j
++)
2245 int aj
= excls
.a
[j
];
2248 check_link(link
, cg_gl
, cg_offset
+a2c
[aj
]);
2253 if (mtop
->bIntermolecularInteractions
)
2255 int i
= ril_intermol
.index
[a
];
2256 while (i
< ril_intermol
.index
[a
+1])
2258 int ftype
= ril_intermol
.il
[i
++];
2259 int nral
= NRAL(ftype
);
2260 /* Skip the ifunc index */
2262 for (int j
= 0; j
< nral
; j
++)
2264 /* Here we assume we have no charge groups;
2265 * this has been checked above.
2267 int aj
= ril_intermol
.il
[i
+ j
];
2268 check_link(link
, cg_gl
, aj
);
2270 i
+= nral_rt(ftype
);
2274 if (link
->index
[cg_gl
+1] - link
->index
[cg_gl
] > 0)
2276 SET_CGINFO_BOND_INTER(cgi_mb
->cginfo
[cg
]);
2281 cg_offset
+= cgs
.nr
;
2283 int nlink_mol
= link
->index
[cg_offset
] - link
->index
[cg_offset
- cgs
.nr
];
2287 fprintf(debug
, "molecule type '%s' %d cgs has %d cg links through bonded interac.\n", *molt
.name
, cgs
.nr
, nlink_mol
);
2290 if (molb
.nmol
> mol
)
2292 /* Copy the data for the rest of the molecules in this block */
2293 link
->nalloc_a
+= (molb
.nmol
- mol
)*nlink_mol
;
2294 srenew(link
->a
, link
->nalloc_a
);
2295 for (; mol
< molb
.nmol
; mol
++)
2297 for (int cg
= 0; cg
< cgs
.nr
; cg
++)
2299 int cg_gl
= cg_offset
+ cg
;
2300 link
->index
[cg_gl
+ 1] =
2301 link
->index
[cg_gl
+ 1 - cgs
.nr
] + nlink_mol
;
2302 for (int j
= link
->index
[cg_gl
]; j
< link
->index
[cg_gl
+1]; j
++)
2304 link
->a
[j
] = link
->a
[j
- nlink_mol
] + cgs
.nr
;
2306 if (link
->index
[cg_gl
+1] - link
->index
[cg_gl
] > 0 &&
2307 cg_gl
- cgi_mb
->cg_start
< cgi_mb
->cg_mod
)
2309 SET_CGINFO_BOND_INTER(cgi_mb
->cginfo
[cg_gl
- cgi_mb
->cg_start
]);
2313 cg_offset
+= cgs
.nr
;
2320 fprintf(debug
, "Of the %d charge groups %d are linked via bonded interactions\n", ncg_mtop(mtop
), ncgi
);
2331 } bonded_distance_t
;
2333 /*! \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 */
2334 static void update_max_bonded_distance(real r2
, int ftype
, int a1
, int a2
,
2335 bonded_distance_t
*bd
)
2346 /*! \brief Set the distance, function type and atom indices for the longest distance between charge-groups of molecule type \p molt for two-body and multi-body bonded interactions */
2347 static void bonded_cg_distance_mol(const gmx_moltype_t
*molt
,
2348 const std::vector
<int> &at2cg
,
2349 gmx_bool bBCheck
, gmx_bool bExcl
, rvec
*cg_cm
,
2350 bonded_distance_t
*bd_2b
,
2351 bonded_distance_t
*bd_mb
)
2353 for (int ftype
= 0; ftype
< F_NRE
; ftype
++)
2355 if (dd_check_ftype(ftype
, bBCheck
, FALSE
, FALSE
))
2357 const auto &il
= molt
->ilist
[ftype
];
2358 int nral
= NRAL(ftype
);
2361 for (int i
= 0; i
< il
.size(); i
+= 1+nral
)
2363 for (int ai
= 0; ai
< nral
; ai
++)
2365 int cgi
= at2cg
[il
.iatoms
[i
+1+ai
]];
2366 for (int aj
= ai
+ 1; aj
< nral
; aj
++)
2368 int cgj
= at2cg
[il
.iatoms
[i
+1+aj
]];
2371 real rij2
= distance2(cg_cm
[cgi
], cg_cm
[cgj
]);
2373 update_max_bonded_distance(rij2
, ftype
,
2376 (nral
== 2) ? bd_2b
: bd_mb
);
2386 const t_blocka
*excls
= &molt
->excls
;
2387 for (int ai
= 0; ai
< excls
->nr
; ai
++)
2389 int cgi
= at2cg
[ai
];
2390 for (int j
= excls
->index
[ai
]; j
< excls
->index
[ai
+1]; j
++)
2392 int cgj
= at2cg
[excls
->a
[j
]];
2395 real rij2
= distance2(cg_cm
[cgi
], cg_cm
[cgj
]);
2397 /* There is no function type for exclusions, use -1 */
2398 update_max_bonded_distance(rij2
, -1, ai
, excls
->a
[j
], bd_2b
);
2405 /*! \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 */
2406 static void bonded_distance_intermol(const InteractionLists
&ilists_intermol
,
2408 const rvec
*x
, int ePBC
, const matrix box
,
2409 bonded_distance_t
*bd_2b
,
2410 bonded_distance_t
*bd_mb
)
2414 set_pbc(&pbc
, ePBC
, box
);
2416 for (int ftype
= 0; ftype
< F_NRE
; ftype
++)
2418 if (dd_check_ftype(ftype
, bBCheck
, FALSE
, FALSE
))
2420 const auto &il
= ilists_intermol
[ftype
];
2421 int nral
= NRAL(ftype
);
2423 /* No nral>1 check here, since intermol interactions always
2424 * have nral>=2 (and the code is also correct for nral=1).
2426 for (int i
= 0; i
< il
.size(); i
+= 1+nral
)
2428 for (int ai
= 0; ai
< nral
; ai
++)
2430 int atom_i
= il
.iatoms
[i
+ 1 + ai
];
2432 for (int aj
= ai
+ 1; aj
< nral
; aj
++)
2437 int atom_j
= il
.iatoms
[i
+ 1 + aj
];
2439 pbc_dx(&pbc
, x
[atom_i
], x
[atom_j
], dx
);
2443 update_max_bonded_distance(rij2
, ftype
,
2445 (nral
== 2) ? bd_2b
: bd_mb
);
2453 //! Returns whether \p molt has at least one virtual site
2454 static bool moltypeHasVsite(const gmx_moltype_t
&molt
)
2456 bool hasVsite
= false;
2457 for (int i
= 0; i
< F_NRE
; i
++)
2459 if ((interaction_function
[i
].flags
& IF_VSITE
) &&
2460 molt
.ilist
[i
].size() > 0)
2469 //! Compute charge group centers of mass for molecule \p molt
2470 static void get_cgcm_mol(const gmx_moltype_t
*molt
,
2471 const gmx_ffparams_t
*ffparams
,
2472 int ePBC
, t_graph
*graph
, const matrix box
,
2473 const rvec
*x
, rvec
*xs
, rvec
*cg_cm
)
2477 if (ePBC
!= epbcNONE
)
2479 mk_mshift(nullptr, graph
, ePBC
, box
, x
);
2481 shift_x(graph
, box
, x
, xs
);
2482 /* By doing an extra mk_mshift the molecules that are broken
2483 * because they were e.g. imported from another software
2484 * will be made whole again. Such are the healing powers
2487 mk_mshift(nullptr, graph
, ePBC
, box
, xs
);
2491 /* We copy the coordinates so the original coordinates remain
2492 * unchanged, just to be 100% sure that we do not affect
2493 * binary reproducibility of simulations.
2495 n
= molt
->cgs
.index
[molt
->cgs
.nr
];
2496 for (i
= 0; i
< n
; i
++)
2498 copy_rvec(x
[i
], xs
[i
]);
2502 if (moltypeHasVsite(*molt
))
2504 /* Convert to old, deprecated format */
2505 t_ilist ilist
[F_NRE
];
2506 for (int ftype
= 0; ftype
< F_NRE
; ftype
++)
2508 if (interaction_function
[ftype
].flags
& IF_VSITE
)
2510 ilist
[ftype
].nr
= molt
->ilist
[ftype
].size();
2511 ilist
[ftype
].iatoms
= const_cast<int *>(molt
->ilist
[ftype
].iatoms
.data());
2515 construct_vsites(nullptr, xs
, 0.0, nullptr,
2516 ffparams
->iparams
.data(), ilist
,
2517 epbcNONE
, TRUE
, nullptr, nullptr);
2520 calc_cgcm(nullptr, 0, molt
->cgs
.nr
, &molt
->cgs
, xs
, cg_cm
);
2523 void dd_bonded_cg_distance(const gmx::MDLogger
&mdlog
,
2524 const gmx_mtop_t
*mtop
,
2525 const t_inputrec
*ir
,
2526 const rvec
*x
, const matrix box
,
2528 real
*r_2b
, real
*r_mb
)
2530 gmx_bool bExclRequired
;
2534 bonded_distance_t bd_2b
= { 0, -1, -1, -1 };
2535 bonded_distance_t bd_mb
= { 0, -1, -1, -1 };
2537 bExclRequired
= inputrecExclForces(ir
);
2542 for (const gmx_molblock_t
&molb
: mtop
->molblock
)
2544 const gmx_moltype_t
&molt
= mtop
->moltype
[molb
.type
];
2545 if (molt
.cgs
.nr
== 1 || molb
.nmol
== 0)
2547 at_offset
+= molb
.nmol
*molt
.atoms
.nr
;
2551 if (ir
->ePBC
!= epbcNONE
)
2553 mk_graph_moltype(molt
, &graph
);
2556 std::vector
<int> at2cg
= make_at2cg(molt
.cgs
);
2557 snew(xs
, molt
.atoms
.nr
);
2558 snew(cg_cm
, molt
.cgs
.nr
);
2559 for (int mol
= 0; mol
< molb
.nmol
; mol
++)
2561 get_cgcm_mol(&molt
, &mtop
->ffparams
, ir
->ePBC
, &graph
, box
,
2562 x
+at_offset
, xs
, cg_cm
);
2564 bonded_distance_t bd_mol_2b
= { 0, -1, -1, -1 };
2565 bonded_distance_t bd_mol_mb
= { 0, -1, -1, -1 };
2567 bonded_cg_distance_mol(&molt
, at2cg
, bBCheck
, bExclRequired
, cg_cm
,
2568 &bd_mol_2b
, &bd_mol_mb
);
2570 /* Process the mol data adding the atom index offset */
2571 update_max_bonded_distance(bd_mol_2b
.r2
, bd_mol_2b
.ftype
,
2572 at_offset
+ bd_mol_2b
.a1
,
2573 at_offset
+ bd_mol_2b
.a2
,
2575 update_max_bonded_distance(bd_mol_mb
.r2
, bd_mol_mb
.ftype
,
2576 at_offset
+ bd_mol_mb
.a1
,
2577 at_offset
+ bd_mol_mb
.a2
,
2580 at_offset
+= molt
.atoms
.nr
;
2584 if (ir
->ePBC
!= epbcNONE
)
2591 if (mtop
->bIntermolecularInteractions
)
2593 if (ncg_mtop(mtop
) < mtop
->natoms
)
2595 gmx_fatal(FARGS
, "The combination of intermolecular interactions, charge groups and domain decomposition is not supported. Use cutoff-scheme=Verlet (which removes the charge groups) or run without domain decomposition.");
2598 GMX_RELEASE_ASSERT(mtop
->intermolecular_ilist
, "We should have an ilist when intermolecular interactions are on");
2600 bonded_distance_intermol(*mtop
->intermolecular_ilist
,
2606 *r_2b
= sqrt(bd_2b
.r2
);
2607 *r_mb
= sqrt(bd_mb
.r2
);
2609 if (*r_2b
> 0 || *r_mb
> 0)
2611 GMX_LOG(mdlog
.info
).appendText("Initial maximum distances in bonded interactions:");
2614 GMX_LOG(mdlog
.info
).appendTextFormatted(
2615 " two-body bonded interactions: %5.3f nm, %s, atoms %d %d",
2616 *r_2b
, (bd_2b
.ftype
>= 0) ? interaction_function
[bd_2b
.ftype
].longname
: "Exclusion",
2617 bd_2b
.a1
+ 1, bd_2b
.a2
+ 1);
2621 GMX_LOG(mdlog
.info
).appendTextFormatted(
2622 " multi-body bonded interactions: %5.3f nm, %s, atoms %d %d",
2623 *r_mb
, interaction_function
[bd_mb
.ftype
].longname
,
2624 bd_mb
.a1
+ 1, bd_mb
.a2
+ 1);