2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006,2007,2008,2009,2010,2011,2012,2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
38 * \brief This file defines functions used by the domdec module
39 * while managing the construction, use and error checking for
40 * topologies local to a DD rank.
42 * \author Berk Hess <hess@kth.se>
43 * \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/mdlib/vsite.h"
66 #include "gromacs/mdtypes/commrec.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/smalloc.h"
81 #include "gromacs/utility/stringutil.h"
83 #include "domdec_constraints.h"
84 #include "domdec_internal.h"
85 #include "domdec_vsite.h"
87 /*! \brief The number of integer item in the local state, used for broadcasting of the state */
88 #define NITEM_DD_INIT_LOCAL_STATE 5
90 struct reverse_ilist_t
92 int *index
; /* Index for each atom into il */
93 int *il
; /* ftype|type|a0|...|an|ftype|... */
94 int numAtomsInMolecule
; /* The number of atoms in this molecule */
104 /*! \brief Struct for thread local work data for local topology generation */
106 t_idef idef
; /**< Partial local topology */
107 int **vsite_pbc
; /**< vsite PBC structure */
108 int *vsite_pbc_nalloc
; /**< Allocation sizes for vsite_pbc */
109 int nbonded
; /**< The number of bondeds in this struct */
110 t_blocka excl
; /**< List of exclusions */
111 int excl_count
; /**< The total exclusion count for \p excl */
114 /*! \brief Struct for the reverse topology: links bonded interactions to atomsx */
115 struct gmx_reverse_top_t
117 //! @cond Doxygen_Suppress
118 gmx_bool bExclRequired
; /**< Do we require all exclusions to be assigned? */
119 int n_excl_at_max
; /**< The maximum number of exclusions one atom can have */
120 gmx_bool bConstr
; /**< Are there constraints in this revserse top? */
121 gmx_bool bSettle
; /**< Are there settles in this revserse top? */
122 gmx_bool bBCheck
; /**< All bonded interactions have to be assigned? */
123 gmx_bool bInterCGInteractions
; /**< Are there bondeds/exclusions between charge-groups? */
124 reverse_ilist_t
*ril_mt
; /**< Reverse ilist for all moltypes */
125 int ril_mt_tot_size
; /**< The size of ril_mt[?].index summed over all entries */
126 int ilsort
; /**< The sorting state of bondeds for free energy */
127 molblock_ind_t
*mbi
; /**< molblock to global atom index for quick lookup of molblocks on atom index */
128 int nmolblock
; /**< The number of molblocks, size of \p mbi */
130 gmx_bool bIntermolecularInteractions
; /**< Do we have intermolecular interactions? */
131 reverse_ilist_t ril_intermol
; /**< Intermolecular reverse ilist */
133 /* Work data structures for multi-threading */
134 int nthread
; /**< The number of threads to be used */
135 thread_work_t
*th_work
; /**< Thread work array for local topology generation */
139 /*! \brief Returns the number of atom entries for il in gmx_reverse_top_t */
140 static int nral_rt(int ftype
)
145 if (interaction_function
[ftype
].flags
& IF_VSITE
)
147 /* With vsites the reverse topology contains
148 * two extra entries for PBC.
156 /*! \brief Return whether interactions of type \p ftype need to be assigned exactly once */
157 static gmx_bool
dd_check_ftype(int ftype
, gmx_bool bBCheck
,
158 gmx_bool bConstr
, gmx_bool bSettle
)
160 return (((interaction_function
[ftype
].flags
& IF_BOND
) &&
161 !(interaction_function
[ftype
].flags
& IF_VSITE
) &&
162 (bBCheck
|| !(interaction_function
[ftype
].flags
& IF_LIMZERO
))) ||
163 (bConstr
&& (ftype
== F_CONSTR
|| ftype
== F_CONSTRNC
)) ||
164 (bSettle
&& ftype
== F_SETTLE
));
167 /*! \brief Print a header on error messages */
168 static void print_error_header(FILE *fplog
, const char *moltypename
, int nprint
)
170 fprintf(fplog
, "\nMolecule type '%s'\n", moltypename
);
171 fprintf(stderr
, "\nMolecule type '%s'\n", moltypename
);
173 "the first %d missing interactions, except for exclusions:\n",
176 "the first %d missing interactions, except for exclusions:\n",
180 /*! \brief Help print error output when interactions are missing */
181 static void print_missing_interactions_mb(FILE *fplog
, t_commrec
*cr
,
182 const gmx_reverse_top_t
*rt
,
183 const char *moltypename
,
184 const reverse_ilist_t
*ril
,
185 int a_start
, int a_end
,
186 int nat_mol
, int nmol
,
190 int nril_mol
= ril
->index
[nat_mol
];
191 snew(assigned
, nmol
*nril_mol
);
193 gmx::ArrayRef
<const int> gatindex
= cr
->dd
->globalAtomIndices
;
194 for (int ftype
= 0; ftype
< F_NRE
; ftype
++)
196 if (dd_check_ftype(ftype
, rt
->bBCheck
, rt
->bConstr
, rt
->bSettle
))
198 int nral
= NRAL(ftype
);
199 const t_ilist
*il
= &idef
->il
[ftype
];
200 const t_iatom
*ia
= il
->iatoms
;
201 for (int i
= 0; i
< il
->nr
; i
+= 1+nral
)
203 int a0
= gatindex
[ia
[1]];
204 /* Check if this interaction is in
205 * the currently checked molblock.
207 if (a0
>= a_start
&& a0
< a_end
)
209 int mol
= (a0
- a_start
)/nat_mol
;
210 int a0_mol
= (a0
- a_start
) - mol
*nat_mol
;
211 int j_mol
= ril
->index
[a0_mol
];
213 while (j_mol
< ril
->index
[a0_mol
+1] && !found
)
215 int j
= mol
*nril_mol
+ j_mol
;
216 int ftype_j
= ril
->il
[j_mol
];
217 /* Here we need to check if this interaction has
218 * not already been assigned, since we could have
219 * multiply defined interactions.
221 if (ftype
== ftype_j
&& ia
[0] == ril
->il
[j_mol
+1] &&
224 /* Check the atoms */
226 for (int a
= 0; a
< nral
; a
++)
228 if (gatindex
[ia
[1+a
]] !=
229 a_start
+ mol
*nat_mol
+ ril
->il
[j_mol
+2+a
])
239 j_mol
+= 2 + nral_rt(ftype_j
);
243 gmx_incons("Some interactions seem to be assigned multiple times");
251 gmx_sumi(nmol
*nril_mol
, assigned
, cr
);
255 for (int mol
= 0; mol
< nmol
; mol
++)
258 while (j_mol
< nril_mol
)
260 int ftype
= ril
->il
[j_mol
];
261 int nral
= NRAL(ftype
);
262 int j
= mol
*nril_mol
+ j_mol
;
263 if (assigned
[j
] == 0 &&
264 !(interaction_function
[ftype
].flags
& IF_VSITE
))
266 if (DDMASTER(cr
->dd
))
270 print_error_header(fplog
, moltypename
, nprint
);
272 fprintf(fplog
, "%20s atoms",
273 interaction_function
[ftype
].longname
);
274 fprintf(stderr
, "%20s atoms",
275 interaction_function
[ftype
].longname
);
277 for (a
= 0; a
< nral
; a
++)
279 fprintf(fplog
, "%5d", ril
->il
[j_mol
+2+a
]+1);
280 fprintf(stderr
, "%5d", ril
->il
[j_mol
+2+a
]+1);
285 fprintf(stderr
, " ");
288 fprintf(fplog
, " global");
289 fprintf(stderr
, " global");
290 for (a
= 0; a
< nral
; a
++)
292 fprintf(fplog
, "%6d",
293 a_start
+mol
*nat_mol
+ril
->il
[j_mol
+2+a
]+1);
294 fprintf(stderr
, "%6d",
295 a_start
+mol
*nat_mol
+ril
->il
[j_mol
+2+a
]+1);
297 fprintf(fplog
, "\n");
298 fprintf(stderr
, "\n");
306 j_mol
+= 2 + nral_rt(ftype
);
313 /*! \brief Help print error output when interactions are missing */
314 static void print_missing_interactions_atoms(FILE *fplog
, t_commrec
*cr
,
315 const gmx_mtop_t
*mtop
,
318 const gmx_reverse_top_t
*rt
= cr
->dd
->reverse_top
;
320 /* Print the atoms in the missing interactions per molblock */
322 for (const gmx_molblock_t
&molb
: mtop
->molblock
)
324 const gmx_moltype_t
&moltype
= mtop
->moltype
[molb
.type
];
326 a_end
= a_start
+ molb
.nmol
*moltype
.atoms
.nr
;
328 print_missing_interactions_mb(fplog
, cr
, rt
,
330 &rt
->ril_mt
[molb
.type
],
331 a_start
, a_end
, moltype
.atoms
.nr
,
337 void dd_print_missing_interactions(FILE *fplog
, t_commrec
*cr
,
339 const gmx_mtop_t
*top_global
,
340 const gmx_localtop_t
*top_local
,
341 const t_state
*state_local
)
343 int ndiff_tot
, cl
[F_NRE
], n
, ndiff
, rest_global
, rest_local
;
352 fprintf(fplog
, "\nNot all bonded interactions have been properly assigned to the domain decomposition cells\n");
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
);
370 fprintf(fplog
, "\nA list of missing interactions:\n");
372 fprintf(stderr
, "\nA list of missing interactions:\n");
373 rest_global
= dd
->nbonded_global
;
374 rest_local
= local_count
;
375 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
377 /* In the reverse and local top all constraints are merged
378 * into F_CONSTR. So in the if statement we skip F_CONSTRNC
379 * and add these constraints when doing F_CONSTR.
381 if (((interaction_function
[ftype
].flags
& IF_BOND
) &&
382 (dd
->reverse_top
->bBCheck
383 || !(interaction_function
[ftype
].flags
& IF_LIMZERO
)))
384 || (dd
->reverse_top
->bConstr
&& ftype
== F_CONSTR
)
385 || (dd
->reverse_top
->bSettle
&& ftype
== F_SETTLE
))
387 n
= gmx_mtop_ftype_count(top_global
, ftype
);
388 if (ftype
== F_CONSTR
)
390 n
+= gmx_mtop_ftype_count(top_global
, F_CONSTRNC
);
392 ndiff
= cl
[ftype
] - n
;
395 sprintf(buf
, "%20s of %6d missing %6d",
396 interaction_function
[ftype
].longname
, n
, -ndiff
);
399 fprintf(fplog
, "%s\n", buf
);
401 fprintf(stderr
, "%s\n", buf
);
404 rest_local
-= cl
[ftype
];
408 ndiff
= rest_local
- rest_global
;
411 sprintf(buf
, "%20s of %6d missing %6d", "exclusions",
412 rest_global
, -ndiff
);
415 fprintf(fplog
, "%s\n", buf
);
417 fprintf(stderr
, "%s\n", buf
);
421 print_missing_interactions_atoms(fplog
, cr
, top_global
, &top_local
->idef
);
422 write_dd_pdb("dd_dump_err", 0, "dump", top_global
, cr
,
423 -1, as_rvec_array(state_local
->x
.data()), state_local
->box
);
425 std::string errorMessage
;
429 errorMessage
= "One or more interactions were assigned to multiple domains of the domain decompostion. Please report this bug.";
433 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
));
435 gmx_fatal_collective(FARGS
, cr
->mpi_comm_mygroup
, MASTER(cr
), errorMessage
.c_str());
438 /*! \brief Return global topology molecule information for global atom index \p i_gl */
439 static void global_atomnr_to_moltype_ind(gmx_reverse_top_t
*rt
, int i_gl
,
440 int *mb
, int *mt
, int *mol
, int *i_mol
)
442 molblock_ind_t
*mbi
= rt
->mbi
;
444 int end
= rt
->nmolblock
; /* exclusive */
447 /* binary search for molblock_ind */
451 if (i_gl
>= mbi
[mid
].a_end
)
455 else if (i_gl
< mbi
[mid
].a_start
)
469 *mol
= (i_gl
- mbi
->a_start
) / mbi
->natoms_mol
;
470 *i_mol
= (i_gl
- mbi
->a_start
) - (*mol
)*mbi
->natoms_mol
;
473 /*! \brief Count the exclusions for all atoms in \p cgs */
474 static void count_excls(const t_block
*cgs
, const t_blocka
*excls
,
475 int *n_excl
, int *n_intercg_excl
, int *n_excl_at_max
)
477 int cg
, at0
, at1
, at
, excl
, atj
;
482 for (cg
= 0; cg
< cgs
->nr
; cg
++)
484 at0
= cgs
->index
[cg
];
485 at1
= cgs
->index
[cg
+1];
486 for (at
= at0
; at
< at1
; at
++)
488 for (excl
= excls
->index
[at
]; excl
< excls
->index
[at
+1]; excl
++)
490 atj
= excls
->a
[excl
];
494 if (atj
< at0
|| atj
>= at1
)
501 *n_excl_at_max
= std::max(*n_excl_at_max
,
502 excls
->index
[at
+1] - excls
->index
[at
]);
507 /*! \brief Run the reverse ilist generation and store it in r_il when \p bAssign = TRUE */
508 static int low_make_reverse_ilist(const t_ilist
*il_mt
, const t_atom
*atom
,
509 const int * const * vsite_pbc
,
511 gmx_bool bConstr
, gmx_bool bSettle
,
513 int *r_index
, int *r_il
,
514 gmx_bool bLinkToAllAtoms
,
517 int ftype
, nral
, i
, j
, nlink
, link
;
525 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
527 if ((interaction_function
[ftype
].flags
& (IF_BOND
| IF_VSITE
)) ||
528 (bConstr
&& (ftype
== F_CONSTR
|| ftype
== F_CONSTRNC
)) ||
529 (bSettle
&& ftype
== F_SETTLE
))
531 bVSite
= (interaction_function
[ftype
].flags
& IF_VSITE
);
534 for (i
= 0; i
< il
->nr
; i
+= 1+nral
)
541 /* We don't need the virtual sites for the cg-links */
551 /* Couple to the first atom in the interaction */
554 for (link
= 0; link
< nlink
; link
++)
559 assert(r_il
); //with bAssign not allowed to be null
561 r_il
[r_index
[a
]+count
[a
]] =
562 (ftype
== F_CONSTRNC
? F_CONSTR
: ftype
);
563 r_il
[r_index
[a
]+count
[a
]+1] = ia
[0];
564 for (j
= 1; j
< 1+nral
; j
++)
566 /* Store the molecular atom number */
567 r_il
[r_index
[a
]+count
[a
]+1+j
] = ia
[j
];
570 if (interaction_function
[ftype
].flags
& IF_VSITE
)
574 /* Add an entry to iatoms for storing
575 * which of the constructing atoms are
578 r_il
[r_index
[a
]+count
[a
]+2+nral
] = 0;
579 for (j
= 2; j
< 1+nral
; j
++)
581 if (atom
[ia
[j
]].ptype
== eptVSite
)
583 r_il
[r_index
[a
]+count
[a
]+2+nral
] |= (2<<j
);
586 /* Store vsite pbc atom in a second extra entry */
587 r_il
[r_index
[a
]+count
[a
]+2+nral
+1] =
588 (vsite_pbc
? vsite_pbc
[ftype
-F_VSITE2
][i
/(1+nral
)] : -2);
593 /* We do not count vsites since they are always
594 * uniquely assigned and can be assigned
595 * to multiple nodes with recursive vsites.
598 !(interaction_function
[ftype
].flags
& IF_LIMZERO
))
603 count
[a
] += 2 + nral_rt(ftype
);
612 /*! \brief Make the reverse ilist: a list of bonded interactions linked to atoms */
613 static int make_reverse_ilist(const t_ilist
*ilist
,
614 const t_atoms
*atoms
,
615 const int * const * vsite_pbc
,
616 gmx_bool bConstr
, gmx_bool bSettle
,
618 gmx_bool bLinkToAllAtoms
,
619 reverse_ilist_t
*ril_mt
)
621 int nat_mt
, *count
, i
, nint_mt
;
623 /* Count the interactions */
626 low_make_reverse_ilist(ilist
, atoms
->atom
, vsite_pbc
,
628 bConstr
, bSettle
, bBCheck
, nullptr, nullptr,
629 bLinkToAllAtoms
, FALSE
);
631 snew(ril_mt
->index
, nat_mt
+1);
632 ril_mt
->index
[0] = 0;
633 for (i
= 0; i
< nat_mt
; i
++)
635 ril_mt
->index
[i
+1] = ril_mt
->index
[i
] + count
[i
];
638 snew(ril_mt
->il
, ril_mt
->index
[nat_mt
]);
640 /* Store the interactions */
642 low_make_reverse_ilist(ilist
, atoms
->atom
, vsite_pbc
,
644 bConstr
, bSettle
, bBCheck
,
645 ril_mt
->index
, ril_mt
->il
,
646 bLinkToAllAtoms
, TRUE
);
650 ril_mt
->numAtomsInMolecule
= atoms
->nr
;
655 /*! \brief Destroys a reverse_ilist_t struct */
656 static void destroy_reverse_ilist(reverse_ilist_t
*ril
)
662 /*! \brief Generate the reverse topology */
663 static gmx_reverse_top_t
*make_reverse_top(const gmx_mtop_t
*mtop
, gmx_bool bFE
,
664 const int * const * const * vsite_pbc_molt
,
665 gmx_bool bConstr
, gmx_bool bSettle
,
666 gmx_bool bBCheck
, int *nint
)
668 gmx_reverse_top_t
*rt
;
674 /* Should we include constraints (for SHAKE) in rt? */
675 rt
->bConstr
= bConstr
;
676 rt
->bSettle
= bSettle
;
677 rt
->bBCheck
= bBCheck
;
679 rt
->bInterCGInteractions
= mtop
->bIntermolecularInteractions
;
680 snew(nint_mt
, mtop
->moltype
.size());
681 snew(rt
->ril_mt
, mtop
->moltype
.size());
682 rt
->ril_mt_tot_size
= 0;
683 for (size_t mt
= 0; mt
< mtop
->moltype
.size(); mt
++)
685 const gmx_moltype_t
&molt
= mtop
->moltype
[mt
];
688 rt
->bInterCGInteractions
= TRUE
;
691 /* Make the atom to interaction list for this molecule type */
693 make_reverse_ilist(molt
.ilist
, &molt
.atoms
,
694 vsite_pbc_molt
? vsite_pbc_molt
[mt
] : nullptr,
695 rt
->bConstr
, rt
->bSettle
, rt
->bBCheck
, FALSE
,
698 rt
->ril_mt_tot_size
+= rt
->ril_mt
[mt
].index
[molt
.atoms
.nr
];
702 fprintf(debug
, "The total size of the atom to interaction index is %d integers\n", rt
->ril_mt_tot_size
);
706 for (const gmx_molblock_t
&molblock
: mtop
->molblock
)
708 *nint
+= molblock
.nmol
*nint_mt
[molblock
.type
];
712 /* Make an intermolecular reverse top, if necessary */
713 rt
->bIntermolecularInteractions
= mtop
->bIntermolecularInteractions
;
714 if (rt
->bIntermolecularInteractions
)
716 t_atoms atoms_global
;
718 rt
->ril_intermol
.index
= nullptr;
719 rt
->ril_intermol
.il
= nullptr;
721 atoms_global
.nr
= mtop
->natoms
;
722 atoms_global
.atom
= nullptr; /* Only used with virtual sites */
725 make_reverse_ilist(mtop
->intermolecular_ilist
, &atoms_global
,
727 rt
->bConstr
, rt
->bSettle
, rt
->bBCheck
, FALSE
,
731 if (bFE
&& gmx_mtop_bondeds_free_energy(mtop
))
733 rt
->ilsort
= ilsortFE_UNSORTED
;
737 rt
->ilsort
= ilsortNO_FE
;
740 /* Make a molblock index for fast searching */
741 snew(rt
->mbi
, mtop
->molblock
.size());
742 rt
->nmolblock
= mtop
->molblock
.size();
744 for (size_t mb
= 0; mb
< mtop
->molblock
.size(); mb
++)
746 const gmx_molblock_t
&molb
= mtop
->molblock
[mb
];
747 int numAtomsPerMol
= mtop
->moltype
[molb
.type
].atoms
.nr
;
748 rt
->mbi
[mb
].a_start
= i
;
749 i
+= molb
.nmol
*numAtomsPerMol
;
750 rt
->mbi
[mb
].a_end
= i
;
751 rt
->mbi
[mb
].natoms_mol
= numAtomsPerMol
;
752 rt
->mbi
[mb
].type
= molb
.type
;
755 rt
->nthread
= gmx_omp_nthreads_get(emntDomdec
);
756 snew(rt
->th_work
, rt
->nthread
);
757 if (vsite_pbc_molt
!= nullptr)
759 for (thread
= 0; thread
< rt
->nthread
; thread
++)
761 snew(rt
->th_work
[thread
].vsite_pbc
, F_VSITEN
-F_VSITE2
+1);
762 snew(rt
->th_work
[thread
].vsite_pbc_nalloc
, F_VSITEN
-F_VSITE2
+1);
769 void dd_make_reverse_top(FILE *fplog
,
770 gmx_domdec_t
*dd
, const gmx_mtop_t
*mtop
,
771 const gmx_vsite_t
*vsite
,
772 const t_inputrec
*ir
, gmx_bool bBCheck
)
776 fprintf(fplog
, "\nLinking all bonded interactions to atoms\n");
779 /* If normal and/or settle constraints act only within charge groups,
780 * we can store them in the reverse top and simply assign them to domains.
781 * Otherwise we need to assign them to multiple domains and set up
782 * the parallel version constraint algorithm(s).
785 dd
->reverse_top
= make_reverse_top(mtop
, ir
->efep
!= efepNO
,
786 vsite
? vsite
->vsite_pbc_molt
: nullptr,
787 !dd
->bInterCGcons
, !dd
->bInterCGsettles
,
788 bBCheck
, &dd
->nbonded_global
);
790 gmx_reverse_top_t
*rt
= dd
->reverse_top
;
792 /* With the Verlet scheme, exclusions are handled in the non-bonded
793 * kernels and only exclusions inside the cut-off lead to exclusion
794 * forces. Since each atom pair is treated at most once in the non-bonded
795 * kernels, it doesn't matter if the exclusions for the same atom pair
796 * appear multiple times in the exclusion list. In contrast, the, old,
797 * group cut-off scheme loops over a list of exclusions, so there each
798 * excluded pair should appear exactly once.
800 rt
->bExclRequired
= (ir
->cutoff_scheme
== ecutsGROUP
&&
801 inputrecExclForces(ir
));
804 dd
->n_intercg_excl
= 0;
805 rt
->n_excl_at_max
= 0;
806 for (const gmx_molblock_t
&molb
: mtop
->molblock
)
808 int n_excl_mol
, n_excl_icg
, n_excl_at_max
;
810 const gmx_moltype_t
&molt
= mtop
->moltype
[molb
.type
];
811 count_excls(&molt
.cgs
, &molt
.excls
,
812 &n_excl_mol
, &n_excl_icg
, &n_excl_at_max
);
813 nexcl
+= molb
.nmol
*n_excl_mol
;
814 dd
->n_intercg_excl
+= molb
.nmol
*n_excl_icg
;
815 rt
->n_excl_at_max
= std::max(rt
->n_excl_at_max
, n_excl_at_max
);
817 if (rt
->bExclRequired
)
819 dd
->nbonded_global
+= nexcl
;
820 if (EEL_FULL(ir
->coulombtype
) && dd
->n_intercg_excl
> 0 && fplog
)
822 fprintf(fplog
, "There are %d inter charge-group exclusions,\n"
823 "will use an extra communication step for exclusion forces for %s\n",
824 dd
->n_intercg_excl
, eel_names
[ir
->coulombtype
]);
828 if (vsite
&& vsite
->n_intercg_vsite
> 0)
832 fprintf(fplog
, "There are %d inter charge-group virtual sites,\n"
833 "will an extra communication step for selected coordinates and forces\n",
834 vsite
->n_intercg_vsite
);
836 init_domdec_vsites(dd
, vsite
->n_intercg_vsite
);
839 if (dd
->bInterCGcons
|| dd
->bInterCGsettles
)
841 init_domdec_constraints(dd
, mtop
);
845 fprintf(fplog
, "\n");
849 /*! \brief Store a vsite interaction at the end of \p il
851 * This routine is very similar to add_ifunc, but vsites interactions
852 * have more work to do than other kinds of interactions, and the
853 * complex way nral (and thus vector contents) depends on ftype
854 * confuses static analysis tools unless we fuse the vsite
855 * atom-indexing organization code with the ifunc-adding code, so that
856 * they can see that nral is the same value. */
858 add_ifunc_for_vsites(t_iatom
*tiatoms
, gmx_ga2la_t
*ga2la
,
859 int nral
, gmx_bool bHomeA
,
860 int a
, int a_gl
, int a_mol
,
861 const t_iatom
*iatoms
,
866 if (il
->nr
+1+nral
> il
->nalloc
)
868 il
->nalloc
= over_alloc_large(il
->nr
+1+nral
);
869 srenew(il
->iatoms
, il
->nalloc
);
871 liatoms
= il
->iatoms
+ il
->nr
;
875 tiatoms
[0] = iatoms
[0];
879 /* We know the local index of the first atom */
884 /* Convert later in make_local_vsites */
885 tiatoms
[1] = -a_gl
- 1;
888 for (int k
= 2; k
< 1+nral
; k
++)
890 int ak_gl
= a_gl
+ iatoms
[k
] - a_mol
;
891 if (!ga2la_get_home(ga2la
, ak_gl
, &tiatoms
[k
]))
893 /* Copy the global index, convert later in make_local_vsites */
894 tiatoms
[k
] = -(ak_gl
+ 1);
896 // Note that ga2la_get_home always sets the third parameter if
899 for (int k
= 0; k
< 1+nral
; k
++)
901 liatoms
[k
] = tiatoms
[k
];
905 /*! \brief Store a bonded interaction at the end of \p il */
906 static inline void add_ifunc(int nral
, t_iatom
*tiatoms
, t_ilist
*il
)
911 if (il
->nr
+1+nral
> il
->nalloc
)
913 il
->nalloc
= over_alloc_large(il
->nr
+1+nral
);
914 srenew(il
->iatoms
, il
->nalloc
);
916 liatoms
= il
->iatoms
+ il
->nr
;
917 for (k
= 0; k
<= nral
; k
++)
919 liatoms
[k
] = tiatoms
[k
];
924 /*! \brief Store a position restraint in idef and iatoms, complex because the parameters are different for each entry */
925 static void add_posres(int mol
, int a_mol
, int numAtomsInMolecule
,
926 const gmx_molblock_t
*molb
,
927 t_iatom
*iatoms
, const t_iparams
*ip_in
,
933 /* This position restraint has not been added yet,
934 * so it's index is the current number of position restraints.
936 n
= idef
->il
[F_POSRES
].nr
/2;
937 if (n
+1 > idef
->iparams_posres_nalloc
)
939 idef
->iparams_posres_nalloc
= over_alloc_dd(n
+1);
940 srenew(idef
->iparams_posres
, idef
->iparams_posres_nalloc
);
942 ip
= &idef
->iparams_posres
[n
];
943 /* Copy the force constants */
944 *ip
= ip_in
[iatoms
[0]];
946 /* Get the position restraint coordinates from the molblock */
947 a_molb
= mol
*numAtomsInMolecule
+ a_mol
;
948 GMX_ASSERT(a_molb
< static_cast<int>(molb
->posres_xA
.size()), "We need a sufficient number of position restraint coordinates");
949 ip
->posres
.pos0A
[XX
] = molb
->posres_xA
[a_molb
][XX
];
950 ip
->posres
.pos0A
[YY
] = molb
->posres_xA
[a_molb
][YY
];
951 ip
->posres
.pos0A
[ZZ
] = molb
->posres_xA
[a_molb
][ZZ
];
952 if (!molb
->posres_xB
.empty())
954 ip
->posres
.pos0B
[XX
] = molb
->posres_xB
[a_molb
][XX
];
955 ip
->posres
.pos0B
[YY
] = molb
->posres_xB
[a_molb
][YY
];
956 ip
->posres
.pos0B
[ZZ
] = molb
->posres_xB
[a_molb
][ZZ
];
960 ip
->posres
.pos0B
[XX
] = ip
->posres
.pos0A
[XX
];
961 ip
->posres
.pos0B
[YY
] = ip
->posres
.pos0A
[YY
];
962 ip
->posres
.pos0B
[ZZ
] = ip
->posres
.pos0A
[ZZ
];
964 /* Set the parameter index for idef->iparams_posre */
968 /*! \brief Store a flat-bottomed position restraint in idef and iatoms, complex because the parameters are different for each entry */
969 static void add_fbposres(int mol
, int a_mol
, int numAtomsInMolecule
,
970 const gmx_molblock_t
*molb
,
971 t_iatom
*iatoms
, const t_iparams
*ip_in
,
977 /* This flat-bottom position restraint has not been added yet,
978 * so it's index is the current number of position restraints.
980 n
= idef
->il
[F_FBPOSRES
].nr
/2;
981 if (n
+1 > idef
->iparams_fbposres_nalloc
)
983 idef
->iparams_fbposres_nalloc
= over_alloc_dd(n
+1);
984 srenew(idef
->iparams_fbposres
, idef
->iparams_fbposres_nalloc
);
986 ip
= &idef
->iparams_fbposres
[n
];
987 /* Copy the force constants */
988 *ip
= ip_in
[iatoms
[0]];
990 /* Get the position restraint coordinats from the molblock */
991 a_molb
= mol
*numAtomsInMolecule
+ a_mol
;
992 GMX_ASSERT(a_molb
< static_cast<int>(molb
->posres_xA
.size()), "We need a sufficient number of position restraint coordinates");
993 /* Take reference positions from A position of normal posres */
994 ip
->fbposres
.pos0
[XX
] = molb
->posres_xA
[a_molb
][XX
];
995 ip
->fbposres
.pos0
[YY
] = molb
->posres_xA
[a_molb
][YY
];
996 ip
->fbposres
.pos0
[ZZ
] = molb
->posres_xA
[a_molb
][ZZ
];
998 /* Note: no B-type for flat-bottom posres */
1000 /* Set the parameter index for idef->iparams_posre */
1004 /*! \brief Store a virtual site interaction, complex because of PBC and recursion */
1005 static void add_vsite(gmx_ga2la_t
*ga2la
, const int *index
, const int *rtil
,
1006 int ftype
, int nral
,
1007 gmx_bool bHomeA
, int a
, int a_gl
, int a_mol
,
1008 const t_iatom
*iatoms
,
1009 t_idef
*idef
, int **vsite_pbc
, int *vsite_pbc_nalloc
)
1011 int k
, vsi
, pbc_a_mol
;
1012 t_iatom tiatoms
[1+MAXATOMLIST
];
1013 int j
, ftype_r
, nral_r
;
1015 /* Add this interaction to the local topology */
1016 add_ifunc_for_vsites(tiatoms
, ga2la
, nral
, bHomeA
, a
, a_gl
, a_mol
, iatoms
, &idef
->il
[ftype
]);
1020 vsi
= idef
->il
[ftype
].nr
/(1+nral
) - 1;
1021 if (vsi
>= vsite_pbc_nalloc
[ftype
-F_VSITE2
])
1023 vsite_pbc_nalloc
[ftype
-F_VSITE2
] = over_alloc_large(vsi
+1);
1024 srenew(vsite_pbc
[ftype
-F_VSITE2
], vsite_pbc_nalloc
[ftype
-F_VSITE2
]);
1028 pbc_a_mol
= iatoms
[1+nral
+1];
1031 /* The pbc flag is one of the following two options:
1032 * -2: vsite and all constructing atoms are within the same cg, no pbc
1033 * -1: vsite and its first constructing atom are in the same cg, do pbc
1035 vsite_pbc
[ftype
-F_VSITE2
][vsi
] = pbc_a_mol
;
1039 /* Set the pbc atom for this vsite so we can make its pbc
1040 * identical to the rest of the atoms in its charge group.
1041 * Since the order of the atoms does not change within a charge
1042 * group, we do not need the global to local atom index.
1044 vsite_pbc
[ftype
-F_VSITE2
][vsi
] = a
+ pbc_a_mol
- iatoms
[1];
1049 /* This vsite is non-home (required for recursion),
1050 * and therefore there is no charge group to match pbc with.
1051 * But we always turn on full_pbc to assure that higher order
1052 * recursion works correctly.
1054 vsite_pbc
[ftype
-F_VSITE2
][vsi
] = -1;
1060 /* Check for recursion */
1061 for (k
= 2; k
< 1+nral
; k
++)
1063 if ((iatoms
[1+nral
] & (2<<k
)) && (tiatoms
[k
] < 0))
1065 /* This construction atoms is a vsite and not a home atom */
1068 fprintf(debug
, "Constructing atom %d of vsite atom %d is a vsite and non-home\n", iatoms
[k
]+1, a_mol
+1);
1070 /* Find the vsite construction */
1072 /* Check all interactions assigned to this atom */
1073 j
= index
[iatoms
[k
]];
1074 while (j
< index
[iatoms
[k
]+1])
1076 ftype_r
= rtil
[j
++];
1077 nral_r
= NRAL(ftype_r
);
1078 if (interaction_function
[ftype_r
].flags
& IF_VSITE
)
1080 /* Add this vsite (recursion) */
1081 add_vsite(ga2la
, index
, rtil
, ftype_r
, nral_r
,
1082 FALSE
, -1, a_gl
+iatoms
[k
]-iatoms
[1], iatoms
[k
],
1083 rtil
+j
, idef
, vsite_pbc
, vsite_pbc_nalloc
);
1084 j
+= 1 + nral_r
+ 2;
1096 /*! \brief Build the index that maps each local atom to its local atom group */
1097 static void makeLocalAtomGroupFromAtom(gmx_domdec_t
*dd
)
1099 const gmx::BlockRanges
&atomGroups
= dd
->atomGroups();
1101 dd
->localAtomGroupFromAtom
.clear();
1103 for (size_t g
= 0; g
< dd
->globalAtomGroupIndices
.size(); g
++)
1105 for (int a
= atomGroups
.index
[g
]; a
< atomGroups
.index
[g
+ 1]; a
++)
1107 dd
->localAtomGroupFromAtom
.push_back(g
);
1112 /*! \brief Returns the squared distance between charge groups \p i and \p j */
1113 static real
dd_dist2(t_pbc
*pbc_null
, rvec
*cg_cm
, const int *la2lc
, int i
, int j
)
1119 pbc_dx_aiuc(pbc_null
, cg_cm
[la2lc
[i
]], cg_cm
[la2lc
[j
]], dx
);
1123 rvec_sub(cg_cm
[la2lc
[i
]], cg_cm
[la2lc
[j
]], dx
);
1129 /*! \brief Append t_blocka block structures 1 to nsrc in src to *dest */
1130 static void combine_blocka(t_blocka
*dest
, const thread_work_t
*src
, int nsrc
)
1134 ni
= src
[nsrc
-1].excl
.nr
;
1136 for (s
= 0; s
< nsrc
; s
++)
1138 na
+= src
[s
].excl
.nra
;
1140 if (ni
+ 1 > dest
->nalloc_index
)
1142 dest
->nalloc_index
= over_alloc_large(ni
+1);
1143 srenew(dest
->index
, dest
->nalloc_index
);
1145 if (dest
->nra
+ na
> dest
->nalloc_a
)
1147 dest
->nalloc_a
= over_alloc_large(dest
->nra
+na
);
1148 srenew(dest
->a
, dest
->nalloc_a
);
1150 for (s
= 1; s
< nsrc
; s
++)
1152 for (i
= dest
->nr
+1; i
< src
[s
].excl
.nr
+1; i
++)
1154 dest
->index
[i
] = dest
->nra
+ src
[s
].excl
.index
[i
];
1156 for (i
= 0; i
< src
[s
].excl
.nra
; i
++)
1158 dest
->a
[dest
->nra
+i
] = src
[s
].excl
.a
[i
];
1160 dest
->nr
= src
[s
].excl
.nr
;
1161 dest
->nra
+= src
[s
].excl
.nra
;
1165 /*! \brief Append t_idef structures 1 to nsrc in src to *dest,
1166 * virtual sites need special attention, as pbc info differs per vsite.
1168 static void combine_idef(t_idef
*dest
, const thread_work_t
*src
, int nsrc
,
1173 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
1178 for (s
= 1; s
< nsrc
; s
++)
1180 n
+= src
[s
].idef
.il
[ftype
].nr
;
1186 ild
= &dest
->il
[ftype
];
1188 if (ild
->nr
+ n
> ild
->nalloc
)
1190 ild
->nalloc
= over_alloc_large(ild
->nr
+n
);
1191 srenew(ild
->iatoms
, ild
->nalloc
);
1195 int nral1
= 0, ftv
= 0;
1197 vpbc
= ((interaction_function
[ftype
].flags
& IF_VSITE
) &&
1198 vsite
->vsite_pbc_loc
!= nullptr);
1201 nral1
= 1 + NRAL(ftype
);
1202 ftv
= ftype
- F_VSITE2
;
1203 if ((ild
->nr
+ n
)/nral1
> vsite
->vsite_pbc_loc_nalloc
[ftv
])
1205 vsite
->vsite_pbc_loc_nalloc
[ftv
] =
1206 over_alloc_large((ild
->nr
+ n
)/nral1
);
1207 srenew(vsite
->vsite_pbc_loc
[ftv
],
1208 vsite
->vsite_pbc_loc_nalloc
[ftv
]);
1212 for (s
= 1; s
< nsrc
; s
++)
1217 ils
= &src
[s
].idef
.il
[ftype
];
1218 for (i
= 0; i
< ils
->nr
; i
++)
1220 ild
->iatoms
[ild
->nr
+i
] = ils
->iatoms
[i
];
1224 for (i
= 0; i
< ils
->nr
; i
+= nral1
)
1226 vsite
->vsite_pbc_loc
[ftv
][(ild
->nr
+i
)/nral1
] =
1227 src
[s
].vsite_pbc
[ftv
][i
/nral1
];
1234 /* Position restraints need an additional treatment */
1235 if (ftype
== F_POSRES
|| ftype
== F_FBPOSRES
)
1237 int nposres
= dest
->il
[ftype
].nr
/2;
1238 // TODO: Simplify this code using std::vector
1239 t_iparams
* &iparams_dest
= (ftype
== F_POSRES
? dest
->iparams_posres
: dest
->iparams_fbposres
);
1240 int &posres_nalloc
= (ftype
== F_POSRES
? dest
->iparams_posres_nalloc
: dest
->iparams_fbposres_nalloc
);
1241 if (nposres
> posres_nalloc
)
1243 posres_nalloc
= over_alloc_large(nposres
);
1244 srenew(iparams_dest
, posres_nalloc
);
1247 /* Set nposres to the number of original position restraints in dest */
1248 for (int s
= 1; s
< nsrc
; s
++)
1250 nposres
-= src
[s
].idef
.il
[ftype
].nr
/2;
1253 for (int s
= 1; s
< nsrc
; s
++)
1255 const t_iparams
*iparams_src
= (ftype
== F_POSRES
? src
[s
].idef
.iparams_posres
: src
[s
].idef
.iparams_fbposres
);
1257 for (int i
= 0; i
< src
[s
].idef
.il
[ftype
].nr
/2; i
++)
1259 /* Correct the index into iparams_posres */
1260 dest
->il
[ftype
].iatoms
[nposres
*2] = nposres
;
1261 /* Copy the position restraint force parameters */
1262 iparams_dest
[nposres
] = iparams_src
[i
];
1271 /*! \brief Check and when available assign bonded interactions for local atom i
1274 check_assign_interactions_atom(int i
, int i_gl
,
1276 int numAtomsInMolecule
,
1277 const int *index
, const int *rtil
,
1278 gmx_bool bInterMolInteractions
,
1279 int ind_start
, int ind_end
,
1280 const gmx_domdec_t
*dd
,
1281 const gmx_domdec_zones_t
*zones
,
1282 const gmx_molblock_t
*molb
,
1283 gmx_bool bRCheckMB
, ivec rcheck
, gmx_bool bRCheck2B
,
1288 const t_iparams
*ip_in
,
1290 int **vsite_pbc
, int *vsite_pbc_nalloc
,
1301 const t_iatom
*iatoms
;
1303 t_iatom tiatoms
[1 + MAXATOMLIST
];
1308 if (ftype
== F_SETTLE
)
1310 /* Settles are only in the reverse top when they
1311 * operate within a charge group. So we can assign
1312 * them without checks. We do this only for performance
1313 * reasons; it could be handled by the code below.
1317 /* Home zone: add this settle to the local topology */
1318 tiatoms
[0] = iatoms
[0];
1320 tiatoms
[2] = i
+ iatoms
[2] - iatoms
[1];
1321 tiatoms
[3] = i
+ iatoms
[3] - iatoms
[1];
1322 add_ifunc(nral
, tiatoms
, &idef
->il
[ftype
]);
1327 else if (interaction_function
[ftype
].flags
& IF_VSITE
)
1329 assert(!bInterMolInteractions
);
1330 /* The vsite construction goes where the vsite itself is */
1333 add_vsite(dd
->ga2la
, index
, rtil
, ftype
, nral
,
1334 TRUE
, i
, i_gl
, i_mol
,
1335 iatoms
, idef
, vsite_pbc
, vsite_pbc_nalloc
);
1344 tiatoms
[0] = iatoms
[0];
1348 assert(!bInterMolInteractions
);
1349 /* Assign single-body interactions to the home zone */
1354 if (ftype
== F_POSRES
)
1356 add_posres(mol
, i_mol
, numAtomsInMolecule
,
1357 molb
, tiatoms
, ip_in
, idef
);
1359 else if (ftype
== F_FBPOSRES
)
1361 add_fbposres(mol
, i_mol
, numAtomsInMolecule
,
1362 molb
, tiatoms
, ip_in
, idef
);
1372 /* This is a two-body interaction, we can assign
1373 * analogous to the non-bonded assignments.
1375 int k_gl
, a_loc
, kz
;
1377 if (!bInterMolInteractions
)
1379 /* Get the global index using the offset in the molecule */
1380 k_gl
= i_gl
+ iatoms
[2] - i_mol
;
1386 if (!ga2la_get(dd
->ga2la
, k_gl
, &a_loc
, &kz
))
1396 /* Check zone interaction assignments */
1397 bUse
= ((iz
< zones
->nizone
&&
1399 kz
>= zones
->izone
[iz
].j0
&&
1400 kz
< zones
->izone
[iz
].j1
) ||
1401 (kz
< zones
->nizone
&&
1403 iz
>= zones
->izone
[kz
].j0
&&
1404 iz
< zones
->izone
[kz
].j1
));
1409 /* If necessary check the cgcm distance */
1411 dd_dist2(pbc_null
, cg_cm
, la2lc
,
1412 tiatoms
[1], tiatoms
[2]) >= rc2
)
1421 /* Assign this multi-body bonded interaction to
1422 * the local node if we have all the atoms involved
1423 * (local or communicated) and the minimum zone shift
1424 * in each dimension is zero, for dimensions
1425 * with 2 DD cells an extra check may be necessary.
1427 ivec k_zero
, k_plus
;
1433 for (k
= 1; k
<= nral
&& bUse
; k
++)
1439 if (!bInterMolInteractions
)
1441 /* Get the global index using the offset in the molecule */
1442 k_gl
= i_gl
+ iatoms
[k
] - i_mol
;
1448 bLocal
= ga2la_get(dd
->ga2la
, k_gl
, &a_loc
, &kz
);
1449 if (!bLocal
|| kz
>= zones
->n
)
1451 /* We do not have this atom of this interaction
1452 * locally, or it comes from more than one cell
1462 for (d
= 0; d
< DIM
; d
++)
1464 if (zones
->shift
[kz
][d
] == 0)
1476 k_zero
[XX
] && k_zero
[YY
] && k_zero
[ZZ
]);
1481 for (d
= 0; (d
< DIM
&& bUse
); d
++)
1483 /* Check if the cg_cm distance falls within
1484 * the cut-off to avoid possible multiple
1485 * assignments of bonded interactions.
1489 dd_dist2(pbc_null
, cg_cm
, la2lc
,
1490 tiatoms
[k_zero
[d
]], tiatoms
[k_plus
[d
]]) >= rc2
)
1499 /* Add this interaction to the local topology */
1500 add_ifunc(nral
, tiatoms
, &idef
->il
[ftype
]);
1501 /* Sum so we can check in global_stat
1502 * if we have everything.
1505 !(interaction_function
[ftype
].flags
& IF_LIMZERO
))
1515 /*! \brief This function looks up and assigns bonded interactions for zone iz.
1517 * With thread parallelizing each thread acts on a different atom range:
1518 * at_start to at_end.
1520 static int make_bondeds_zone(gmx_domdec_t
*dd
,
1521 const gmx_domdec_zones_t
*zones
,
1522 const std::vector
<gmx_molblock_t
> &molb
,
1523 gmx_bool bRCheckMB
, ivec rcheck
, gmx_bool bRCheck2B
,
1525 int *la2lc
, t_pbc
*pbc_null
, rvec
*cg_cm
,
1526 const t_iparams
*ip_in
,
1529 int *vsite_pbc_nalloc
,
1531 int at_start
, int at_end
)
1533 int mb
, mt
, mol
, i_mol
;
1536 gmx_reverse_top_t
*rt
;
1539 rt
= dd
->reverse_top
;
1541 bBCheck
= rt
->bBCheck
;
1545 for (int i
= at_start
; i
< at_end
; i
++)
1547 /* Get the global atom number */
1548 const int i_gl
= dd
->globalAtomIndices
[i
];
1549 global_atomnr_to_moltype_ind(rt
, i_gl
, &mb
, &mt
, &mol
, &i_mol
);
1550 /* Check all intramolecular interactions assigned to this atom */
1551 index
= rt
->ril_mt
[mt
].index
;
1552 rtil
= rt
->ril_mt
[mt
].il
;
1554 check_assign_interactions_atom(i
, i_gl
, mol
, i_mol
,
1555 rt
->ril_mt
[mt
].numAtomsInMolecule
,
1557 index
[i_mol
], index
[i_mol
+1],
1560 bRCheckMB
, rcheck
, bRCheck2B
, rc2
,
1565 idef
, vsite_pbc
, vsite_pbc_nalloc
,
1571 if (rt
->bIntermolecularInteractions
)
1573 /* Check all intermolecular interactions assigned to this atom */
1574 index
= rt
->ril_intermol
.index
;
1575 rtil
= rt
->ril_intermol
.il
;
1577 check_assign_interactions_atom(i
, i_gl
, mol
, i_mol
,
1578 rt
->ril_mt
[mt
].numAtomsInMolecule
,
1580 index
[i_gl
], index
[i_gl
+ 1],
1583 bRCheckMB
, rcheck
, bRCheck2B
, rc2
,
1588 idef
, vsite_pbc
, vsite_pbc_nalloc
,
1595 return nbonded_local
;
1598 /*! \brief Set the exclusion data for i-zone \p iz for the case of no exclusions */
1599 static void set_no_exclusions_zone(gmx_domdec_t
*dd
, gmx_domdec_zones_t
*zones
,
1600 int iz
, t_blocka
*lexcls
)
1602 int a0
= dd
->atomGroups().index
[zones
->cg_range
[iz
]];
1603 int a1
= dd
->atomGroups().index
[zones
->cg_range
[iz
+ 1]];
1605 for (int a
= a0
+ 1; a
< a1
+ 1; a
++)
1607 lexcls
->index
[a
] = lexcls
->nra
;
1611 /*! \brief Set the exclusion data for i-zone \p iz
1613 * This is a legacy version for the group scheme of the same routine below.
1614 * Here charge groups and distance checks to ensure unique exclusions
1617 static int make_exclusions_zone_cg(gmx_domdec_t
*dd
, gmx_domdec_zones_t
*zones
,
1618 const std::vector
<gmx_moltype_t
> &moltype
,
1619 gmx_bool bRCheck
, real rc2
,
1620 int *la2lc
, t_pbc
*pbc_null
, rvec
*cg_cm
,
1624 int cg_start
, int cg_end
)
1628 const t_blocka
*excls
;
1634 int jla0
= dd
->atomGroups().index
[zones
->izone
[iz
].jcg0
];
1635 int jla1
= dd
->atomGroups().index
[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
= dd
->atomGroups().index
[cg_end
];
1642 int n
= lexcls
->nra
;
1644 for (int cg
= cg_start
; cg
< cg_end
; cg
++)
1646 if (n
+ (cg_end
- cg_start
)*n_excl_at_max
> lexcls
->nalloc_a
)
1648 lexcls
->nalloc_a
= over_alloc_large(n
+ (cg_end
- cg_start
)*n_excl_at_max
);
1649 srenew(lexcls
->a
, lexcls
->nalloc_a
);
1651 int la0
= dd
->atomGroups().index
[cg
];
1652 int la1
= dd
->atomGroups().index
[cg
+ 1];
1653 if (GET_CGINFO_EXCL_INTER(cginfo
[cg
]) ||
1654 !GET_CGINFO_EXCL_INTRA(cginfo
[cg
]))
1656 /* Copy the exclusions from the global top */
1657 for (int la
= la0
; la
< la1
; la
++)
1659 lexcls
->index
[la
] = n
;
1660 int a_gl
= dd
->globalAtomIndices
[la
];
1662 global_atomnr_to_moltype_ind(dd
->reverse_top
, a_gl
, &mb
, &mt
, &mol
, &a_mol
);
1663 excls
= &moltype
[mt
].excls
;
1664 for (int j
= excls
->index
[a_mol
]; j
< excls
->index
[a_mol
+1]; j
++)
1666 int aj_mol
= excls
->a
[j
];
1667 /* This computation of jla is only correct intra-cg */
1668 int jla
= la
+ aj_mol
- a_mol
;
1669 if (jla
>= la0
&& jla
< la1
)
1671 /* This is an intra-cg exclusion. We can skip
1672 * the global indexing and distance checking.
1674 /* Intra-cg exclusions are only required
1675 * for the home zone.
1679 lexcls
->a
[n
++] = jla
;
1680 /* Check to avoid double counts */
1689 /* This is a inter-cg exclusion */
1690 /* Since exclusions are pair interactions,
1691 * just like non-bonded interactions,
1692 * they can be assigned properly up
1693 * to the DD cutoff (not cutoff_min as
1694 * for the other bonded interactions).
1696 if (ga2la_get(ga2la
, a_gl
+aj_mol
-a_mol
, &jla
, &cell
))
1698 if (iz
== 0 && cell
== 0)
1700 lexcls
->a
[n
++] = jla
;
1701 /* Check to avoid double counts */
1707 else if (jla
>= jla0
&& jla
< jla1
&&
1709 dd_dist2(pbc_null
, cg_cm
, la2lc
, la
, jla
) < rc2
))
1711 /* jla > la, since jla0 > la */
1712 lexcls
->a
[n
++] = jla
;
1722 /* There are no inter-cg excls and this cg is self-excluded.
1723 * These exclusions are only required for zone 0,
1724 * since other zones do not see themselves.
1728 for (int la
= la0
; la
< la1
; la
++)
1730 lexcls
->index
[la
] = n
;
1731 for (int j
= la0
; j
< la1
; j
++)
1736 count
+= ((la1
- la0
)*(la1
- la0
- 1))/2;
1740 /* We don't need exclusions for this cg */
1741 for (int la
= la0
; la
< la1
; la
++)
1743 lexcls
->index
[la
] = n
;
1749 lexcls
->index
[lexcls
->nr
] = n
;
1755 /*! \brief Set the exclusion data for i-zone \p iz */
1756 static void make_exclusions_zone(gmx_domdec_t
*dd
,
1757 gmx_domdec_zones_t
*zones
,
1758 const std::vector
<gmx_moltype_t
> &moltype
,
1762 int at_start
, int at_end
)
1765 int n_excl_at_max
, n
, at
;
1769 int jla0
= dd
->atomGroups().index
[zones
->izone
[iz
].jcg0
];
1770 int jla1
= dd
->atomGroups().index
[zones
->izone
[iz
].jcg1
];
1772 n_excl_at_max
= dd
->reverse_top
->n_excl_at_max
;
1774 /* We set the end index, but note that we might not start at zero here */
1775 lexcls
->nr
= at_end
;
1778 for (at
= at_start
; at
< at_end
; at
++)
1780 if (n
+ 1000 > lexcls
->nalloc_a
)
1782 lexcls
->nalloc_a
= over_alloc_large(n
+ 1000);
1783 srenew(lexcls
->a
, lexcls
->nalloc_a
);
1785 if (GET_CGINFO_EXCL_INTER(cginfo
[at
]))
1787 int a_gl
, mb
, mt
, mol
, a_mol
, j
;
1788 const t_blocka
*excls
;
1790 if (n
+ n_excl_at_max
> lexcls
->nalloc_a
)
1792 lexcls
->nalloc_a
= over_alloc_large(n
+ n_excl_at_max
);
1793 srenew(lexcls
->a
, lexcls
->nalloc_a
);
1796 /* Copy the exclusions from the global top */
1797 lexcls
->index
[at
] = n
;
1798 a_gl
= dd
->globalAtomIndices
[at
];
1799 global_atomnr_to_moltype_ind(dd
->reverse_top
, a_gl
,
1800 &mb
, &mt
, &mol
, &a_mol
);
1801 excls
= &moltype
[mt
].excls
;
1802 for (j
= excls
->index
[a_mol
]; j
< excls
->index
[a_mol
+ 1]; j
++)
1804 int aj_mol
, at_j
, cell
;
1806 aj_mol
= excls
->a
[j
];
1808 if (ga2la_get(ga2la
, a_gl
+ aj_mol
- a_mol
, &at_j
, &cell
))
1810 /* This check is not necessary, but it can reduce
1811 * the number of exclusions in the list, which in turn
1812 * can speed up the pair list construction a bit.
1814 if (at_j
>= jla0
&& at_j
< jla1
)
1816 lexcls
->a
[n
++] = at_j
;
1823 /* We don't need exclusions for this atom */
1824 lexcls
->index
[at
] = n
;
1828 lexcls
->index
[lexcls
->nr
] = n
;
1833 /*! \brief Ensure we have enough space in \p ba for \p nindex_max indices */
1834 static void check_alloc_index(t_blocka
*ba
, int nindex_max
)
1836 if (nindex_max
+1 > ba
->nalloc_index
)
1838 ba
->nalloc_index
= over_alloc_dd(nindex_max
+1);
1839 srenew(ba
->index
, ba
->nalloc_index
);
1843 /*! \brief Ensure that we have enough space for exclusion storate in \p lexcls */
1844 static void check_exclusions_alloc(gmx_domdec_t
*dd
, gmx_domdec_zones_t
*zones
,
1847 int nr
= dd
->atomGroups().index
[zones
->izone
[zones
->nizone
- 1].cg1
];
1849 check_alloc_index(lexcls
, nr
);
1851 for (int thread
= 1; thread
< dd
->reverse_top
->nthread
; thread
++)
1853 check_alloc_index(&dd
->reverse_top
->th_work
[thread
].excl
, nr
);
1857 /*! \brief Set the total count indexes for the local exclusions, needed by several functions */
1858 static void finish_local_exclusions(gmx_domdec_t
*dd
, gmx_domdec_zones_t
*zones
,
1861 lexcls
->nr
= dd
->atomGroups().index
[zones
->izone
[zones
->nizone
- 1].cg1
];
1863 if (dd
->n_intercg_excl
== 0)
1865 /* There are no exclusions involving non-home charge groups,
1866 * but we need to set the indices for neighborsearching.
1868 int la0
= dd
->atomGroups().index
[zones
->izone
[0].cg1
];
1869 for (int la
= la0
; la
< lexcls
->nr
; la
++)
1871 lexcls
->index
[la
] = lexcls
->nra
;
1874 /* nr is only used to loop over the exclusions for Ewald and RF,
1875 * so we can set it to the number of home atoms for efficiency.
1877 lexcls
->nr
= dd
->atomGroups().index
[zones
->izone
[0].cg1
];
1881 /*! \brief Clear a t_idef struct */
1882 static void clear_idef(t_idef
*idef
)
1886 /* Clear the counts */
1887 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
1889 idef
->il
[ftype
].nr
= 0;
1893 /*! \brief Generate and store all required local bonded interactions in \p idef and local exclusions in \p lexcls */
1894 static int make_local_bondeds_excls(gmx_domdec_t
*dd
,
1895 gmx_domdec_zones_t
*zones
,
1896 const gmx_mtop_t
*mtop
,
1898 gmx_bool bRCheckMB
, ivec rcheck
, gmx_bool bRCheck2B
,
1900 int *la2lc
, t_pbc
*pbc_null
, rvec
*cg_cm
,
1901 t_idef
*idef
, gmx_vsite_t
*vsite
,
1902 t_blocka
*lexcls
, int *excl_count
)
1904 int nzone_bondeds
, nzone_excl
;
1905 int izone
, cg0
, cg1
;
1909 gmx_reverse_top_t
*rt
;
1911 if (dd
->reverse_top
->bInterCGInteractions
)
1913 nzone_bondeds
= zones
->n
;
1917 /* Only single charge group (or atom) molecules, so interactions don't
1918 * cross zone boundaries and we only need to assign in the home zone.
1923 if (dd
->n_intercg_excl
> 0)
1925 /* We only use exclusions from i-zones to i- and j-zones */
1926 nzone_excl
= zones
->nizone
;
1930 /* There are no inter-cg exclusions and only zone 0 sees itself */
1934 check_exclusions_alloc(dd
, zones
, lexcls
);
1936 rt
= dd
->reverse_top
;
1940 /* Clear the counts */
1948 for (izone
= 0; izone
< nzone_bondeds
; izone
++)
1950 cg0
= zones
->cg_range
[izone
];
1951 cg1
= zones
->cg_range
[izone
+ 1];
1953 #pragma omp parallel for num_threads(rt->nthread) schedule(static)
1954 for (thread
= 0; thread
< rt
->nthread
; thread
++)
1961 int *vsite_pbc_nalloc
;
1964 cg0t
= cg0
+ ((cg1
- cg0
)* thread
)/rt
->nthread
;
1965 cg1t
= cg0
+ ((cg1
- cg0
)*(thread
+1))/rt
->nthread
;
1973 idef_t
= &rt
->th_work
[thread
].idef
;
1977 if (vsite
&& vsite
->bHaveChargeGroups
&& vsite
->n_intercg_vsite
> 0)
1981 vsite_pbc
= vsite
->vsite_pbc_loc
;
1982 vsite_pbc_nalloc
= vsite
->vsite_pbc_loc_nalloc
;
1986 vsite_pbc
= rt
->th_work
[thread
].vsite_pbc
;
1987 vsite_pbc_nalloc
= rt
->th_work
[thread
].vsite_pbc_nalloc
;
1992 vsite_pbc
= nullptr;
1993 vsite_pbc_nalloc
= nullptr;
1996 rt
->th_work
[thread
].nbonded
=
1997 make_bondeds_zone(dd
, zones
,
1999 bRCheckMB
, rcheck
, bRCheck2B
, rc2
,
2000 la2lc
, pbc_null
, cg_cm
, idef
->iparams
,
2002 vsite_pbc
, vsite_pbc_nalloc
,
2004 dd
->atomGroups().index
[cg0t
],
2005 dd
->atomGroups().index
[cg1t
]);
2007 if (izone
< nzone_excl
)
2015 excl_t
= &rt
->th_work
[thread
].excl
;
2020 int numAtomGroups
= dd
->globalAtomGroupIndices
.size();
2021 if (dd
->atomGroups().index
[numAtomGroups
] == numAtomGroups
&&
2024 /* No charge groups and no distance check required */
2025 make_exclusions_zone(dd
, zones
,
2026 mtop
->moltype
, cginfo
,
2033 rt
->th_work
[thread
].excl_count
=
2034 make_exclusions_zone_cg(dd
, zones
,
2035 mtop
->moltype
, bRCheck2B
, rc2
,
2036 la2lc
, pbc_null
, cg_cm
, cginfo
,
2043 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
2046 if (rt
->nthread
> 1)
2048 combine_idef(idef
, rt
->th_work
, rt
->nthread
, vsite
);
2051 for (thread
= 0; thread
< rt
->nthread
; thread
++)
2053 nbonded_local
+= rt
->th_work
[thread
].nbonded
;
2056 if (izone
< nzone_excl
)
2058 if (rt
->nthread
> 1)
2060 combine_blocka(lexcls
, rt
->th_work
, rt
->nthread
);
2063 for (thread
= 0; thread
< rt
->nthread
; thread
++)
2065 *excl_count
+= rt
->th_work
[thread
].excl_count
;
2070 /* Some zones might not have exclusions, but some code still needs to
2071 * loop over the index, so we set the indices here.
2073 for (izone
= nzone_excl
; izone
< zones
->nizone
; izone
++)
2075 set_no_exclusions_zone(dd
, zones
, izone
, lexcls
);
2078 finish_local_exclusions(dd
, zones
, lexcls
);
2081 fprintf(debug
, "We have %d exclusions, check count %d\n",
2082 lexcls
->nra
, *excl_count
);
2085 return nbonded_local
;
2088 void dd_make_local_cgs(gmx_domdec_t
*dd
, t_block
*lcgs
)
2090 lcgs
->nr
= dd
->globalAtomGroupIndices
.size();
2091 lcgs
->index
= dd
->atomGroups_
.index
.data();
2094 void dd_make_local_top(gmx_domdec_t
*dd
, gmx_domdec_zones_t
*zones
,
2095 int npbcdim
, matrix box
,
2096 rvec cellsize_min
, ivec npulse
,
2100 const gmx_mtop_t
*mtop
, gmx_localtop_t
*ltop
)
2102 gmx_bool bRCheckMB
, bRCheck2B
;
2106 t_pbc pbc
, *pbc_null
= nullptr;
2110 fprintf(debug
, "Making local topology\n");
2113 dd_make_local_cgs(dd
, <op
->cgs
);
2118 if (dd
->reverse_top
->bInterCGInteractions
)
2120 /* We need to check to which cell bondeds should be assigned */
2121 rc
= dd_cutoff_twobody(dd
);
2124 fprintf(debug
, "Two-body bonded cut-off distance is %g\n", rc
);
2127 /* Should we check cg_cm distances when assigning bonded interactions? */
2128 for (d
= 0; d
< DIM
; d
++)
2131 /* Only need to check for dimensions where the part of the box
2132 * that is not communicated is smaller than the cut-off.
2134 if (d
< npbcdim
&& dd
->nc
[d
] > 1 &&
2135 (dd
->nc
[d
] - npulse
[d
])*cellsize_min
[d
] < 2*rc
)
2142 /* Check for interactions between two atoms,
2143 * where we can allow interactions up to the cut-off,
2144 * instead of up to the smallest cell dimension.
2151 "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %d\n",
2152 d
, cellsize_min
[d
], d
, rcheck
[d
], bRCheck2B
);
2155 if (bRCheckMB
|| bRCheck2B
)
2157 makeLocalAtomGroupFromAtom(dd
);
2160 pbc_null
= set_pbc_dd(&pbc
, fr
->ePBC
, dd
->nc
, TRUE
, box
);
2170 make_local_bondeds_excls(dd
, zones
, mtop
, fr
->cginfo
,
2171 bRCheckMB
, rcheck
, bRCheck2B
, rc
,
2172 dd
->localAtomGroupFromAtom
.data(),
2173 pbc_null
, cgcm_or_x
,
2175 <op
->excls
, &nexcl
);
2177 /* The ilist is not sorted yet,
2178 * we can only do this when we have the charge arrays.
2180 ltop
->idef
.ilsort
= ilsortUNKNOWN
;
2182 if (dd
->reverse_top
->bExclRequired
)
2184 dd
->nbonded_local
+= nexcl
;
2187 ltop
->atomtypes
= mtop
->atomtypes
;
2190 void dd_sort_local_top(gmx_domdec_t
*dd
, const t_mdatoms
*mdatoms
,
2191 gmx_localtop_t
*ltop
)
2193 if (dd
->reverse_top
->ilsort
== ilsortNO_FE
)
2195 ltop
->idef
.ilsort
= ilsortNO_FE
;
2199 gmx_sort_ilist_fe(<op
->idef
, mdatoms
->chargeA
, mdatoms
->chargeB
);
2203 gmx_localtop_t
*dd_init_local_top(const gmx_mtop_t
*top_global
)
2205 gmx_localtop_t
*top
;
2210 top
->idef
.ntypes
= top_global
->ffparams
.ntypes
;
2211 top
->idef
.atnr
= top_global
->ffparams
.atnr
;
2212 top
->idef
.functype
= top_global
->ffparams
.functype
;
2213 top
->idef
.iparams
= top_global
->ffparams
.iparams
;
2214 top
->idef
.fudgeQQ
= top_global
->ffparams
.fudgeQQ
;
2215 top
->idef
.cmap_grid
= top_global
->ffparams
.cmap_grid
;
2217 for (i
= 0; i
< F_NRE
; i
++)
2219 top
->idef
.il
[i
].iatoms
= nullptr;
2220 top
->idef
.il
[i
].nalloc
= 0;
2222 top
->idef
.ilsort
= ilsortUNKNOWN
;
2227 void dd_init_local_state(gmx_domdec_t
*dd
,
2228 t_state
*state_global
, t_state
*state_local
)
2230 int buf
[NITEM_DD_INIT_LOCAL_STATE
];
2234 buf
[0] = state_global
->flags
;
2235 buf
[1] = state_global
->ngtc
;
2236 buf
[2] = state_global
->nnhpres
;
2237 buf
[3] = state_global
->nhchainlength
;
2238 buf
[4] = state_global
->dfhist
? state_global
->dfhist
->nlambda
: 0;
2240 dd_bcast(dd
, NITEM_DD_INIT_LOCAL_STATE
*sizeof(int), buf
);
2242 init_gtc_state(state_local
, buf
[1], buf
[2], buf
[3]);
2243 init_dfhist_state(state_local
, buf
[4]);
2244 state_local
->flags
= buf
[0];
2247 /*! \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 */
2248 static void check_link(t_blocka
*link
, int cg_gl
, int cg_gl_j
)
2254 for (k
= link
->index
[cg_gl
]; k
< link
->index
[cg_gl
+1]; k
++)
2256 GMX_RELEASE_ASSERT(link
->a
, "Inconsistent NULL pointer while making charge-group links");
2257 if (link
->a
[k
] == cg_gl_j
)
2264 GMX_RELEASE_ASSERT(link
->a
|| link
->index
[cg_gl
+1]+1 > link
->nalloc_a
,
2265 "Inconsistent allocation of link");
2266 /* Add this charge group link */
2267 if (link
->index
[cg_gl
+1]+1 > link
->nalloc_a
)
2269 link
->nalloc_a
= over_alloc_large(link
->index
[cg_gl
+1]+1);
2270 srenew(link
->a
, link
->nalloc_a
);
2272 link
->a
[link
->index
[cg_gl
+1]] = cg_gl_j
;
2273 link
->index
[cg_gl
+1]++;
2277 /*! \brief Return a vector of the charge group index for all atoms */
2278 static std::vector
<int> make_at2cg(const t_block
&cgs
)
2280 std::vector
<int> at2cg(cgs
.index
[cgs
.nr
]);
2281 for (int cg
= 0; cg
< cgs
.nr
; cg
++)
2283 for (int a
= cgs
.index
[cg
]; a
< cgs
.index
[cg
+ 1]; a
++)
2292 t_blocka
*make_charge_group_links(const gmx_mtop_t
*mtop
, gmx_domdec_t
*dd
,
2293 cginfo_mb_t
*cginfo_mb
)
2295 gmx_bool bExclRequired
;
2296 reverse_ilist_t ril
, ril_intermol
;
2298 cginfo_mb_t
*cgi_mb
;
2300 /* For each charge group make a list of other charge groups
2301 * in the system that a linked to it via bonded interactions
2302 * which are also stored in reverse_top.
2305 bExclRequired
= dd
->reverse_top
->bExclRequired
;
2307 if (mtop
->bIntermolecularInteractions
)
2309 if (ncg_mtop(mtop
) < mtop
->natoms
)
2311 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.");
2316 atoms
.nr
= mtop
->natoms
;
2317 atoms
.atom
= nullptr;
2319 make_reverse_ilist(mtop
->intermolecular_ilist
, &atoms
,
2320 nullptr, FALSE
, FALSE
, FALSE
, TRUE
, &ril_intermol
);
2324 snew(link
->index
, ncg_mtop(mtop
)+1);
2331 for (size_t mb
= 0; mb
< mtop
->molblock
.size(); mb
++)
2333 const gmx_molblock_t
&molb
= mtop
->molblock
[mb
];
2338 const gmx_moltype_t
&molt
= mtop
->moltype
[molb
.type
];
2339 const t_block
&cgs
= molt
.cgs
;
2340 const t_blocka
&excls
= molt
.excls
;
2341 std::vector
<int> a2c
= make_at2cg(cgs
);
2342 /* Make a reverse ilist in which the interactions are linked
2343 * to all atoms, not only the first atom as in gmx_reverse_top.
2344 * The constraints are discarded here.
2346 make_reverse_ilist(molt
.ilist
, &molt
.atoms
,
2347 nullptr, FALSE
, FALSE
, FALSE
, TRUE
, &ril
);
2349 cgi_mb
= &cginfo_mb
[mb
];
2352 for (mol
= 0; mol
< (mtop
->bIntermolecularInteractions
? molb
.nmol
: 1); mol
++)
2354 for (int cg
= 0; cg
< cgs
.nr
; cg
++)
2356 int cg_gl
= cg_offset
+ cg
;
2357 link
->index
[cg_gl
+1] = link
->index
[cg_gl
];
2358 for (int a
= cgs
.index
[cg
]; a
< cgs
.index
[cg
+ 1]; a
++)
2360 int i
= ril
.index
[a
];
2361 while (i
< ril
.index
[a
+1])
2363 int ftype
= ril
.il
[i
++];
2364 int nral
= NRAL(ftype
);
2365 /* Skip the ifunc index */
2367 for (int j
= 0; j
< nral
; j
++)
2369 int aj
= ril
.il
[i
+ j
];
2372 check_link(link
, cg_gl
, cg_offset
+a2c
[aj
]);
2375 i
+= nral_rt(ftype
);
2379 /* Exclusions always go both ways */
2380 for (int j
= excls
.index
[a
]; j
< excls
.index
[a
+ 1]; j
++)
2382 int aj
= excls
.a
[j
];
2385 check_link(link
, cg_gl
, cg_offset
+a2c
[aj
]);
2390 if (mtop
->bIntermolecularInteractions
)
2392 int i
= ril_intermol
.index
[a
];
2393 while (i
< ril_intermol
.index
[a
+1])
2395 int ftype
= ril_intermol
.il
[i
++];
2396 int nral
= NRAL(ftype
);
2397 /* Skip the ifunc index */
2399 for (int j
= 0; j
< nral
; j
++)
2401 /* Here we assume we have no charge groups;
2402 * this has been checked above.
2404 int aj
= ril_intermol
.il
[i
+ j
];
2405 check_link(link
, cg_gl
, aj
);
2407 i
+= nral_rt(ftype
);
2411 if (link
->index
[cg_gl
+1] - link
->index
[cg_gl
] > 0)
2413 SET_CGINFO_BOND_INTER(cgi_mb
->cginfo
[cg
]);
2418 cg_offset
+= cgs
.nr
;
2420 int nlink_mol
= link
->index
[cg_offset
] - link
->index
[cg_offset
- cgs
.nr
];
2422 destroy_reverse_ilist(&ril
);
2426 fprintf(debug
, "molecule type '%s' %d cgs has %d cg links through bonded interac.\n", *molt
.name
, cgs
.nr
, nlink_mol
);
2429 if (molb
.nmol
> mol
)
2431 /* Copy the data for the rest of the molecules in this block */
2432 link
->nalloc_a
+= (molb
.nmol
- mol
)*nlink_mol
;
2433 srenew(link
->a
, link
->nalloc_a
);
2434 for (; mol
< molb
.nmol
; mol
++)
2436 for (int cg
= 0; cg
< cgs
.nr
; cg
++)
2438 int cg_gl
= cg_offset
+ cg
;
2439 link
->index
[cg_gl
+ 1] =
2440 link
->index
[cg_gl
+ 1 - cgs
.nr
] + nlink_mol
;
2441 for (int j
= link
->index
[cg_gl
]; j
< link
->index
[cg_gl
+1]; j
++)
2443 link
->a
[j
] = link
->a
[j
- nlink_mol
] + cgs
.nr
;
2445 if (link
->index
[cg_gl
+1] - link
->index
[cg_gl
] > 0 &&
2446 cg_gl
- cgi_mb
->cg_start
< cgi_mb
->cg_mod
)
2448 SET_CGINFO_BOND_INTER(cgi_mb
->cginfo
[cg_gl
- cgi_mb
->cg_start
]);
2452 cg_offset
+= cgs
.nr
;
2457 if (mtop
->bIntermolecularInteractions
)
2459 destroy_reverse_ilist(&ril_intermol
);
2464 fprintf(debug
, "Of the %d charge groups %d are linked via bonded interactions\n", ncg_mtop(mtop
), ncgi
);
2475 } bonded_distance_t
;
2477 /*! \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 */
2478 static void update_max_bonded_distance(real r2
, int ftype
, int a1
, int a2
,
2479 bonded_distance_t
*bd
)
2490 /*! \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 */
2491 static void bonded_cg_distance_mol(const gmx_moltype_t
*molt
,
2492 const std::vector
<int> &at2cg
,
2493 gmx_bool bBCheck
, gmx_bool bExcl
, rvec
*cg_cm
,
2494 bonded_distance_t
*bd_2b
,
2495 bonded_distance_t
*bd_mb
)
2497 for (int ftype
= 0; ftype
< F_NRE
; ftype
++)
2499 if (dd_check_ftype(ftype
, bBCheck
, FALSE
, FALSE
))
2501 const t_ilist
*il
= &molt
->ilist
[ftype
];
2502 int nral
= NRAL(ftype
);
2505 for (int i
= 0; i
< il
->nr
; i
+= 1+nral
)
2507 for (int ai
= 0; ai
< nral
; ai
++)
2509 int cgi
= at2cg
[il
->iatoms
[i
+1+ai
]];
2510 for (int aj
= ai
+ 1; aj
< nral
; aj
++)
2512 int cgj
= at2cg
[il
->iatoms
[i
+1+aj
]];
2515 real rij2
= distance2(cg_cm
[cgi
], cg_cm
[cgj
]);
2517 update_max_bonded_distance(rij2
, ftype
,
2520 (nral
== 2) ? bd_2b
: bd_mb
);
2530 const t_blocka
*excls
= &molt
->excls
;
2531 for (int ai
= 0; ai
< excls
->nr
; ai
++)
2533 int cgi
= at2cg
[ai
];
2534 for (int j
= excls
->index
[ai
]; j
< excls
->index
[ai
+1]; j
++)
2536 int cgj
= at2cg
[excls
->a
[j
]];
2539 real rij2
= distance2(cg_cm
[cgi
], cg_cm
[cgj
]);
2541 /* There is no function type for exclusions, use -1 */
2542 update_max_bonded_distance(rij2
, -1, ai
, excls
->a
[j
], bd_2b
);
2549 /*! \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 */
2550 static void bonded_distance_intermol(const t_ilist
*ilists_intermol
,
2552 const rvec
*x
, int ePBC
, const matrix box
,
2553 bonded_distance_t
*bd_2b
,
2554 bonded_distance_t
*bd_mb
)
2558 set_pbc(&pbc
, ePBC
, box
);
2560 for (int ftype
= 0; ftype
< F_NRE
; ftype
++)
2562 if (dd_check_ftype(ftype
, bBCheck
, FALSE
, FALSE
))
2564 const t_ilist
*il
= &ilists_intermol
[ftype
];
2565 int nral
= NRAL(ftype
);
2567 /* No nral>1 check here, since intermol interactions always
2568 * have nral>=2 (and the code is also correct for nral=1).
2570 for (int i
= 0; i
< il
->nr
; i
+= 1+nral
)
2572 for (int ai
= 0; ai
< nral
; ai
++)
2574 int atom_i
= il
->iatoms
[i
+ 1 + ai
];
2576 for (int aj
= ai
+ 1; aj
< nral
; aj
++)
2581 int atom_j
= il
->iatoms
[i
+ 1 + aj
];
2583 pbc_dx(&pbc
, x
[atom_i
], x
[atom_j
], dx
);
2587 update_max_bonded_distance(rij2
, ftype
,
2589 (nral
== 2) ? bd_2b
: bd_mb
);
2597 //! Returns whether \p molt has at least one virtual site
2598 static bool moltypeHasVsite(const gmx_moltype_t
&molt
)
2600 bool hasVsite
= false;
2601 for (int i
= 0; i
< F_NRE
; i
++)
2603 if ((interaction_function
[i
].flags
& IF_VSITE
) &&
2604 molt
.ilist
[i
].nr
> 0)
2613 //! Compute charge group centers of mass for molecule \p molt
2614 static void get_cgcm_mol(const gmx_moltype_t
*molt
,
2615 const gmx_ffparams_t
*ffparams
,
2616 int ePBC
, t_graph
*graph
, const matrix box
,
2617 const rvec
*x
, rvec
*xs
, rvec
*cg_cm
)
2621 if (ePBC
!= epbcNONE
)
2623 mk_mshift(nullptr, graph
, ePBC
, box
, x
);
2625 shift_x(graph
, box
, x
, xs
);
2626 /* By doing an extra mk_mshift the molecules that are broken
2627 * because they were e.g. imported from another software
2628 * will be made whole again. Such are the healing powers
2631 mk_mshift(nullptr, graph
, ePBC
, box
, xs
);
2635 /* We copy the coordinates so the original coordinates remain
2636 * unchanged, just to be 100% sure that we do not affect
2637 * binary reproducibility of simulations.
2639 n
= molt
->cgs
.index
[molt
->cgs
.nr
];
2640 for (i
= 0; i
< n
; i
++)
2642 copy_rvec(x
[i
], xs
[i
]);
2646 if (moltypeHasVsite(*molt
))
2648 construct_vsites(nullptr, xs
, 0.0, nullptr,
2649 ffparams
->iparams
, molt
->ilist
,
2650 epbcNONE
, TRUE
, nullptr, nullptr);
2653 calc_cgcm(nullptr, 0, molt
->cgs
.nr
, &molt
->cgs
, xs
, cg_cm
);
2656 void dd_bonded_cg_distance(FILE *fplog
,
2657 const gmx_mtop_t
*mtop
,
2658 const t_inputrec
*ir
,
2659 const rvec
*x
, const matrix box
,
2661 real
*r_2b
, real
*r_mb
)
2663 gmx_bool bExclRequired
;
2667 bonded_distance_t bd_2b
= { 0, -1, -1, -1 };
2668 bonded_distance_t bd_mb
= { 0, -1, -1, -1 };
2670 bExclRequired
= inputrecExclForces(ir
);
2675 for (const gmx_molblock_t
&molb
: mtop
->molblock
)
2677 const gmx_moltype_t
&molt
= mtop
->moltype
[molb
.type
];
2678 if (molt
.cgs
.nr
== 1 || molb
.nmol
== 0)
2680 at_offset
+= molb
.nmol
*molt
.atoms
.nr
;
2684 if (ir
->ePBC
!= epbcNONE
)
2686 mk_graph_ilist(nullptr, molt
.ilist
, 0, molt
.atoms
.nr
, FALSE
, FALSE
,
2690 std::vector
<int> at2cg
= make_at2cg(molt
.cgs
);
2691 snew(xs
, molt
.atoms
.nr
);
2692 snew(cg_cm
, molt
.cgs
.nr
);
2693 for (int mol
= 0; mol
< molb
.nmol
; mol
++)
2695 get_cgcm_mol(&molt
, &mtop
->ffparams
, ir
->ePBC
, &graph
, box
,
2696 x
+at_offset
, xs
, cg_cm
);
2698 bonded_distance_t bd_mol_2b
= { 0, -1, -1, -1 };
2699 bonded_distance_t bd_mol_mb
= { 0, -1, -1, -1 };
2701 bonded_cg_distance_mol(&molt
, at2cg
, bBCheck
, bExclRequired
, cg_cm
,
2702 &bd_mol_2b
, &bd_mol_mb
);
2704 /* Process the mol data adding the atom index offset */
2705 update_max_bonded_distance(bd_mol_2b
.r2
, bd_mol_2b
.ftype
,
2706 at_offset
+ bd_mol_2b
.a1
,
2707 at_offset
+ bd_mol_2b
.a2
,
2709 update_max_bonded_distance(bd_mol_mb
.r2
, bd_mol_mb
.ftype
,
2710 at_offset
+ bd_mol_mb
.a1
,
2711 at_offset
+ bd_mol_mb
.a2
,
2714 at_offset
+= molt
.atoms
.nr
;
2718 if (ir
->ePBC
!= epbcNONE
)
2725 if (mtop
->bIntermolecularInteractions
)
2727 if (ncg_mtop(mtop
) < mtop
->natoms
)
2729 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.");
2732 bonded_distance_intermol(mtop
->intermolecular_ilist
,
2738 *r_2b
= sqrt(bd_2b
.r2
);
2739 *r_mb
= sqrt(bd_mb
.r2
);
2741 if (fplog
&& (*r_2b
> 0 || *r_mb
> 0))
2744 "Initial maximum inter charge-group distances:\n");
2748 " two-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2749 *r_2b
, (bd_2b
.ftype
>= 0) ? interaction_function
[bd_2b
.ftype
].longname
: "Exclusion",
2750 bd_2b
.a1
+ 1, bd_2b
.a2
+ 1);
2755 " multi-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2756 *r_mb
, interaction_function
[bd_mb
.ftype
].longname
,
2757 bd_mb
.a1
+ 1, bd_mb
.a2
+ 1);