2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006,2007,2008,2009,2010,2011,2012,2013,2014,2015, 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
54 #include "gromacs/domdec/domdec.h"
55 #include "gromacs/domdec/domdec_network.h"
56 #include "gromacs/domdec/ga2la.h"
57 #include "gromacs/gmxlib/chargegroup.h"
58 #include "gromacs/gmxlib/gmx_omp_nthreads.h"
59 #include "gromacs/gmxlib/network.h"
60 #include "gromacs/legacyheaders/types/commrec.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/mdlib/force.h"
63 #include "gromacs/mdlib/forcerec.h"
64 #include "gromacs/mdlib/vsite.h"
65 #include "gromacs/mdtypes/inputrec.h"
66 #include "gromacs/mdtypes/md_enums.h"
67 #include "gromacs/mdtypes/state.h"
68 #include "gromacs/pbcutil/mshift.h"
69 #include "gromacs/pbcutil/pbc.h"
70 #include "gromacs/topology/mtop_util.h"
71 #include "gromacs/topology/topsort.h"
72 #include "gromacs/utility/cstringutil.h"
73 #include "gromacs/utility/exceptions.h"
74 #include "gromacs/utility/fatalerror.h"
75 #include "gromacs/utility/gmxassert.h"
76 #include "gromacs/utility/smalloc.h"
78 #include "domdec_constraints.h"
79 #include "domdec_internal.h"
80 #include "domdec_vsite.h"
82 /*! \brief The number of integer item in the local state, used for broadcasting of the state */
83 #define NITEM_DD_INIT_LOCAL_STATE 5
86 int *index
; /* Index for each atom into il */
87 int *il
; /* ftype|type|a0|...|an|ftype|... */
97 /*! \brief Struct for thread local work data for local topology generation */
99 t_idef idef
; /**< Partial local topology */
100 int **vsite_pbc
; /**< vsite PBC structure */
101 int *vsite_pbc_nalloc
; /**< Allocation sizes for vsite_pbc */
102 int nbonded
; /**< The number of bondeds in this struct */
103 t_blocka excl
; /**< List of exclusions */
104 int excl_count
; /**< The total exclusion count for \p excl */
107 /*! \brief Struct for the reverse topology: links bonded interactions to atomsx */
108 struct gmx_reverse_top_t
110 //! @cond Doxygen_Suppress
111 gmx_bool bExclRequired
; /**< Do we require all exclusions to be assigned? */
112 int n_excl_at_max
; /**< The maximum number of exclusions one atom can have */
113 gmx_bool bConstr
; /**< Are there constraints in this revserse top? */
114 gmx_bool bSettle
; /**< Are there settles in this revserse top? */
115 gmx_bool bBCheck
; /**< All bonded interactions have to be assigned? */
116 gmx_bool bInterCGInteractions
; /**< Are there bondeds/exclusions between charge-groups? */
117 reverse_ilist_t
*ril_mt
; /**< Reverse ilist for all moltypes */
118 int ril_mt_tot_size
; /**< The size of ril_mt[?].index summed over all entries */
119 int ilsort
; /**< The sorting state of bondeds for free energy */
120 molblock_ind_t
*mbi
; /**< molblock to global atom index for quick lookup of molblocks on atom index */
121 int nmolblock
; /**< The number of molblocks, size of \p mbi */
123 gmx_bool bIntermolecularInteractions
; /**< Do we have intermolecular interactions? */
124 reverse_ilist_t ril_intermol
; /**< Intermolecular reverse ilist */
126 /* Work data structures for multi-threading */
127 int nthread
; /**< The number of threads to be used */
128 thread_work_t
*th_work
; /**< Thread work array for local topology generation */
130 /* Pointers only used for an error message */
131 gmx_mtop_t
*err_top_global
; /**< Pointer to the global top, only used for error reporting */
132 gmx_localtop_t
*err_top_local
; /**< Pointer to the local top, only used for error reporting */
136 /*! \brief Returns the number of atom entries for il in gmx_reverse_top_t */
137 static int nral_rt(int ftype
)
142 if (interaction_function
[ftype
].flags
& IF_VSITE
)
144 /* With vsites the reverse topology contains
145 * two extra entries for PBC.
153 /*! \brief Return whether interactions of type \p ftype need to be assigned exactly once */
154 static gmx_bool
dd_check_ftype(int ftype
, gmx_bool bBCheck
,
155 gmx_bool bConstr
, gmx_bool bSettle
)
157 return (((interaction_function
[ftype
].flags
& IF_BOND
) &&
158 !(interaction_function
[ftype
].flags
& IF_VSITE
) &&
159 (bBCheck
|| !(interaction_function
[ftype
].flags
& IF_LIMZERO
))) ||
160 (bConstr
&& (ftype
== F_CONSTR
|| ftype
== F_CONSTRNC
)) ||
161 (bSettle
&& ftype
== F_SETTLE
));
164 /*! \brief Print a header on error messages */
165 static void print_error_header(FILE *fplog
, char *moltypename
, int nprint
)
167 fprintf(fplog
, "\nMolecule type '%s'\n", moltypename
);
168 fprintf(stderr
, "\nMolecule type '%s'\n", moltypename
);
170 "the first %d missing interactions, except for exclusions:\n",
173 "the first %d missing interactions, except for exclusions:\n",
177 /*! \brief Help print error output when interactions are missing */
178 static void print_missing_interactions_mb(FILE *fplog
, t_commrec
*cr
,
179 gmx_reverse_top_t
*rt
,
181 reverse_ilist_t
*ril
,
182 int a_start
, int a_end
,
183 int nat_mol
, int nmol
,
186 int nril_mol
, *assigned
, *gatindex
;
187 int ftype
, ftype_j
, nral
, i
, j_mol
, j
, a0
, a0_mol
, mol
, a
;
193 nril_mol
= ril
->index
[nat_mol
];
194 snew(assigned
, nmol
*nril_mol
);
196 gatindex
= cr
->dd
->gatindex
;
197 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
199 if (dd_check_ftype(ftype
, rt
->bBCheck
, rt
->bConstr
, rt
->bSettle
))
202 il
= &idef
->il
[ftype
];
204 for (i
= 0; i
< il
->nr
; i
+= 1+nral
)
206 a0
= gatindex
[ia
[1]];
207 /* Check if this interaction is in
208 * the currently checked molblock.
210 if (a0
>= a_start
&& a0
< a_end
)
212 mol
= (a0
- a_start
)/nat_mol
;
213 a0_mol
= (a0
- a_start
) - mol
*nat_mol
;
214 j_mol
= ril
->index
[a0_mol
];
216 while (j_mol
< ril
->index
[a0_mol
+1] && !bFound
)
218 j
= mol
*nril_mol
+ j_mol
;
219 ftype_j
= ril
->il
[j_mol
];
220 /* Here we need to check if this interaction has
221 * not already been assigned, since we could have
222 * multiply defined interactions.
224 if (ftype
== ftype_j
&& ia
[0] == ril
->il
[j_mol
+1] &&
227 /* Check the atoms */
229 for (a
= 0; a
< nral
; a
++)
231 if (gatindex
[ia
[1+a
]] !=
232 a_start
+ mol
*nat_mol
+ ril
->il
[j_mol
+2+a
])
242 j_mol
+= 2 + nral_rt(ftype_j
);
246 gmx_incons("Some interactions seem to be assigned multiple times");
254 gmx_sumi(nmol
*nril_mol
, assigned
, cr
);
258 for (mol
= 0; mol
< nmol
; mol
++)
261 while (j_mol
< nril_mol
)
263 ftype
= ril
->il
[j_mol
];
265 j
= mol
*nril_mol
+ j_mol
;
266 if (assigned
[j
] == 0 &&
267 !(interaction_function
[ftype
].flags
& IF_VSITE
))
269 if (DDMASTER(cr
->dd
))
273 print_error_header(fplog
, moltypename
, nprint
);
275 fprintf(fplog
, "%20s atoms",
276 interaction_function
[ftype
].longname
);
277 fprintf(stderr
, "%20s atoms",
278 interaction_function
[ftype
].longname
);
279 for (a
= 0; a
< nral
; a
++)
281 fprintf(fplog
, "%5d", ril
->il
[j_mol
+2+a
]+1);
282 fprintf(stderr
, "%5d", ril
->il
[j_mol
+2+a
]+1);
287 fprintf(stderr
, " ");
290 fprintf(fplog
, " global");
291 fprintf(stderr
, " global");
292 for (a
= 0; a
< nral
; a
++)
294 fprintf(fplog
, "%6d",
295 a_start
+mol
*nat_mol
+ril
->il
[j_mol
+2+a
]+1);
296 fprintf(stderr
, "%6d",
297 a_start
+mol
*nat_mol
+ril
->il
[j_mol
+2+a
]+1);
299 fprintf(fplog
, "\n");
300 fprintf(stderr
, "\n");
308 j_mol
+= 2 + nral_rt(ftype
);
315 /*! \brief Help print error output when interactions are missing */
316 static void print_missing_interactions_atoms(FILE *fplog
, t_commrec
*cr
,
317 gmx_mtop_t
*mtop
, t_idef
*idef
)
319 int mb
, a_start
, a_end
;
320 gmx_molblock_t
*molb
;
321 gmx_reverse_top_t
*rt
;
323 rt
= cr
->dd
->reverse_top
;
325 /* Print the atoms in the missing interactions per molblock */
327 for (mb
= 0; mb
< mtop
->nmolblock
; mb
++)
329 molb
= &mtop
->molblock
[mb
];
331 a_end
= a_start
+ molb
->nmol
*molb
->natoms_mol
;
333 print_missing_interactions_mb(fplog
, cr
, rt
,
334 *(mtop
->moltype
[molb
->type
].name
),
335 &rt
->ril_mt
[molb
->type
],
336 a_start
, a_end
, molb
->natoms_mol
,
342 void dd_print_missing_interactions(FILE *fplog
, t_commrec
*cr
, int local_count
, gmx_mtop_t
*top_global
, t_state
*state_local
)
344 int ndiff_tot
, cl
[F_NRE
], n
, ndiff
, rest_global
, rest_local
;
348 gmx_mtop_t
*err_top_global
;
349 gmx_localtop_t
*err_top_local
;
353 err_top_global
= dd
->reverse_top
->err_top_global
;
354 err_top_local
= dd
->reverse_top
->err_top_local
;
358 fprintf(fplog
, "\nNot all bonded interactions have been properly assigned to the domain decomposition cells\n");
362 ndiff_tot
= local_count
- dd
->nbonded_global
;
364 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
367 cl
[ftype
] = err_top_local
->idef
.il
[ftype
].nr
/(1+nral
);
370 gmx_sumi(F_NRE
, cl
, cr
);
376 fprintf(fplog
, "\nA list of missing interactions:\n");
378 fprintf(stderr
, "\nA list of missing interactions:\n");
379 rest_global
= dd
->nbonded_global
;
380 rest_local
= local_count
;
381 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
383 /* In the reverse and local top all constraints are merged
384 * into F_CONSTR. So in the if statement we skip F_CONSTRNC
385 * and add these constraints when doing F_CONSTR.
387 if (((interaction_function
[ftype
].flags
& IF_BOND
) &&
388 (dd
->reverse_top
->bBCheck
389 || !(interaction_function
[ftype
].flags
& IF_LIMZERO
)))
390 || (dd
->reverse_top
->bConstr
&& ftype
== F_CONSTR
)
391 || (dd
->reverse_top
->bSettle
&& ftype
== F_SETTLE
))
393 n
= gmx_mtop_ftype_count(err_top_global
, ftype
);
394 if (ftype
== F_CONSTR
)
396 n
+= gmx_mtop_ftype_count(err_top_global
, F_CONSTRNC
);
398 ndiff
= cl
[ftype
] - n
;
401 sprintf(buf
, "%20s of %6d missing %6d",
402 interaction_function
[ftype
].longname
, n
, -ndiff
);
405 fprintf(fplog
, "%s\n", buf
);
407 fprintf(stderr
, "%s\n", buf
);
410 rest_local
-= cl
[ftype
];
414 ndiff
= rest_local
- rest_global
;
417 sprintf(buf
, "%20s of %6d missing %6d", "exclusions",
418 rest_global
, -ndiff
);
421 fprintf(fplog
, "%s\n", buf
);
423 fprintf(stderr
, "%s\n", buf
);
427 print_missing_interactions_atoms(fplog
, cr
, err_top_global
,
428 &err_top_local
->idef
);
429 write_dd_pdb("dd_dump_err", 0, "dump", top_global
, cr
,
430 -1, state_local
->x
, state_local
->box
);
435 gmx_incons("One or more interactions were multiple assigned in the domain decompostion");
439 gmx_fatal(FARGS
, "%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(cr
->dd
), dd_cutoff_twobody(cr
->dd
));
444 /*! \brief Return global topology molecule information for global atom index \p i_gl */
445 static void global_atomnr_to_moltype_ind(gmx_reverse_top_t
*rt
, int i_gl
,
446 int *mb
, int *mt
, int *mol
, int *i_mol
)
448 molblock_ind_t
*mbi
= rt
->mbi
;
450 int end
= rt
->nmolblock
; /* exclusive */
453 /* binary search for molblock_ind */
457 if (i_gl
>= mbi
[mid
].a_end
)
461 else if (i_gl
< mbi
[mid
].a_start
)
475 *mol
= (i_gl
- mbi
->a_start
) / mbi
->natoms_mol
;
476 *i_mol
= (i_gl
- mbi
->a_start
) - (*mol
)*mbi
->natoms_mol
;
479 /*! \brief Count the exclusions for all atoms in \p cgs */
480 static void count_excls(const t_block
*cgs
, const t_blocka
*excls
,
481 int *n_excl
, int *n_intercg_excl
, int *n_excl_at_max
)
483 int cg
, at0
, at1
, at
, excl
, atj
;
488 for (cg
= 0; cg
< cgs
->nr
; cg
++)
490 at0
= cgs
->index
[cg
];
491 at1
= cgs
->index
[cg
+1];
492 for (at
= at0
; at
< at1
; at
++)
494 for (excl
= excls
->index
[at
]; excl
< excls
->index
[at
+1]; excl
++)
496 atj
= excls
->a
[excl
];
500 if (atj
< at0
|| atj
>= at1
)
507 *n_excl_at_max
= std::max(*n_excl_at_max
,
508 excls
->index
[at
+1] - excls
->index
[at
]);
513 /*! \brief Run the reverse ilist generation and store it when \p bAssign = TRUE */
514 static int low_make_reverse_ilist(const t_ilist
*il_mt
,
516 int **vsite_pbc
, /* should be const */
518 gmx_bool bConstr
, gmx_bool bSettle
,
520 int *r_index
, int *r_il
,
521 gmx_bool bLinkToAllAtoms
,
524 int ftype
, nral
, i
, j
, nlink
, link
;
532 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
534 if ((interaction_function
[ftype
].flags
& (IF_BOND
| IF_VSITE
)) ||
535 (bConstr
&& (ftype
== F_CONSTR
|| ftype
== F_CONSTRNC
)) ||
536 (bSettle
&& ftype
== F_SETTLE
))
538 bVSite
= (interaction_function
[ftype
].flags
& IF_VSITE
);
541 for (i
= 0; i
< il
->nr
; i
+= 1+nral
)
548 /* We don't need the virtual sites for the cg-links */
558 /* Couple to the first atom in the interaction */
561 for (link
= 0; link
< nlink
; link
++)
566 r_il
[r_index
[a
]+count
[a
]] =
567 (ftype
== F_CONSTRNC
? F_CONSTR
: ftype
);
568 r_il
[r_index
[a
]+count
[a
]+1] = ia
[0];
569 for (j
= 1; j
< 1+nral
; j
++)
571 /* Store the molecular atom number */
572 r_il
[r_index
[a
]+count
[a
]+1+j
] = ia
[j
];
575 if (interaction_function
[ftype
].flags
& IF_VSITE
)
579 /* Add an entry to iatoms for storing
580 * which of the constructing atoms are
583 r_il
[r_index
[a
]+count
[a
]+2+nral
] = 0;
584 for (j
= 2; j
< 1+nral
; j
++)
586 if (atom
[ia
[j
]].ptype
== eptVSite
)
588 r_il
[r_index
[a
]+count
[a
]+2+nral
] |= (2<<j
);
591 /* Store vsite pbc atom in a second extra entry */
592 r_il
[r_index
[a
]+count
[a
]+2+nral
+1] =
593 (vsite_pbc
? vsite_pbc
[ftype
-F_VSITE2
][i
/(1+nral
)] : -2);
598 /* We do not count vsites since they are always
599 * uniquely assigned and can be assigned
600 * to multiple nodes with recursive vsites.
603 !(interaction_function
[ftype
].flags
& IF_LIMZERO
))
608 count
[a
] += 2 + nral_rt(ftype
);
617 /*! \brief Make the reverse ilist: a list of bonded interactions linked to atoms */
618 static int make_reverse_ilist(const t_ilist
*ilist
,
619 const t_atoms
*atoms
,
620 int **vsite_pbc
, /* should be const (C issue) */
621 gmx_bool bConstr
, gmx_bool bSettle
,
623 gmx_bool bLinkToAllAtoms
,
624 reverse_ilist_t
*ril_mt
)
626 int nat_mt
, *count
, i
, nint_mt
;
628 /* Count the interactions */
631 low_make_reverse_ilist(ilist
, atoms
->atom
, vsite_pbc
,
633 bConstr
, bSettle
, bBCheck
, NULL
, NULL
,
634 bLinkToAllAtoms
, FALSE
);
636 snew(ril_mt
->index
, nat_mt
+1);
637 ril_mt
->index
[0] = 0;
638 for (i
= 0; i
< nat_mt
; i
++)
640 ril_mt
->index
[i
+1] = ril_mt
->index
[i
] + count
[i
];
643 snew(ril_mt
->il
, ril_mt
->index
[nat_mt
]);
645 /* Store the interactions */
647 low_make_reverse_ilist(ilist
, atoms
->atom
, vsite_pbc
,
649 bConstr
, bSettle
, bBCheck
,
650 ril_mt
->index
, ril_mt
->il
,
651 bLinkToAllAtoms
, TRUE
);
658 /*! \brief Destroys a reverse_ilist_t struct */
659 static void destroy_reverse_ilist(reverse_ilist_t
*ril
)
665 /*! \brief Generate the reverse topology */
666 static gmx_reverse_top_t
*make_reverse_top(gmx_mtop_t
*mtop
, gmx_bool bFE
,
667 int ***vsite_pbc_molt
,
668 gmx_bool bConstr
, gmx_bool bSettle
,
669 gmx_bool bBCheck
, int *nint
)
672 gmx_reverse_top_t
*rt
;
679 /* Should we include constraints (for SHAKE) in rt? */
680 rt
->bConstr
= bConstr
;
681 rt
->bSettle
= bSettle
;
682 rt
->bBCheck
= bBCheck
;
684 rt
->bInterCGInteractions
= mtop
->bIntermolecularInteractions
;
685 snew(nint_mt
, mtop
->nmoltype
);
686 snew(rt
->ril_mt
, mtop
->nmoltype
);
687 rt
->ril_mt_tot_size
= 0;
688 for (mt
= 0; mt
< mtop
->nmoltype
; mt
++)
690 molt
= &mtop
->moltype
[mt
];
691 if (molt
->cgs
.nr
> 1)
693 rt
->bInterCGInteractions
= TRUE
;
696 /* Make the atom to interaction list for this molecule type */
698 make_reverse_ilist(molt
->ilist
, &molt
->atoms
,
699 vsite_pbc_molt
? vsite_pbc_molt
[mt
] : NULL
,
700 rt
->bConstr
, rt
->bSettle
, rt
->bBCheck
, FALSE
,
703 rt
->ril_mt_tot_size
+= rt
->ril_mt
[mt
].index
[molt
->atoms
.nr
];
707 fprintf(debug
, "The total size of the atom to interaction index is %d integers\n", rt
->ril_mt_tot_size
);
711 for (mb
= 0; mb
< mtop
->nmolblock
; mb
++)
713 *nint
+= mtop
->molblock
[mb
].nmol
*nint_mt
[mtop
->molblock
[mb
].type
];
717 /* Make an intermolecular reverse top, if necessary */
718 rt
->bIntermolecularInteractions
= mtop
->bIntermolecularInteractions
;
719 if (rt
->bIntermolecularInteractions
)
721 t_atoms atoms_global
;
723 rt
->ril_intermol
.index
= NULL
;
724 rt
->ril_intermol
.il
= NULL
;
726 atoms_global
.nr
= mtop
->natoms
;
727 atoms_global
.atom
= NULL
; /* Only used with virtual sites */
730 make_reverse_ilist(mtop
->intermolecular_ilist
, &atoms_global
,
732 rt
->bConstr
, rt
->bSettle
, rt
->bBCheck
, FALSE
,
736 if (bFE
&& gmx_mtop_bondeds_free_energy(mtop
))
738 rt
->ilsort
= ilsortFE_UNSORTED
;
742 rt
->ilsort
= ilsortNO_FE
;
745 /* Make a molblock index for fast searching */
746 snew(rt
->mbi
, mtop
->nmolblock
);
747 rt
->nmolblock
= mtop
->nmolblock
;
749 for (mb
= 0; mb
< mtop
->nmolblock
; mb
++)
751 rt
->mbi
[mb
].a_start
= i
;
752 i
+= mtop
->molblock
[mb
].nmol
*mtop
->molblock
[mb
].natoms_mol
;
753 rt
->mbi
[mb
].a_end
= i
;
754 rt
->mbi
[mb
].natoms_mol
= mtop
->molblock
[mb
].natoms_mol
;
755 rt
->mbi
[mb
].type
= mtop
->molblock
[mb
].type
;
758 rt
->nthread
= gmx_omp_nthreads_get(emntDomdec
);
759 snew(rt
->th_work
, rt
->nthread
);
760 if (vsite_pbc_molt
!= NULL
)
762 for (thread
= 0; thread
< rt
->nthread
; thread
++)
764 snew(rt
->th_work
[thread
].vsite_pbc
, F_VSITEN
-F_VSITE2
+1);
765 snew(rt
->th_work
[thread
].vsite_pbc_nalloc
, F_VSITEN
-F_VSITE2
+1);
772 void dd_make_reverse_top(FILE *fplog
,
773 gmx_domdec_t
*dd
, gmx_mtop_t
*mtop
,
775 t_inputrec
*ir
, gmx_bool bBCheck
)
779 fprintf(fplog
, "\nLinking all bonded interactions to atoms\n");
782 /* If normal and/or settle constraints act only within charge groups,
783 * we can store them in the reverse top and simply assign them to domains.
784 * Otherwise we need to assign them to multiple domains and set up
785 * the parallel version constraint algorithm(s).
788 dd
->reverse_top
= make_reverse_top(mtop
, ir
->efep
!= efepNO
,
789 vsite
? vsite
->vsite_pbc_molt
: NULL
,
790 !dd
->bInterCGcons
, !dd
->bInterCGsettles
,
791 bBCheck
, &dd
->nbonded_global
);
793 gmx_reverse_top_t
*rt
= dd
->reverse_top
;
795 if (rt
->ril_mt_tot_size
>= 200000 &&
797 mtop
->nmolblock
== 1 && mtop
->molblock
[0].nmol
== 1)
799 /* mtop comes from a pre Gromacs 4 tpr file */
800 const char *note
= "NOTE: The tpr file used for this simulation is in an old format, for less memory usage and possibly more performance create a new tpr file with an up to date version of grompp";
803 fprintf(fplog
, "\n%s\n\n", note
);
807 fprintf(stderr
, "\n%s\n\n", note
);
811 /* With the Verlet scheme, exclusions are handled in the non-bonded
812 * kernels and only exclusions inside the cut-off lead to exclusion
813 * forces. Since each atom pair is treated at most once in the non-bonded
814 * kernels, it doesn't matter if the exclusions for the same atom pair
815 * appear multiple times in the exclusion list. In contrast, the, old,
816 * group cut-off scheme loops over a list of exclusions, so there each
817 * excluded pair should appear exactly once.
819 rt
->bExclRequired
= (ir
->cutoff_scheme
== ecutsGROUP
&&
820 inputrecExclForces(ir
));
825 dd
->n_intercg_excl
= 0;
826 rt
->n_excl_at_max
= 0;
827 for (mb
= 0; mb
< mtop
->nmolblock
; mb
++)
829 gmx_molblock_t
*molb
;
831 int n_excl_mol
, n_excl_icg
, n_excl_at_max
;
833 molb
= &mtop
->molblock
[mb
];
834 molt
= &mtop
->moltype
[molb
->type
];
835 count_excls(&molt
->cgs
, &molt
->excls
,
836 &n_excl_mol
, &n_excl_icg
, &n_excl_at_max
);
837 nexcl
+= molb
->nmol
*n_excl_mol
;
838 dd
->n_intercg_excl
+= molb
->nmol
*n_excl_icg
;
839 rt
->n_excl_at_max
= std::max(rt
->n_excl_at_max
, n_excl_at_max
);
841 if (rt
->bExclRequired
)
843 dd
->nbonded_global
+= nexcl
;
844 if (EEL_FULL(ir
->coulombtype
) && dd
->n_intercg_excl
> 0 && fplog
)
846 fprintf(fplog
, "There are %d inter charge-group exclusions,\n"
847 "will use an extra communication step for exclusion forces for %s\n",
848 dd
->n_intercg_excl
, eel_names
[ir
->coulombtype
]);
852 if (vsite
&& vsite
->n_intercg_vsite
> 0)
856 fprintf(fplog
, "There are %d inter charge-group virtual sites,\n"
857 "will an extra communication step for selected coordinates and forces\n",
858 vsite
->n_intercg_vsite
);
860 init_domdec_vsites(dd
, vsite
->n_intercg_vsite
);
863 if (dd
->bInterCGcons
|| dd
->bInterCGsettles
)
865 init_domdec_constraints(dd
, mtop
);
869 fprintf(fplog
, "\n");
873 /*! \brief Store a vsite interaction at the end of \p il
875 * This routine is very similar to add_ifunc, but vsites interactions
876 * have more work to do than other kinds of interactions, and the
877 * complex way nral (and thus vector contents) depends on ftype
878 * confuses static analysis tools unless we fuse the vsite
879 * atom-indexing organization code with the ifunc-adding code, so that
880 * they can see that nral is the same value. */
881 static gmx_inline
void
882 add_ifunc_for_vsites(t_iatom
*tiatoms
, gmx_ga2la_t
*ga2la
,
883 int nral
, gmx_bool bHomeA
,
884 int a
, int a_gl
, int a_mol
,
885 const t_iatom
*iatoms
,
890 if (il
->nr
+1+nral
> il
->nalloc
)
892 il
->nalloc
= over_alloc_large(il
->nr
+1+nral
);
893 srenew(il
->iatoms
, il
->nalloc
);
895 liatoms
= il
->iatoms
+ il
->nr
;
899 tiatoms
[0] = iatoms
[0];
903 /* We know the local index of the first atom */
908 /* Convert later in make_local_vsites */
909 tiatoms
[1] = -a_gl
- 1;
912 for (int k
= 2; k
< 1+nral
; k
++)
914 int ak_gl
= a_gl
+ iatoms
[k
] - a_mol
;
915 if (!ga2la_get_home(ga2la
, ak_gl
, &tiatoms
[k
]))
917 /* Copy the global index, convert later in make_local_vsites */
918 tiatoms
[k
] = -(ak_gl
+ 1);
920 // Note that ga2la_get_home always sets the third parameter if
923 for (int k
= 0; k
< 1+nral
; k
++)
925 liatoms
[k
] = tiatoms
[k
];
929 /*! \brief Store a bonded interaction at the end of \p il */
930 static gmx_inline
void add_ifunc(int nral
, t_iatom
*tiatoms
, t_ilist
*il
)
935 if (il
->nr
+1+nral
> il
->nalloc
)
937 il
->nalloc
= over_alloc_large(il
->nr
+1+nral
);
938 srenew(il
->iatoms
, il
->nalloc
);
940 liatoms
= il
->iatoms
+ il
->nr
;
941 for (k
= 0; k
<= nral
; k
++)
943 liatoms
[k
] = tiatoms
[k
];
948 /*! \brief Store a position restraint in idef and iatoms, complex because the parameters are different for each entry */
949 static void add_posres(int mol
, int a_mol
, const gmx_molblock_t
*molb
,
950 t_iatom
*iatoms
, const t_iparams
*ip_in
,
956 /* This position restraint has not been added yet,
957 * so it's index is the current number of position restraints.
959 n
= idef
->il
[F_POSRES
].nr
/2;
960 if (n
+1 > idef
->iparams_posres_nalloc
)
962 idef
->iparams_posres_nalloc
= over_alloc_dd(n
+1);
963 srenew(idef
->iparams_posres
, idef
->iparams_posres_nalloc
);
965 ip
= &idef
->iparams_posres
[n
];
966 /* Copy the force constants */
967 *ip
= ip_in
[iatoms
[0]];
969 /* Get the position restraint coordinates from the molblock */
970 a_molb
= mol
*molb
->natoms_mol
+ a_mol
;
971 if (a_molb
>= molb
->nposres_xA
)
973 gmx_incons("Not enough position restraint coordinates");
975 ip
->posres
.pos0A
[XX
] = molb
->posres_xA
[a_molb
][XX
];
976 ip
->posres
.pos0A
[YY
] = molb
->posres_xA
[a_molb
][YY
];
977 ip
->posres
.pos0A
[ZZ
] = molb
->posres_xA
[a_molb
][ZZ
];
978 if (molb
->nposres_xB
> 0)
980 ip
->posres
.pos0B
[XX
] = molb
->posres_xB
[a_molb
][XX
];
981 ip
->posres
.pos0B
[YY
] = molb
->posres_xB
[a_molb
][YY
];
982 ip
->posres
.pos0B
[ZZ
] = molb
->posres_xB
[a_molb
][ZZ
];
986 ip
->posres
.pos0B
[XX
] = ip
->posres
.pos0A
[XX
];
987 ip
->posres
.pos0B
[YY
] = ip
->posres
.pos0A
[YY
];
988 ip
->posres
.pos0B
[ZZ
] = ip
->posres
.pos0A
[ZZ
];
990 /* Set the parameter index for idef->iparams_posre */
994 /*! \brief Store a flat-bottomed position restraint in idef and iatoms, complex because the parameters are different for each entry */
995 static void add_fbposres(int mol
, int a_mol
, const gmx_molblock_t
*molb
,
996 t_iatom
*iatoms
, const t_iparams
*ip_in
,
1002 /* This flat-bottom position restraint has not been added yet,
1003 * so it's index is the current number of position restraints.
1005 n
= idef
->il
[F_FBPOSRES
].nr
/2;
1006 if (n
+1 > idef
->iparams_fbposres_nalloc
)
1008 idef
->iparams_fbposres_nalloc
= over_alloc_dd(n
+1);
1009 srenew(idef
->iparams_fbposres
, idef
->iparams_fbposres_nalloc
);
1011 ip
= &idef
->iparams_fbposres
[n
];
1012 /* Copy the force constants */
1013 *ip
= ip_in
[iatoms
[0]];
1015 /* Get the position restriant coordinats from the molblock */
1016 a_molb
= mol
*molb
->natoms_mol
+ a_mol
;
1017 if (a_molb
>= molb
->nposres_xA
)
1019 gmx_incons("Not enough position restraint coordinates");
1021 /* Take reference positions from A position of normal posres */
1022 ip
->fbposres
.pos0
[XX
] = molb
->posres_xA
[a_molb
][XX
];
1023 ip
->fbposres
.pos0
[YY
] = molb
->posres_xA
[a_molb
][YY
];
1024 ip
->fbposres
.pos0
[ZZ
] = molb
->posres_xA
[a_molb
][ZZ
];
1026 /* Note: no B-type for flat-bottom posres */
1028 /* Set the parameter index for idef->iparams_posre */
1032 /*! \brief Store a virtual site interaction, complex because of PBC and recursion */
1033 static void add_vsite(gmx_ga2la_t
*ga2la
, const int *index
, const int *rtil
,
1034 int ftype
, int nral
,
1035 gmx_bool bHomeA
, int a
, int a_gl
, int a_mol
,
1036 const t_iatom
*iatoms
,
1037 t_idef
*idef
, int **vsite_pbc
, int *vsite_pbc_nalloc
)
1039 int k
, vsi
, pbc_a_mol
;
1040 t_iatom tiatoms
[1+MAXATOMLIST
];
1041 int j
, ftype_r
, nral_r
;
1043 /* Add this interaction to the local topology */
1044 add_ifunc_for_vsites(tiatoms
, ga2la
, nral
, bHomeA
, a
, a_gl
, a_mol
, iatoms
, &idef
->il
[ftype
]);
1048 vsi
= idef
->il
[ftype
].nr
/(1+nral
) - 1;
1049 if (vsi
>= vsite_pbc_nalloc
[ftype
-F_VSITE2
])
1051 vsite_pbc_nalloc
[ftype
-F_VSITE2
] = over_alloc_large(vsi
+1);
1052 srenew(vsite_pbc
[ftype
-F_VSITE2
], vsite_pbc_nalloc
[ftype
-F_VSITE2
]);
1056 pbc_a_mol
= iatoms
[1+nral
+1];
1059 /* The pbc flag is one of the following two options:
1060 * -2: vsite and all constructing atoms are within the same cg, no pbc
1061 * -1: vsite and its first constructing atom are in the same cg, do pbc
1063 vsite_pbc
[ftype
-F_VSITE2
][vsi
] = pbc_a_mol
;
1067 /* Set the pbc atom for this vsite so we can make its pbc
1068 * identical to the rest of the atoms in its charge group.
1069 * Since the order of the atoms does not change within a charge
1070 * group, we do not need the global to local atom index.
1072 vsite_pbc
[ftype
-F_VSITE2
][vsi
] = a
+ pbc_a_mol
- iatoms
[1];
1077 /* This vsite is non-home (required for recursion),
1078 * and therefore there is no charge group to match pbc with.
1079 * But we always turn on full_pbc to assure that higher order
1080 * recursion works correctly.
1082 vsite_pbc
[ftype
-F_VSITE2
][vsi
] = -1;
1088 /* Check for recursion */
1089 for (k
= 2; k
< 1+nral
; k
++)
1091 if ((iatoms
[1+nral
] & (2<<k
)) && (tiatoms
[k
] < 0))
1093 /* This construction atoms is a vsite and not a home atom */
1096 fprintf(debug
, "Constructing atom %d of vsite atom %d is a vsite and non-home\n", iatoms
[k
]+1, a_mol
+1);
1098 /* Find the vsite construction */
1100 /* Check all interactions assigned to this atom */
1101 j
= index
[iatoms
[k
]];
1102 while (j
< index
[iatoms
[k
]+1])
1104 ftype_r
= rtil
[j
++];
1105 nral_r
= NRAL(ftype_r
);
1106 if (interaction_function
[ftype_r
].flags
& IF_VSITE
)
1108 /* Add this vsite (recursion) */
1109 add_vsite(ga2la
, index
, rtil
, ftype_r
, nral_r
,
1110 FALSE
, -1, a_gl
+iatoms
[k
]-iatoms
[1], iatoms
[k
],
1111 rtil
+j
, idef
, vsite_pbc
, vsite_pbc_nalloc
);
1112 j
+= 1 + nral_r
+ 2;
1124 /*! \brief Update the local atom to local charge group index */
1125 static void make_la2lc(gmx_domdec_t
*dd
)
1127 int *cgindex
, *la2lc
, cg
, a
;
1129 cgindex
= dd
->cgindex
;
1131 if (dd
->nat_tot
> dd
->la2lc_nalloc
)
1133 dd
->la2lc_nalloc
= over_alloc_dd(dd
->nat_tot
);
1134 snew(dd
->la2lc
, dd
->la2lc_nalloc
);
1138 /* Make the local atom to local cg index */
1139 for (cg
= 0; cg
< dd
->ncg_tot
; cg
++)
1141 for (a
= cgindex
[cg
]; a
< cgindex
[cg
+1]; a
++)
1148 /*! \brief Returns the squared distance between charge groups \p i and \p j */
1149 static real
dd_dist2(t_pbc
*pbc_null
, rvec
*cg_cm
, const int *la2lc
, int i
, int j
)
1155 pbc_dx_aiuc(pbc_null
, cg_cm
[la2lc
[i
]], cg_cm
[la2lc
[j
]], dx
);
1159 rvec_sub(cg_cm
[la2lc
[i
]], cg_cm
[la2lc
[j
]], dx
);
1165 /*! \brief Append t_blocka block structures 1 to nsrc in src to *dest */
1166 static void combine_blocka(t_blocka
*dest
, const thread_work_t
*src
, int nsrc
)
1170 ni
= src
[nsrc
-1].excl
.nr
;
1172 for (s
= 0; s
< nsrc
; s
++)
1174 na
+= src
[s
].excl
.nra
;
1176 if (ni
+ 1 > dest
->nalloc_index
)
1178 dest
->nalloc_index
= over_alloc_large(ni
+1);
1179 srenew(dest
->index
, dest
->nalloc_index
);
1181 if (dest
->nra
+ na
> dest
->nalloc_a
)
1183 dest
->nalloc_a
= over_alloc_large(dest
->nra
+na
);
1184 srenew(dest
->a
, dest
->nalloc_a
);
1186 for (s
= 1; s
< nsrc
; s
++)
1188 for (i
= dest
->nr
+1; i
< src
[s
].excl
.nr
+1; i
++)
1190 dest
->index
[i
] = dest
->nra
+ src
[s
].excl
.index
[i
];
1192 for (i
= 0; i
< src
[s
].excl
.nra
; i
++)
1194 dest
->a
[dest
->nra
+i
] = src
[s
].excl
.a
[i
];
1196 dest
->nr
= src
[s
].excl
.nr
;
1197 dest
->nra
+= src
[s
].excl
.nra
;
1201 /*! \brief Append t_idef structures 1 to nsrc in src to *dest,
1202 * virtual sites need special attention, as pbc info differs per vsite.
1204 static void combine_idef(t_idef
*dest
, const thread_work_t
*src
, int nsrc
,
1209 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
1214 for (s
= 1; s
< nsrc
; s
++)
1216 n
+= src
[s
].idef
.il
[ftype
].nr
;
1222 ild
= &dest
->il
[ftype
];
1224 if (ild
->nr
+ n
> ild
->nalloc
)
1226 ild
->nalloc
= over_alloc_large(ild
->nr
+n
);
1227 srenew(ild
->iatoms
, ild
->nalloc
);
1231 int nral1
= 0, ftv
= 0;
1233 vpbc
= ((interaction_function
[ftype
].flags
& IF_VSITE
) &&
1234 vsite
->vsite_pbc_loc
!= NULL
);
1237 nral1
= 1 + NRAL(ftype
);
1238 ftv
= ftype
- F_VSITE2
;
1239 if ((ild
->nr
+ n
)/nral1
> vsite
->vsite_pbc_loc_nalloc
[ftv
])
1241 vsite
->vsite_pbc_loc_nalloc
[ftv
] =
1242 over_alloc_large((ild
->nr
+ n
)/nral1
);
1243 srenew(vsite
->vsite_pbc_loc
[ftv
],
1244 vsite
->vsite_pbc_loc_nalloc
[ftv
]);
1248 for (s
= 1; s
< nsrc
; s
++)
1253 ils
= &src
[s
].idef
.il
[ftype
];
1254 for (i
= 0; i
< ils
->nr
; i
++)
1256 ild
->iatoms
[ild
->nr
+i
] = ils
->iatoms
[i
];
1260 for (i
= 0; i
< ils
->nr
; i
+= nral1
)
1262 vsite
->vsite_pbc_loc
[ftv
][(ild
->nr
+i
)/nral1
] =
1263 src
[s
].vsite_pbc
[ftv
][i
/nral1
];
1272 /* Position restraints need an additional treatment */
1273 if (dest
->il
[F_POSRES
].nr
> 0)
1277 n
= dest
->il
[F_POSRES
].nr
/2;
1278 if (n
> dest
->iparams_posres_nalloc
)
1280 dest
->iparams_posres_nalloc
= over_alloc_large(n
);
1281 srenew(dest
->iparams_posres
, dest
->iparams_posres_nalloc
);
1283 /* Set n to the number of original position restraints in dest */
1284 for (s
= 1; s
< nsrc
; s
++)
1286 n
-= src
[s
].idef
.il
[F_POSRES
].nr
/2;
1288 for (s
= 1; s
< nsrc
; s
++)
1290 for (i
= 0; i
< src
[s
].idef
.il
[F_POSRES
].nr
/2; i
++)
1292 /* Correct the index into iparams_posres */
1293 dest
->il
[F_POSRES
].iatoms
[n
*2] = n
;
1294 /* Copy the position restraint force parameters */
1295 dest
->iparams_posres
[n
] = src
[s
].idef
.iparams_posres
[i
];
1302 /*! \brief Check and when available assign bonded interactions for local atom i
1304 static gmx_inline
void
1305 check_assign_interactions_atom(int i
, int i_gl
,
1307 const int *index
, const int *rtil
,
1308 gmx_bool bInterMolInteractions
,
1309 int ind_start
, int ind_end
,
1310 const gmx_domdec_t
*dd
,
1311 const gmx_domdec_zones_t
*zones
,
1312 const gmx_molblock_t
*molb
,
1313 gmx_bool bRCheckMB
, ivec rcheck
, gmx_bool bRCheck2B
,
1318 const t_iparams
*ip_in
,
1320 int **vsite_pbc
, int *vsite_pbc_nalloc
,
1331 const t_iatom
*iatoms
;
1333 t_iatom tiatoms
[1 + MAXATOMLIST
];
1338 if (ftype
== F_SETTLE
)
1340 /* Settles are only in the reverse top when they
1341 * operate within a charge group. So we can assign
1342 * them without checks. We do this only for performance
1343 * reasons; it could be handled by the code below.
1347 /* Home zone: add this settle to the local topology */
1348 tiatoms
[0] = iatoms
[0];
1350 tiatoms
[2] = i
+ iatoms
[2] - iatoms
[1];
1351 tiatoms
[3] = i
+ iatoms
[3] - iatoms
[1];
1352 add_ifunc(nral
, tiatoms
, &idef
->il
[ftype
]);
1357 else if (interaction_function
[ftype
].flags
& IF_VSITE
)
1359 assert(!bInterMolInteractions
);
1360 /* The vsite construction goes where the vsite itself is */
1363 add_vsite(dd
->ga2la
, index
, rtil
, ftype
, nral
,
1364 TRUE
, i
, i_gl
, i_mol
,
1365 iatoms
, idef
, vsite_pbc
, vsite_pbc_nalloc
);
1374 tiatoms
[0] = iatoms
[0];
1378 assert(!bInterMolInteractions
);
1379 /* Assign single-body interactions to the home zone */
1384 if (ftype
== F_POSRES
)
1386 add_posres(mol
, i_mol
, molb
, tiatoms
, ip_in
,
1389 else if (ftype
== F_FBPOSRES
)
1391 add_fbposres(mol
, i_mol
, molb
, tiatoms
, ip_in
,
1402 /* This is a two-body interaction, we can assign
1403 * analogous to the non-bonded assignments.
1405 int k_gl
, a_loc
, kz
;
1407 if (!bInterMolInteractions
)
1409 /* Get the global index using the offset in the molecule */
1410 k_gl
= i_gl
+ iatoms
[2] - i_mol
;
1416 if (!ga2la_get(dd
->ga2la
, k_gl
, &a_loc
, &kz
))
1426 /* Check zone interaction assignments */
1427 bUse
= ((iz
< zones
->nizone
&&
1429 kz
>= zones
->izone
[iz
].j0
&&
1430 kz
< zones
->izone
[iz
].j1
) ||
1431 (kz
< zones
->nizone
&&
1433 iz
>= zones
->izone
[kz
].j0
&&
1434 iz
< zones
->izone
[kz
].j1
));
1439 /* If necessary check the cgcm distance */
1441 dd_dist2(pbc_null
, cg_cm
, la2lc
,
1442 tiatoms
[1], tiatoms
[2]) >= rc2
)
1451 /* Assign this multi-body bonded interaction to
1452 * the local node if we have all the atoms involved
1453 * (local or communicated) and the minimum zone shift
1454 * in each dimension is zero, for dimensions
1455 * with 2 DD cells an extra check may be necessary.
1457 ivec k_zero
, k_plus
;
1463 for (k
= 1; k
<= nral
&& bUse
; k
++)
1469 if (!bInterMolInteractions
)
1471 /* Get the global index using the offset in the molecule */
1472 k_gl
= i_gl
+ iatoms
[k
] - i_mol
;
1478 bLocal
= ga2la_get(dd
->ga2la
, k_gl
, &a_loc
, &kz
);
1479 if (!bLocal
|| kz
>= zones
->n
)
1481 /* We do not have this atom of this interaction
1482 * locally, or it comes from more than one cell
1492 for (d
= 0; d
< DIM
; d
++)
1494 if (zones
->shift
[kz
][d
] == 0)
1506 k_zero
[XX
] && k_zero
[YY
] && k_zero
[ZZ
]);
1511 for (d
= 0; (d
< DIM
&& bUse
); d
++)
1513 /* Check if the cg_cm distance falls within
1514 * the cut-off to avoid possible multiple
1515 * assignments of bonded interactions.
1519 dd_dist2(pbc_null
, cg_cm
, la2lc
,
1520 tiatoms
[k_zero
[d
]], tiatoms
[k_plus
[d
]]) >= rc2
)
1529 /* Add this interaction to the local topology */
1530 add_ifunc(nral
, tiatoms
, &idef
->il
[ftype
]);
1531 /* Sum so we can check in global_stat
1532 * if we have everything.
1535 !(interaction_function
[ftype
].flags
& IF_LIMZERO
))
1545 /*! \brief This function looks up and assigns bonded interactions for zone iz.
1547 * With thread parallelizing each thread acts on a different atom range:
1548 * at_start to at_end.
1550 static int make_bondeds_zone(gmx_domdec_t
*dd
,
1551 const gmx_domdec_zones_t
*zones
,
1552 const gmx_molblock_t
*molb
,
1553 gmx_bool bRCheckMB
, ivec rcheck
, gmx_bool bRCheck2B
,
1555 int *la2lc
, t_pbc
*pbc_null
, rvec
*cg_cm
,
1556 const t_iparams
*ip_in
,
1559 int *vsite_pbc_nalloc
,
1561 int at_start
, int at_end
)
1563 int i
, i_gl
, mb
, mt
, mol
, i_mol
;
1566 gmx_reverse_top_t
*rt
;
1569 rt
= dd
->reverse_top
;
1571 bBCheck
= rt
->bBCheck
;
1575 for (i
= at_start
; i
< at_end
; i
++)
1577 /* Get the global atom number */
1578 i_gl
= dd
->gatindex
[i
];
1579 global_atomnr_to_moltype_ind(rt
, i_gl
, &mb
, &mt
, &mol
, &i_mol
);
1580 /* Check all intramolecular interactions assigned to this atom */
1581 index
= rt
->ril_mt
[mt
].index
;
1582 rtil
= rt
->ril_mt
[mt
].il
;
1584 check_assign_interactions_atom(i
, i_gl
, mol
, i_mol
,
1586 index
[i_mol
], index
[i_mol
+1],
1589 bRCheckMB
, rcheck
, bRCheck2B
, rc2
,
1594 idef
, vsite_pbc
, vsite_pbc_nalloc
,
1600 if (rt
->bIntermolecularInteractions
)
1602 /* Check all intermolecular interactions assigned to this atom */
1603 index
= rt
->ril_intermol
.index
;
1604 rtil
= rt
->ril_intermol
.il
;
1606 check_assign_interactions_atom(i
, i_gl
, mol
, i_mol
,
1608 index
[i_gl
], index
[i_gl
+ 1],
1611 bRCheckMB
, rcheck
, bRCheck2B
, rc2
,
1616 idef
, vsite_pbc
, vsite_pbc_nalloc
,
1623 return nbonded_local
;
1626 /*! \brief Set the exclusion data for i-zone \p iz for the case of no exclusions */
1627 static void set_no_exclusions_zone(gmx_domdec_t
*dd
, gmx_domdec_zones_t
*zones
,
1628 int iz
, t_blocka
*lexcls
)
1632 a0
= dd
->cgindex
[zones
->cg_range
[iz
]];
1633 a1
= dd
->cgindex
[zones
->cg_range
[iz
+1]];
1635 for (a
= a0
+1; a
< a1
+1; a
++)
1637 lexcls
->index
[a
] = lexcls
->nra
;
1641 /*! \brief Set the exclusion data for i-zone \p iz
1643 * This is a legacy version for the group scheme of the same routine below.
1644 * Here charge groups and distance checks to ensure unique exclusions
1647 static int make_exclusions_zone_cg(gmx_domdec_t
*dd
, gmx_domdec_zones_t
*zones
,
1648 const gmx_moltype_t
*moltype
,
1649 gmx_bool bRCheck
, real rc2
,
1650 int *la2lc
, t_pbc
*pbc_null
, rvec
*cg_cm
,
1654 int cg_start
, int cg_end
)
1656 int n_excl_at_max
, n
, count
, jla0
, jla1
, jla
;
1657 int cg
, la0
, la1
, la
, a_gl
, mb
, mt
, mol
, a_mol
, j
, aj_mol
;
1658 const t_blocka
*excls
;
1664 jla0
= dd
->cgindex
[zones
->izone
[iz
].jcg0
];
1665 jla1
= dd
->cgindex
[zones
->izone
[iz
].jcg1
];
1667 n_excl_at_max
= dd
->reverse_top
->n_excl_at_max
;
1669 /* We set the end index, but note that we might not start at zero here */
1670 lexcls
->nr
= dd
->cgindex
[cg_end
];
1674 for (cg
= cg_start
; cg
< cg_end
; cg
++)
1676 if (n
+ (cg_end
- cg_start
)*n_excl_at_max
> lexcls
->nalloc_a
)
1678 lexcls
->nalloc_a
= over_alloc_large(n
+ (cg_end
- cg_start
)*n_excl_at_max
);
1679 srenew(lexcls
->a
, lexcls
->nalloc_a
);
1681 la0
= dd
->cgindex
[cg
];
1682 la1
= dd
->cgindex
[cg
+1];
1683 if (GET_CGINFO_EXCL_INTER(cginfo
[cg
]) ||
1684 !GET_CGINFO_EXCL_INTRA(cginfo
[cg
]))
1686 /* Copy the exclusions from the global top */
1687 for (la
= la0
; la
< la1
; la
++)
1689 lexcls
->index
[la
] = n
;
1690 a_gl
= dd
->gatindex
[la
];
1691 global_atomnr_to_moltype_ind(dd
->reverse_top
, a_gl
, &mb
, &mt
, &mol
, &a_mol
);
1692 excls
= &moltype
[mt
].excls
;
1693 for (j
= excls
->index
[a_mol
]; j
< excls
->index
[a_mol
+1]; j
++)
1695 aj_mol
= excls
->a
[j
];
1696 /* This computation of jla is only correct intra-cg */
1697 jla
= la
+ aj_mol
- a_mol
;
1698 if (jla
>= la0
&& jla
< la1
)
1700 /* This is an intra-cg exclusion. We can skip
1701 * the global indexing and distance checking.
1703 /* Intra-cg exclusions are only required
1704 * for the home zone.
1708 lexcls
->a
[n
++] = jla
;
1709 /* Check to avoid double counts */
1718 /* This is a inter-cg exclusion */
1719 /* Since exclusions are pair interactions,
1720 * just like non-bonded interactions,
1721 * they can be assigned properly up
1722 * to the DD cutoff (not cutoff_min as
1723 * for the other bonded interactions).
1725 if (ga2la_get(ga2la
, a_gl
+aj_mol
-a_mol
, &jla
, &cell
))
1727 if (iz
== 0 && cell
== 0)
1729 lexcls
->a
[n
++] = jla
;
1730 /* Check to avoid double counts */
1736 else if (jla
>= jla0
&& jla
< jla1
&&
1738 dd_dist2(pbc_null
, cg_cm
, la2lc
, la
, jla
) < rc2
))
1740 /* jla > la, since jla0 > la */
1741 lexcls
->a
[n
++] = jla
;
1751 /* There are no inter-cg excls and this cg is self-excluded.
1752 * These exclusions are only required for zone 0,
1753 * since other zones do not see themselves.
1757 for (la
= la0
; la
< la1
; la
++)
1759 lexcls
->index
[la
] = n
;
1760 for (j
= la0
; j
< la1
; j
++)
1765 count
+= ((la1
- la0
)*(la1
- la0
- 1))/2;
1769 /* We don't need exclusions for this cg */
1770 for (la
= la0
; la
< la1
; la
++)
1772 lexcls
->index
[la
] = n
;
1778 lexcls
->index
[lexcls
->nr
] = n
;
1784 /*! \brief Set the exclusion data for i-zone \p iz */
1785 static void make_exclusions_zone(gmx_domdec_t
*dd
,
1786 gmx_domdec_zones_t
*zones
,
1787 const gmx_moltype_t
*moltype
,
1791 int at_start
, int at_end
)
1795 int n_excl_at_max
, n
, at
;
1799 jla0
= dd
->cgindex
[zones
->izone
[iz
].jcg0
];
1800 jla1
= dd
->cgindex
[zones
->izone
[iz
].jcg1
];
1802 n_excl_at_max
= dd
->reverse_top
->n_excl_at_max
;
1804 /* We set the end index, but note that we might not start at zero here */
1805 lexcls
->nr
= at_end
;
1808 for (at
= at_start
; at
< at_end
; at
++)
1810 if (n
+ 1000 > lexcls
->nalloc_a
)
1812 lexcls
->nalloc_a
= over_alloc_large(n
+ 1000);
1813 srenew(lexcls
->a
, lexcls
->nalloc_a
);
1815 if (GET_CGINFO_EXCL_INTER(cginfo
[at
]))
1817 int a_gl
, mb
, mt
, mol
, a_mol
, j
;
1818 const t_blocka
*excls
;
1820 if (n
+ n_excl_at_max
> lexcls
->nalloc_a
)
1822 lexcls
->nalloc_a
= over_alloc_large(n
+ n_excl_at_max
);
1823 srenew(lexcls
->a
, lexcls
->nalloc_a
);
1826 /* Copy the exclusions from the global top */
1827 lexcls
->index
[at
] = n
;
1828 a_gl
= dd
->gatindex
[at
];
1829 global_atomnr_to_moltype_ind(dd
->reverse_top
, a_gl
,
1830 &mb
, &mt
, &mol
, &a_mol
);
1831 excls
= &moltype
[mt
].excls
;
1832 for (j
= excls
->index
[a_mol
]; j
< excls
->index
[a_mol
+ 1]; j
++)
1834 int aj_mol
, at_j
, cell
;
1836 aj_mol
= excls
->a
[j
];
1838 if (ga2la_get(ga2la
, a_gl
+ aj_mol
- a_mol
, &at_j
, &cell
))
1840 /* This check is not necessary, but it can reduce
1841 * the number of exclusions in the list, which in turn
1842 * can speed up the pair list construction a bit.
1844 if (at_j
>= jla0
&& at_j
< jla1
)
1846 lexcls
->a
[n
++] = at_j
;
1853 /* We don't need exclusions for this atom */
1854 lexcls
->index
[at
] = n
;
1858 lexcls
->index
[lexcls
->nr
] = n
;
1863 /*! \brief Ensure we have enough space in \p ba for \p nindex_max indices */
1864 static void check_alloc_index(t_blocka
*ba
, int nindex_max
)
1866 if (nindex_max
+1 > ba
->nalloc_index
)
1868 ba
->nalloc_index
= over_alloc_dd(nindex_max
+1);
1869 srenew(ba
->index
, ba
->nalloc_index
);
1873 /*! \brief Ensure that we have enough space for exclusion storate in \p lexcls */
1874 static void check_exclusions_alloc(gmx_domdec_t
*dd
, gmx_domdec_zones_t
*zones
,
1880 nr
= dd
->cgindex
[zones
->izone
[zones
->nizone
-1].cg1
];
1882 check_alloc_index(lexcls
, nr
);
1884 for (thread
= 1; thread
< dd
->reverse_top
->nthread
; thread
++)
1886 check_alloc_index(&dd
->reverse_top
->th_work
[thread
].excl
, nr
);
1890 /*! \brief Set the total count indexes for the local exclusions, needed by several functions */
1891 static void finish_local_exclusions(gmx_domdec_t
*dd
, gmx_domdec_zones_t
*zones
,
1896 lexcls
->nr
= dd
->cgindex
[zones
->izone
[zones
->nizone
-1].cg1
];
1898 if (dd
->n_intercg_excl
== 0)
1900 /* There are no exclusions involving non-home charge groups,
1901 * but we need to set the indices for neighborsearching.
1903 la0
= dd
->cgindex
[zones
->izone
[0].cg1
];
1904 for (la
= la0
; la
< lexcls
->nr
; la
++)
1906 lexcls
->index
[la
] = lexcls
->nra
;
1909 /* nr is only used to loop over the exclusions for Ewald and RF,
1910 * so we can set it to the number of home atoms for efficiency.
1912 lexcls
->nr
= dd
->cgindex
[zones
->izone
[0].cg1
];
1916 /*! \brief Clear a t_idef struct */
1917 static void clear_idef(t_idef
*idef
)
1921 /* Clear the counts */
1922 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
1924 idef
->il
[ftype
].nr
= 0;
1928 /*! \brief Generate and store all required local bonded interactions in \p idef and local exclusions in \p lexcls */
1929 static int make_local_bondeds_excls(gmx_domdec_t
*dd
,
1930 gmx_domdec_zones_t
*zones
,
1931 const gmx_mtop_t
*mtop
,
1933 gmx_bool bRCheckMB
, ivec rcheck
, gmx_bool bRCheck2B
,
1935 int *la2lc
, t_pbc
*pbc_null
, rvec
*cg_cm
,
1936 t_idef
*idef
, gmx_vsite_t
*vsite
,
1937 t_blocka
*lexcls
, int *excl_count
)
1939 int nzone_bondeds
, nzone_excl
;
1940 int izone
, cg0
, cg1
;
1944 gmx_reverse_top_t
*rt
;
1946 if (dd
->reverse_top
->bInterCGInteractions
)
1948 nzone_bondeds
= zones
->n
;
1952 /* Only single charge group (or atom) molecules, so interactions don't
1953 * cross zone boundaries and we only need to assign in the home zone.
1958 if (dd
->n_intercg_excl
> 0)
1960 /* We only use exclusions from i-zones to i- and j-zones */
1961 nzone_excl
= zones
->nizone
;
1965 /* There are no inter-cg exclusions and only zone 0 sees itself */
1969 check_exclusions_alloc(dd
, zones
, lexcls
);
1971 rt
= dd
->reverse_top
;
1975 /* Clear the counts */
1983 for (izone
= 0; izone
< nzone_bondeds
; izone
++)
1985 cg0
= zones
->cg_range
[izone
];
1986 cg1
= zones
->cg_range
[izone
+ 1];
1988 #pragma omp parallel for num_threads(rt->nthread) schedule(static)
1989 for (thread
= 0; thread
< rt
->nthread
; thread
++)
1996 int *vsite_pbc_nalloc
;
1999 cg0t
= cg0
+ ((cg1
- cg0
)* thread
)/rt
->nthread
;
2000 cg1t
= cg0
+ ((cg1
- cg0
)*(thread
+1))/rt
->nthread
;
2008 idef_t
= &rt
->th_work
[thread
].idef
;
2012 if (vsite
&& vsite
->bHaveChargeGroups
&& vsite
->n_intercg_vsite
> 0)
2016 vsite_pbc
= vsite
->vsite_pbc_loc
;
2017 vsite_pbc_nalloc
= vsite
->vsite_pbc_loc_nalloc
;
2021 vsite_pbc
= rt
->th_work
[thread
].vsite_pbc
;
2022 vsite_pbc_nalloc
= rt
->th_work
[thread
].vsite_pbc_nalloc
;
2028 vsite_pbc_nalloc
= NULL
;
2031 rt
->th_work
[thread
].nbonded
=
2032 make_bondeds_zone(dd
, zones
,
2034 bRCheckMB
, rcheck
, bRCheck2B
, rc2
,
2035 la2lc
, pbc_null
, cg_cm
, idef
->iparams
,
2037 vsite_pbc
, vsite_pbc_nalloc
,
2039 dd
->cgindex
[cg0t
], dd
->cgindex
[cg1t
]);
2041 if (izone
< nzone_excl
)
2049 excl_t
= &rt
->th_work
[thread
].excl
;
2054 if (dd
->cgindex
[dd
->ncg_tot
] == dd
->ncg_tot
&&
2057 /* No charge groups and no distance check required */
2058 make_exclusions_zone(dd
, zones
,
2059 mtop
->moltype
, cginfo
,
2066 rt
->th_work
[thread
].excl_count
=
2067 make_exclusions_zone_cg(dd
, zones
,
2068 mtop
->moltype
, bRCheck2B
, rc2
,
2069 la2lc
, pbc_null
, cg_cm
, cginfo
,
2076 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
2079 if (rt
->nthread
> 1)
2081 combine_idef(idef
, rt
->th_work
, rt
->nthread
, vsite
);
2084 for (thread
= 0; thread
< rt
->nthread
; thread
++)
2086 nbonded_local
+= rt
->th_work
[thread
].nbonded
;
2089 if (izone
< nzone_excl
)
2091 if (rt
->nthread
> 1)
2093 combine_blocka(lexcls
, rt
->th_work
, rt
->nthread
);
2096 for (thread
= 0; thread
< rt
->nthread
; thread
++)
2098 *excl_count
+= rt
->th_work
[thread
].excl_count
;
2103 /* Some zones might not have exclusions, but some code still needs to
2104 * loop over the index, so we set the indices here.
2106 for (izone
= nzone_excl
; izone
< zones
->nizone
; izone
++)
2108 set_no_exclusions_zone(dd
, zones
, izone
, lexcls
);
2111 finish_local_exclusions(dd
, zones
, lexcls
);
2114 fprintf(debug
, "We have %d exclusions, check count %d\n",
2115 lexcls
->nra
, *excl_count
);
2118 return nbonded_local
;
2121 void dd_make_local_cgs(gmx_domdec_t
*dd
, t_block
*lcgs
)
2123 lcgs
->nr
= dd
->ncg_tot
;
2124 lcgs
->index
= dd
->cgindex
;
2127 void dd_make_local_top(gmx_domdec_t
*dd
, gmx_domdec_zones_t
*zones
,
2128 int npbcdim
, matrix box
,
2129 rvec cellsize_min
, ivec npulse
,
2133 gmx_mtop_t
*mtop
, gmx_localtop_t
*ltop
)
2135 gmx_bool bRCheckMB
, bRCheck2B
;
2139 t_pbc pbc
, *pbc_null
= NULL
;
2143 fprintf(debug
, "Making local topology\n");
2146 dd_make_local_cgs(dd
, <op
->cgs
);
2151 if (dd
->reverse_top
->bInterCGInteractions
)
2153 /* We need to check to which cell bondeds should be assigned */
2154 rc
= dd_cutoff_twobody(dd
);
2157 fprintf(debug
, "Two-body bonded cut-off distance is %g\n", rc
);
2160 /* Should we check cg_cm distances when assigning bonded interactions? */
2161 for (d
= 0; d
< DIM
; d
++)
2164 /* Only need to check for dimensions where the part of the box
2165 * that is not communicated is smaller than the cut-off.
2167 if (d
< npbcdim
&& dd
->nc
[d
] > 1 &&
2168 (dd
->nc
[d
] - npulse
[d
])*cellsize_min
[d
] < 2*rc
)
2175 /* Check for interactions between two atoms,
2176 * where we can allow interactions up to the cut-off,
2177 * instead of up to the smallest cell dimension.
2184 "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %d\n",
2185 d
, cellsize_min
[d
], d
, rcheck
[d
], bRCheck2B
);
2188 if (bRCheckMB
|| bRCheck2B
)
2193 set_pbc_dd(&pbc
, fr
->ePBC
, dd
, TRUE
, box
);
2204 make_local_bondeds_excls(dd
, zones
, mtop
, fr
->cginfo
,
2205 bRCheckMB
, rcheck
, bRCheck2B
, rc
,
2207 pbc_null
, cgcm_or_x
,
2209 <op
->excls
, &nexcl
);
2211 /* The ilist is not sorted yet,
2212 * we can only do this when we have the charge arrays.
2214 ltop
->idef
.ilsort
= ilsortUNKNOWN
;
2216 if (dd
->reverse_top
->bExclRequired
)
2218 dd
->nbonded_local
+= nexcl
;
2220 forcerec_set_excl_load(fr
, ltop
);
2223 ltop
->atomtypes
= mtop
->atomtypes
;
2225 /* For an error message only */
2226 dd
->reverse_top
->err_top_global
= mtop
;
2227 dd
->reverse_top
->err_top_local
= ltop
;
2230 void dd_sort_local_top(gmx_domdec_t
*dd
, t_mdatoms
*mdatoms
,
2231 gmx_localtop_t
*ltop
)
2233 if (dd
->reverse_top
->ilsort
== ilsortNO_FE
)
2235 ltop
->idef
.ilsort
= ilsortNO_FE
;
2239 gmx_sort_ilist_fe(<op
->idef
, mdatoms
->chargeA
, mdatoms
->chargeB
);
2243 gmx_localtop_t
*dd_init_local_top(gmx_mtop_t
*top_global
)
2245 gmx_localtop_t
*top
;
2250 top
->idef
.ntypes
= top_global
->ffparams
.ntypes
;
2251 top
->idef
.atnr
= top_global
->ffparams
.atnr
;
2252 top
->idef
.functype
= top_global
->ffparams
.functype
;
2253 top
->idef
.iparams
= top_global
->ffparams
.iparams
;
2254 top
->idef
.fudgeQQ
= top_global
->ffparams
.fudgeQQ
;
2255 top
->idef
.cmap_grid
= top_global
->ffparams
.cmap_grid
;
2257 for (i
= 0; i
< F_NRE
; i
++)
2259 top
->idef
.il
[i
].iatoms
= NULL
;
2260 top
->idef
.il
[i
].nalloc
= 0;
2262 top
->idef
.ilsort
= ilsortUNKNOWN
;
2267 void dd_init_local_state(gmx_domdec_t
*dd
,
2268 t_state
*state_global
, t_state
*state_local
)
2270 int buf
[NITEM_DD_INIT_LOCAL_STATE
];
2274 buf
[0] = state_global
->flags
;
2275 buf
[1] = state_global
->ngtc
;
2276 buf
[2] = state_global
->nnhpres
;
2277 buf
[3] = state_global
->nhchainlength
;
2278 buf
[4] = state_global
->dfhist
.nlambda
;
2280 dd_bcast(dd
, NITEM_DD_INIT_LOCAL_STATE
*sizeof(int), buf
);
2282 init_state(state_local
, 0, buf
[1], buf
[2], buf
[3], buf
[4]);
2283 state_local
->flags
= buf
[0];
2286 /*! \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 */
2287 static void check_link(t_blocka
*link
, int cg_gl
, int cg_gl_j
)
2293 for (k
= link
->index
[cg_gl
]; k
< link
->index
[cg_gl
+1]; k
++)
2295 GMX_RELEASE_ASSERT(link
->a
, "Inconsistent NULL pointer while making charge-group links");
2296 if (link
->a
[k
] == cg_gl_j
)
2303 GMX_RELEASE_ASSERT(link
->a
|| link
->index
[cg_gl
+1]+1 > link
->nalloc_a
,
2304 "Inconsistent allocation of link");
2305 /* Add this charge group link */
2306 if (link
->index
[cg_gl
+1]+1 > link
->nalloc_a
)
2308 link
->nalloc_a
= over_alloc_large(link
->index
[cg_gl
+1]+1);
2309 srenew(link
->a
, link
->nalloc_a
);
2311 link
->a
[link
->index
[cg_gl
+1]] = cg_gl_j
;
2312 link
->index
[cg_gl
+1]++;
2316 /*! \brief Allocate and return an array of the charge group index for all atoms */
2317 static int *make_at2cg(t_block
*cgs
)
2321 snew(at2cg
, cgs
->index
[cgs
->nr
]);
2322 for (cg
= 0; cg
< cgs
->nr
; cg
++)
2324 for (a
= cgs
->index
[cg
]; a
< cgs
->index
[cg
+1]; a
++)
2333 t_blocka
*make_charge_group_links(gmx_mtop_t
*mtop
, gmx_domdec_t
*dd
,
2334 cginfo_mb_t
*cginfo_mb
)
2336 gmx_bool bExclRequired
;
2337 int mb
, cg_offset
, cg
, cg_gl
, a
, aj
, i
, j
, ftype
, nral
, nlink_mol
, mol
, ncgi
;
2338 gmx_molblock_t
*molb
;
2339 gmx_moltype_t
*molt
;
2343 reverse_ilist_t ril
, ril_intermol
;
2345 cginfo_mb_t
*cgi_mb
;
2347 /* For each charge group make a list of other charge groups
2348 * in the system that a linked to it via bonded interactions
2349 * which are also stored in reverse_top.
2352 bExclRequired
= dd
->reverse_top
->bExclRequired
;
2354 if (mtop
->bIntermolecularInteractions
)
2356 if (ncg_mtop(mtop
) < mtop
->natoms
)
2358 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.");
2363 atoms
.nr
= mtop
->natoms
;
2366 make_reverse_ilist(mtop
->intermolecular_ilist
, &atoms
,
2367 NULL
, FALSE
, FALSE
, FALSE
, TRUE
, &ril_intermol
);
2371 snew(link
->index
, ncg_mtop(mtop
)+1);
2378 for (mb
= 0; mb
< mtop
->nmolblock
; mb
++)
2380 molb
= &mtop
->molblock
[mb
];
2381 if (molb
->nmol
== 0)
2385 molt
= &mtop
->moltype
[molb
->type
];
2387 excls
= &molt
->excls
;
2388 a2c
= make_at2cg(cgs
);
2389 /* Make a reverse ilist in which the interactions are linked
2390 * to all atoms, not only the first atom as in gmx_reverse_top.
2391 * The constraints are discarded here.
2393 make_reverse_ilist(molt
->ilist
, &molt
->atoms
,
2394 NULL
, FALSE
, FALSE
, FALSE
, TRUE
, &ril
);
2396 cgi_mb
= &cginfo_mb
[mb
];
2398 for (mol
= 0; mol
< (mtop
->bIntermolecularInteractions
? molb
->nmol
: 1); mol
++)
2400 for (cg
= 0; cg
< cgs
->nr
; cg
++)
2402 cg_gl
= cg_offset
+ cg
;
2403 link
->index
[cg_gl
+1] = link
->index
[cg_gl
];
2404 for (a
= cgs
->index
[cg
]; a
< cgs
->index
[cg
+1]; a
++)
2407 while (i
< ril
.index
[a
+1])
2409 ftype
= ril
.il
[i
++];
2411 /* Skip the ifunc index */
2413 for (j
= 0; j
< nral
; j
++)
2418 check_link(link
, cg_gl
, cg_offset
+a2c
[aj
]);
2421 i
+= nral_rt(ftype
);
2425 /* Exclusions always go both ways */
2426 for (j
= excls
->index
[a
]; j
< excls
->index
[a
+1]; j
++)
2431 check_link(link
, cg_gl
, cg_offset
+a2c
[aj
]);
2436 if (mtop
->bIntermolecularInteractions
)
2438 i
= ril_intermol
.index
[a
];
2439 while (i
< ril_intermol
.index
[a
+1])
2441 ftype
= ril_intermol
.il
[i
++];
2443 /* Skip the ifunc index */
2445 for (j
= 0; j
< nral
; j
++)
2447 /* Here we assume we have no charge groups;
2448 * this has been checked above.
2450 aj
= ril_intermol
.il
[i
+j
];
2451 check_link(link
, cg_gl
, aj
);
2453 i
+= nral_rt(ftype
);
2457 if (link
->index
[cg_gl
+1] - link
->index
[cg_gl
] > 0)
2459 SET_CGINFO_BOND_INTER(cgi_mb
->cginfo
[cg
]);
2464 cg_offset
+= cgs
->nr
;
2466 nlink_mol
= link
->index
[cg_offset
] - link
->index
[cg_offset
-cgs
->nr
];
2468 destroy_reverse_ilist(&ril
);
2473 fprintf(debug
, "molecule type '%s' %d cgs has %d cg links through bonded interac.\n", *molt
->name
, cgs
->nr
, nlink_mol
);
2476 if (molb
->nmol
> mol
)
2478 /* Copy the data for the rest of the molecules in this block */
2479 link
->nalloc_a
+= (molb
->nmol
- mol
)*nlink_mol
;
2480 srenew(link
->a
, link
->nalloc_a
);
2481 for (; mol
< molb
->nmol
; mol
++)
2483 for (cg
= 0; cg
< cgs
->nr
; cg
++)
2485 cg_gl
= cg_offset
+ cg
;
2486 link
->index
[cg_gl
+1] =
2487 link
->index
[cg_gl
+1-cgs
->nr
] + nlink_mol
;
2488 for (j
= link
->index
[cg_gl
]; j
< link
->index
[cg_gl
+1]; j
++)
2490 link
->a
[j
] = link
->a
[j
-nlink_mol
] + cgs
->nr
;
2492 if (link
->index
[cg_gl
+1] - link
->index
[cg_gl
] > 0 &&
2493 cg_gl
- cgi_mb
->cg_start
< cgi_mb
->cg_mod
)
2495 SET_CGINFO_BOND_INTER(cgi_mb
->cginfo
[cg_gl
- cgi_mb
->cg_start
]);
2499 cg_offset
+= cgs
->nr
;
2504 if (mtop
->bIntermolecularInteractions
)
2506 destroy_reverse_ilist(&ril_intermol
);
2511 fprintf(debug
, "Of the %d charge groups %d are linked via bonded interactions\n", ncg_mtop(mtop
), ncgi
);
2522 } bonded_distance_t
;
2524 /*! \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 */
2525 static void update_max_bonded_distance(real r2
, int ftype
, int a1
, int a2
,
2526 bonded_distance_t
*bd
)
2537 /*! \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 */
2538 static void bonded_cg_distance_mol(gmx_moltype_t
*molt
, int *at2cg
,
2539 gmx_bool bBCheck
, gmx_bool bExcl
, rvec
*cg_cm
,
2540 bonded_distance_t
*bd_2b
,
2541 bonded_distance_t
*bd_mb
)
2543 for (int ftype
= 0; ftype
< F_NRE
; ftype
++)
2545 if (dd_check_ftype(ftype
, bBCheck
, FALSE
, FALSE
))
2547 const t_ilist
*il
= &molt
->ilist
[ftype
];
2548 int nral
= NRAL(ftype
);
2551 for (int i
= 0; i
< il
->nr
; i
+= 1+nral
)
2553 for (int ai
= 0; ai
< nral
; ai
++)
2555 int cgi
= at2cg
[il
->iatoms
[i
+1+ai
]];
2556 for (int aj
= ai
+ 1; aj
< nral
; aj
++)
2558 int cgj
= at2cg
[il
->iatoms
[i
+1+aj
]];
2561 real rij2
= distance2(cg_cm
[cgi
], cg_cm
[cgj
]);
2563 update_max_bonded_distance(rij2
, ftype
,
2566 (nral
== 2) ? bd_2b
: bd_mb
);
2576 const t_blocka
*excls
= &molt
->excls
;
2577 for (int ai
= 0; ai
< excls
->nr
; ai
++)
2579 int cgi
= at2cg
[ai
];
2580 for (int j
= excls
->index
[ai
]; j
< excls
->index
[ai
+1]; j
++)
2582 int cgj
= at2cg
[excls
->a
[j
]];
2585 real rij2
= distance2(cg_cm
[cgi
], cg_cm
[cgj
]);
2587 /* There is no function type for exclusions, use -1 */
2588 update_max_bonded_distance(rij2
, -1, ai
, excls
->a
[j
], bd_2b
);
2595 /*! \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 */
2596 static void bonded_distance_intermol(const t_ilist
*ilists_intermol
,
2598 rvec
*x
, int ePBC
, matrix box
,
2599 bonded_distance_t
*bd_2b
,
2600 bonded_distance_t
*bd_mb
)
2604 set_pbc(&pbc
, ePBC
, box
);
2606 for (int ftype
= 0; ftype
< F_NRE
; ftype
++)
2608 if (dd_check_ftype(ftype
, bBCheck
, FALSE
, FALSE
))
2610 const t_ilist
*il
= &ilists_intermol
[ftype
];
2611 int nral
= NRAL(ftype
);
2613 /* No nral>1 check here, since intermol interactions always
2614 * have nral>=2 (and the code is also correct for nral=1).
2616 for (int i
= 0; i
< il
->nr
; i
+= 1+nral
)
2618 for (int ai
= 0; ai
< nral
; ai
++)
2620 int atom_i
= il
->iatoms
[i
+ 1 + ai
];
2622 for (int aj
= ai
+ 1; aj
< nral
; aj
++)
2627 int atom_j
= il
->iatoms
[i
+ 1 + aj
];
2629 pbc_dx(&pbc
, x
[atom_i
], x
[atom_j
], dx
);
2633 update_max_bonded_distance(rij2
, ftype
,
2635 (nral
== 2) ? bd_2b
: bd_mb
);
2643 //! Compute charge group centers of mass for molecule \p molt
2644 static void get_cgcm_mol(gmx_moltype_t
*molt
, gmx_ffparams_t
*ffparams
,
2645 int ePBC
, t_graph
*graph
, matrix box
,
2647 rvec
*x
, rvec
*xs
, rvec
*cg_cm
)
2651 if (ePBC
!= epbcNONE
)
2653 mk_mshift(NULL
, graph
, ePBC
, box
, x
);
2655 shift_x(graph
, box
, x
, xs
);
2656 /* By doing an extra mk_mshift the molecules that are broken
2657 * because they were e.g. imported from another software
2658 * will be made whole again. Such are the healing powers
2661 mk_mshift(NULL
, graph
, ePBC
, box
, xs
);
2665 /* We copy the coordinates so the original coordinates remain
2666 * unchanged, just to be 100% sure that we do not affect
2667 * binary reproducibility of simulations.
2669 n
= molt
->cgs
.index
[molt
->cgs
.nr
];
2670 for (i
= 0; i
< n
; i
++)
2672 copy_rvec(x
[i
], xs
[i
]);
2678 construct_vsites(vsite
, xs
, 0.0, NULL
,
2679 ffparams
->iparams
, molt
->ilist
,
2680 epbcNONE
, TRUE
, NULL
, NULL
);
2683 calc_cgcm(NULL
, 0, molt
->cgs
.nr
, &molt
->cgs
, xs
, cg_cm
);
2686 //! Returns whether \p molt has a virtual site
2687 static int have_vsite_molt(gmx_moltype_t
*molt
)
2693 for (i
= 0; i
< F_NRE
; i
++)
2695 if ((interaction_function
[i
].flags
& IF_VSITE
) &&
2696 molt
->ilist
[i
].nr
> 0)
2705 void dd_bonded_cg_distance(FILE *fplog
,
2707 t_inputrec
*ir
, rvec
*x
, matrix box
,
2709 real
*r_2b
, real
*r_mb
)
2711 gmx_bool bExclRequired
;
2712 int mb
, at_offset
, *at2cg
, mol
;
2715 gmx_molblock_t
*molb
;
2716 gmx_moltype_t
*molt
;
2718 bonded_distance_t bd_2b
= { 0, -1, -1, -1 };
2719 bonded_distance_t bd_mb
= { 0, -1, -1, -1 };
2721 bExclRequired
= inputrecExclForces(ir
);
2723 vsite
= init_vsite(mtop
, NULL
, TRUE
);
2728 for (mb
= 0; mb
< mtop
->nmolblock
; mb
++)
2730 molb
= &mtop
->molblock
[mb
];
2731 molt
= &mtop
->moltype
[molb
->type
];
2732 if (molt
->cgs
.nr
== 1 || molb
->nmol
== 0)
2734 at_offset
+= molb
->nmol
*molt
->atoms
.nr
;
2738 if (ir
->ePBC
!= epbcNONE
)
2740 mk_graph_ilist(NULL
, molt
->ilist
, 0, molt
->atoms
.nr
, FALSE
, FALSE
,
2744 at2cg
= make_at2cg(&molt
->cgs
);
2745 snew(xs
, molt
->atoms
.nr
);
2746 snew(cg_cm
, molt
->cgs
.nr
);
2747 for (mol
= 0; mol
< molb
->nmol
; mol
++)
2749 get_cgcm_mol(molt
, &mtop
->ffparams
, ir
->ePBC
, &graph
, box
,
2750 have_vsite_molt(molt
) ? vsite
: NULL
,
2751 x
+at_offset
, xs
, cg_cm
);
2753 bonded_distance_t bd_mol_2b
= { 0, -1, -1, -1 };
2754 bonded_distance_t bd_mol_mb
= { 0, -1, -1, -1 };
2756 bonded_cg_distance_mol(molt
, at2cg
, bBCheck
, bExclRequired
, cg_cm
,
2757 &bd_mol_2b
, &bd_mol_mb
);
2759 /* Process the mol data adding the atom index offset */
2760 update_max_bonded_distance(bd_mol_2b
.r2
, bd_mol_2b
.ftype
,
2761 at_offset
+ bd_mol_2b
.a1
,
2762 at_offset
+ bd_mol_2b
.a2
,
2764 update_max_bonded_distance(bd_mol_mb
.r2
, bd_mol_mb
.ftype
,
2765 at_offset
+ bd_mol_mb
.a1
,
2766 at_offset
+ bd_mol_mb
.a2
,
2769 at_offset
+= molt
->atoms
.nr
;
2774 if (ir
->ePBC
!= epbcNONE
)
2781 /* We should have a vsite free routine, but here we can simply free */
2784 if (mtop
->bIntermolecularInteractions
)
2786 if (ncg_mtop(mtop
) < mtop
->natoms
)
2788 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.");
2791 bonded_distance_intermol(mtop
->intermolecular_ilist
,
2797 *r_2b
= sqrt(bd_2b
.r2
);
2798 *r_mb
= sqrt(bd_mb
.r2
);
2800 if (fplog
&& (*r_2b
> 0 || *r_mb
> 0))
2803 "Initial maximum inter charge-group distances:\n");
2807 " two-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2808 *r_2b
, (bd_2b
.ftype
>= 0) ? interaction_function
[bd_2b
.ftype
].longname
: "Exclusion",
2809 bd_2b
.a1
+ 1, bd_2b
.a2
+ 1);
2814 " multi-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2815 *r_mb
, interaction_function
[bd_mb
.ftype
].longname
,
2816 bd_mb
.a1
+ 1, bd_mb
.a2
+ 1);