2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006,2007,2008,2009,2010,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 implements functions for domdec to use
39 * while managing inter-atomic constraints.
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_domdec
47 #include "domdec_constraints.h"
53 #include "gromacs/domdec/dlbtiming.h"
54 #include "gromacs/domdec/domdec.h"
55 #include "gromacs/domdec/domdec_struct.h"
56 #include "gromacs/domdec/ga2la.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/mdlib/constr.h"
59 #include "gromacs/mdlib/gmx_omp_nthreads.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/mdtypes/forcerec.h" // only for GET_CGINFO_*
62 #include "gromacs/pbcutil/ishift.h"
63 #include "gromacs/topology/ifunc.h"
64 #include "gromacs/topology/mtop_lookup.h"
65 #include "gromacs/utility/exceptions.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "gromacs/utility/gmxassert.h"
68 #include "gromacs/utility/smalloc.h"
70 #include "domdec_specatomcomm.h"
71 #include "domdec_vsite.h"
74 /*! \brief Struct used during constraint setup with domain decomposition */
75 struct gmx_domdec_constraints_t
{
76 //! @cond Doxygen_Suppress
77 int *molb_con_offset
; /**< Offset in the constraint array for each molblock */
78 int *molb_ncon_mol
; /**< The number of constraints per molecule for each molblock */
80 int ncon
; /**< The fully local and conneced constraints */
81 /* The global constraint number, only required for clearing gc_req */
82 int *con_gl
; /**< Global constraint indices for local constraints */
83 int *con_nlocat
; /**< Number of local atoms (2/1/0) for each constraint */
84 int con_nalloc
; /**< Allocation size for \p con_gl and \p con_nlocat */
86 char *gc_req
; /**< Boolean that tells if a global constraint index has been requested; note: size global #constraints */
87 gmx_hash_t
*ga2la
; /**< Global to local communicated constraint atom only index */
89 /* Multi-threading stuff */
90 int nthread
; /**< Number of threads used for DD constraint setup */
91 t_ilist
*ils
; /**< Constraint ilist working arrays, size \p nthread */
95 void dd_move_x_constraints(gmx_domdec_t
*dd
, matrix box
,
96 rvec
*x0
, rvec
*x1
, gmx_bool bX1IsCoord
)
98 if (dd
->constraint_comm
)
100 dd_move_x_specat(dd
, dd
->constraint_comm
, box
, x0
, x1
, bX1IsCoord
);
102 ddReopenBalanceRegionCpu(dd
);
106 int *dd_constraints_nlocalatoms(gmx_domdec_t
*dd
)
110 return dd
->constraints
->con_nlocat
;
118 void dd_clear_local_constraint_indices(gmx_domdec_t
*dd
)
120 gmx_domdec_constraints_t
*dc
;
123 dc
= dd
->constraints
;
125 for (i
= 0; i
< dc
->ncon
; i
++)
127 dc
->gc_req
[dc
->con_gl
[i
]] = 0;
130 if (dd
->constraint_comm
)
132 gmx_hash_clear_and_optimize(dc
->ga2la
);
136 void dd_clear_local_vsite_indices(gmx_domdec_t
*dd
)
140 gmx_hash_clear_and_optimize(dd
->ga2la_vsite
);
144 /*! \brief Walks over the constraints out from the local atoms into the non-local atoms and adds them to a list */
145 static void walk_out(int con
, int con_offset
, int a
, int offset
, int nrec
,
146 int ncon1
, const t_iatom
*ia1
, const t_iatom
*ia2
,
147 const t_blocka
*at2con
,
148 const gmx_ga2la_t
*ga2la
, gmx_bool bHomeConnect
,
149 gmx_domdec_constraints_t
*dc
,
150 gmx_domdec_specat_comm_t
*dcc
,
154 int a1_gl
, a2_gl
, a_loc
, i
, coni
, b
;
157 if (dc
->gc_req
[con_offset
+con
] == 0)
159 /* Add this non-home constraint to the list */
160 if (dc
->ncon
+1 > dc
->con_nalloc
)
162 dc
->con_nalloc
= over_alloc_large(dc
->ncon
+1);
163 srenew(dc
->con_gl
, dc
->con_nalloc
);
164 srenew(dc
->con_nlocat
, dc
->con_nalloc
);
166 dc
->con_gl
[dc
->ncon
] = con_offset
+ con
;
167 dc
->con_nlocat
[dc
->ncon
] = (bHomeConnect
? 1 : 0);
168 dc
->gc_req
[con_offset
+con
] = 1;
169 if (il_local
->nr
+ 3 > il_local
->nalloc
)
171 il_local
->nalloc
= over_alloc_dd(il_local
->nr
+3);
172 srenew(il_local
->iatoms
, il_local
->nalloc
);
174 iap
= constr_iatomptr(ncon1
, ia1
, ia2
, con
);
175 il_local
->iatoms
[il_local
->nr
++] = iap
[0];
176 a1_gl
= offset
+ iap
[1];
177 a2_gl
= offset
+ iap
[2];
178 /* The following indexing code can probably be optizimed */
179 if (ga2la_get_home(ga2la
, a1_gl
, &a_loc
))
181 il_local
->iatoms
[il_local
->nr
++] = a_loc
;
185 /* We set this index later */
186 il_local
->iatoms
[il_local
->nr
++] = -a1_gl
- 1;
188 if (ga2la_get_home(ga2la
, a2_gl
, &a_loc
))
190 il_local
->iatoms
[il_local
->nr
++] = a_loc
;
194 /* We set this index later */
195 il_local
->iatoms
[il_local
->nr
++] = -a2_gl
- 1;
199 /* Check to not ask for the same atom more than once */
200 if (gmx_hash_get_minone(dc
->ga2la
, offset
+a
) == -1)
203 /* Add this non-home atom to the list */
204 if (ireq
->n
+1 > ireq
->nalloc
)
206 ireq
->nalloc
= over_alloc_large(ireq
->n
+1);
207 srenew(ireq
->ind
, ireq
->nalloc
);
209 ireq
->ind
[ireq
->n
++] = offset
+ a
;
210 /* Temporarily mark with -2, we get the index later */
211 gmx_hash_set(dc
->ga2la
, offset
+a
, -2);
216 for (i
= at2con
->index
[a
]; i
< at2con
->index
[a
+1]; i
++)
222 iap
= constr_iatomptr(ncon1
, ia1
, ia2
, coni
);
231 if (!ga2la_get_home(ga2la
, offset
+b
, &a_loc
))
233 walk_out(coni
, con_offset
, b
, offset
, nrec
-1,
234 ncon1
, ia1
, ia2
, at2con
,
235 ga2la
, FALSE
, dc
, dcc
, il_local
, ireq
);
242 /*! \brief Looks up SETTLE constraints for a range of charge-groups */
243 static void atoms_to_settles(gmx_domdec_t
*dd
,
244 const gmx_mtop_t
*mtop
,
246 const int **at2settle_mt
,
247 int cg_start
, int cg_end
,
251 gmx_ga2la_t
*ga2la
= dd
->ga2la
;
252 int nral
= NRAL(F_SETTLE
);
255 for (int cg
= cg_start
; cg
< cg_end
; cg
++)
257 if (GET_CGINFO_SETTLE(cginfo
[cg
]))
259 for (int a
= dd
->cgindex
[cg
]; a
< dd
->cgindex
[cg
+1]; a
++)
261 int a_gl
= dd
->gatindex
[a
];
263 mtopGetMolblockIndex(mtop
, a_gl
, &mb
, nullptr, &a_mol
);
265 const gmx_molblock_t
*molb
= &mtop
->molblock
[mb
];
266 int settle
= at2settle_mt
[molb
->type
][a_mol
];
270 int offset
= a_gl
- a_mol
;
272 t_iatom
*ia1
= mtop
->moltype
[molb
->type
].ilist
[F_SETTLE
].iatoms
;
274 int a_gls
[3], a_locs
[3];
275 gmx_bool bAssign
= FALSE
;
277 for (int sa
= 0; sa
< nral
; sa
++)
279 int a_glsa
= offset
+ ia1
[settle
*(1+nral
)+1+sa
];
281 if (ga2la_get_home(ga2la
, a_glsa
, &a_locs
[sa
]))
283 if (nlocal
== 0 && a_gl
== a_glsa
)
293 if (ils_local
->nr
+1+nral
> ils_local
->nalloc
)
295 ils_local
->nalloc
= over_alloc_dd(ils_local
->nr
+1+nral
);
296 srenew(ils_local
->iatoms
, ils_local
->nalloc
);
299 ils_local
->iatoms
[ils_local
->nr
++] = ia1
[settle
*4];
301 for (int sa
= 0; sa
< nral
; sa
++)
303 if (ga2la_get_home(ga2la
, a_gls
[sa
], &a_locs
[sa
]))
305 ils_local
->iatoms
[ils_local
->nr
++] = a_locs
[sa
];
309 ils_local
->iatoms
[ils_local
->nr
++] = -a_gls
[sa
] - 1;
310 /* Add this non-home atom to the list */
311 if (ireq
->n
+1 > ireq
->nalloc
)
313 ireq
->nalloc
= over_alloc_large(ireq
->n
+1);
314 srenew(ireq
->ind
, ireq
->nalloc
);
316 ireq
->ind
[ireq
->n
++] = a_gls
[sa
];
317 /* A check on double atom requests is
318 * not required for settle.
329 /*! \brief Looks up constraint for the local atoms */
330 static void atoms_to_constraints(gmx_domdec_t
*dd
,
331 const gmx_mtop_t
*mtop
,
333 const t_blocka
*at2con_mt
, int nrec
,
337 const t_blocka
*at2con
;
339 t_iatom
*ia1
, *ia2
, *iap
;
340 int a_loc
, b_lo
, offset
, b_mol
, i
, con
, con_offset
;
342 gmx_domdec_constraints_t
*dc
= dd
->constraints
;
343 gmx_domdec_specat_comm_t
*dcc
= dd
->constraint_comm
;
345 gmx_ga2la_t
*ga2la
= dd
->ga2la
;
349 for (int cg
= 0; cg
< dd
->ncg_home
; cg
++)
351 if (GET_CGINFO_CONSTR(cginfo
[cg
]))
353 for (int a
= dd
->cgindex
[cg
]; a
< dd
->cgindex
[cg
+1]; a
++)
355 int a_gl
= dd
->gatindex
[a
];
357 mtopGetMolblockIndex(mtop
, a_gl
, &mb
, &molnr
, &a_mol
);
359 const gmx_molblock_t
*molb
= &mtop
->molblock
[mb
];
361 ncon1
= mtop
->moltype
[molb
->type
].ilist
[F_CONSTR
].nr
/NRAL(F_SETTLE
);
363 ia1
= mtop
->moltype
[molb
->type
].ilist
[F_CONSTR
].iatoms
;
364 ia2
= mtop
->moltype
[molb
->type
].ilist
[F_CONSTRNC
].iatoms
;
366 /* Calculate the global constraint number offset for the molecule.
367 * This is only required for the global index to make sure
368 * that we use each constraint only once.
371 dc
->molb_con_offset
[mb
] + molnr
*dc
->molb_ncon_mol
[mb
];
373 /* The global atom number offset for this molecule */
374 offset
= a_gl
- a_mol
;
375 at2con
= &at2con_mt
[molb
->type
];
376 for (i
= at2con
->index
[a_mol
]; i
< at2con
->index
[a_mol
+1]; i
++)
379 iap
= constr_iatomptr(ncon1
, ia1
, ia2
, con
);
388 if (ga2la_get_home(ga2la
, offset
+b_mol
, &a_loc
))
390 /* Add this fully home constraint at the first atom */
393 if (dc
->ncon
+1 > dc
->con_nalloc
)
395 dc
->con_nalloc
= over_alloc_large(dc
->ncon
+1);
396 srenew(dc
->con_gl
, dc
->con_nalloc
);
397 srenew(dc
->con_nlocat
, dc
->con_nalloc
);
399 dc
->con_gl
[dc
->ncon
] = con_offset
+ con
;
400 dc
->con_nlocat
[dc
->ncon
] = 2;
401 if (ilc_local
->nr
+ 3 > ilc_local
->nalloc
)
403 ilc_local
->nalloc
= over_alloc_dd(ilc_local
->nr
+ 3);
404 srenew(ilc_local
->iatoms
, ilc_local
->nalloc
);
407 ilc_local
->iatoms
[ilc_local
->nr
++] = iap
[0];
408 ilc_local
->iatoms
[ilc_local
->nr
++] = (a_gl
== iap
[1] ? a
: b_lo
);
409 ilc_local
->iatoms
[ilc_local
->nr
++] = (a_gl
== iap
[1] ? b_lo
: a
);
416 /* We need the nrec constraints coupled to this constraint,
417 * so we need to walk out of the home cell by nrec+1 atoms,
418 * since already atom bg is not locally present.
419 * Therefore we call walk_out with nrec recursions to go
420 * after this first call.
422 walk_out(con
, con_offset
, b_mol
, offset
, nrec
,
423 ncon1
, ia1
, ia2
, at2con
,
424 dd
->ga2la
, TRUE
, dc
, dcc
, ilc_local
, ireq
);
434 "Constraints: home %3d border %3d atoms: %3d\n",
435 nhome
, dc
->ncon
-nhome
,
436 dd
->constraint_comm
? ireq
->n
: 0);
440 int dd_make_local_constraints(gmx_domdec_t
*dd
, int at_start
,
441 const struct gmx_mtop_t
*mtop
,
443 gmx_constr_t constr
, int nrec
,
446 gmx_domdec_constraints_t
*dc
;
447 t_ilist
*ilc_local
, *ils_local
;
449 const t_blocka
*at2con_mt
;
450 const int **at2settle_mt
;
451 gmx_hash_t
*ga2la_specat
;
455 // This code should not be called unless this condition is true,
456 // because that's the only time init_domdec_constraints is
458 GMX_RELEASE_ASSERT(dd
->bInterCGcons
|| dd
->bInterCGsettles
, "dd_make_local_constraints called when there are no local constraints");
459 // ... and init_domdec_constraints always sets
460 // dd->constraint_comm...
461 GMX_RELEASE_ASSERT(dd
->constraint_comm
, "Invalid use of dd_make_local_constraints before construction of constraint_comm");
462 // ... which static analysis needs to be reassured about, because
463 // otherwise, when dd->bInterCGsettles is
464 // true. dd->constraint_comm is unilaterally dereferenced before
465 // the call to atoms_to_settles.
467 dc
= dd
->constraints
;
469 ilc_local
= &il_local
[F_CONSTR
];
470 ils_local
= &il_local
[F_SETTLE
];
474 if (dd
->constraint_comm
)
476 at2con_mt
= atom2constraints_moltype(constr
);
477 ireq
= &dd
->constraint_comm
->ireq
[0];
482 // Currently unreachable
487 if (dd
->bInterCGsettles
)
489 at2settle_mt
= atom2settle_moltype(constr
);
494 /* Settle works inside charge groups, we assigned them already */
495 at2settle_mt
= nullptr;
498 if (at2settle_mt
== nullptr)
500 atoms_to_constraints(dd
, mtop
, cginfo
, at2con_mt
, nrec
,
508 /* Do the constraints, if present, on the first thread.
509 * Do the settles on all other threads.
511 t0_set
= ((at2con_mt
!= nullptr && dc
->nthread
> 1) ? 1 : 0);
513 #pragma omp parallel for num_threads(dc->nthread) schedule(static)
514 for (thread
= 0; thread
< dc
->nthread
; thread
++)
518 if (at2con_mt
&& thread
== 0)
520 atoms_to_constraints(dd
, mtop
, cginfo
, at2con_mt
, nrec
,
524 if (thread
>= t0_set
)
530 /* Distribute the settle check+assignments over
531 * dc->nthread or dc->nthread-1 threads.
533 cg0
= (dd
->ncg_home
*(thread
-t0_set
))/(dc
->nthread
-t0_set
);
534 cg1
= (dd
->ncg_home
*(thread
-t0_set
+1))/(dc
->nthread
-t0_set
);
536 if (thread
== t0_set
)
542 ilst
= &dc
->ils
[thread
];
546 ireqt
= &dd
->constraint_comm
->ireq
[thread
];
552 atoms_to_settles(dd
, mtop
, cginfo
, at2settle_mt
,
557 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
560 /* Combine the generate settles and requested indices */
561 for (thread
= 1; thread
< dc
->nthread
; thread
++)
569 ilst
= &dc
->ils
[thread
];
570 if (ils_local
->nr
+ ilst
->nr
> ils_local
->nalloc
)
572 ils_local
->nalloc
= over_alloc_large(ils_local
->nr
+ ilst
->nr
);
573 srenew(ils_local
->iatoms
, ils_local
->nalloc
);
575 for (ia
= 0; ia
< ilst
->nr
; ia
++)
577 ils_local
->iatoms
[ils_local
->nr
+ia
] = ilst
->iatoms
[ia
];
579 ils_local
->nr
+= ilst
->nr
;
582 ireqt
= &dd
->constraint_comm
->ireq
[thread
];
583 if (ireq
->n
+ireqt
->n
> ireq
->nalloc
)
585 ireq
->nalloc
= over_alloc_large(ireq
->n
+ireqt
->n
);
586 srenew(ireq
->ind
, ireq
->nalloc
);
588 for (ia
= 0; ia
< ireqt
->n
; ia
++)
590 ireq
->ind
[ireq
->n
+ia
] = ireqt
->ind
[ia
];
597 fprintf(debug
, "Settles: total %3d\n", ils_local
->nr
/4);
601 if (dd
->constraint_comm
)
606 setup_specat_communication(dd
, ireq
, dd
->constraint_comm
,
607 dd
->constraints
->ga2la
,
609 "constraint", " or lincs-order");
611 /* Fill in the missing indices */
612 ga2la_specat
= dd
->constraints
->ga2la
;
614 nral1
= 1 + NRAL(F_CONSTR
);
615 for (i
= 0; i
< ilc_local
->nr
; i
+= nral1
)
617 iap
= ilc_local
->iatoms
+ i
;
618 for (j
= 1; j
< nral1
; j
++)
622 iap
[j
] = gmx_hash_get_minone(ga2la_specat
, -iap
[j
]-1);
627 nral1
= 1 + NRAL(F_SETTLE
);
628 for (i
= 0; i
< ils_local
->nr
; i
+= nral1
)
630 iap
= ils_local
->iatoms
+ i
;
631 for (j
= 1; j
< nral1
; j
++)
635 iap
[j
] = gmx_hash_get_minone(ga2la_specat
, -iap
[j
]-1);
642 // Currently unreachable
649 void init_domdec_constraints(gmx_domdec_t
*dd
,
650 const gmx_mtop_t
*mtop
)
652 gmx_domdec_constraints_t
*dc
;
653 const gmx_molblock_t
*molb
;
657 fprintf(debug
, "Begin init_domdec_constraints\n");
660 snew(dd
->constraints
, 1);
661 dc
= dd
->constraints
;
663 snew(dc
->molb_con_offset
, mtop
->molblock
.size());
664 snew(dc
->molb_ncon_mol
, mtop
->molblock
.size());
667 for (size_t mb
= 0; mb
< mtop
->molblock
.size(); mb
++)
669 molb
= &mtop
->molblock
[mb
];
670 dc
->molb_con_offset
[mb
] = ncon
;
671 dc
->molb_ncon_mol
[mb
] =
672 mtop
->moltype
[molb
->type
].ilist
[F_CONSTR
].nr
/3 +
673 mtop
->moltype
[molb
->type
].ilist
[F_CONSTRNC
].nr
/3;
674 ncon
+= molb
->nmol
*dc
->molb_ncon_mol
[mb
];
679 snew(dc
->gc_req
, ncon
);
680 for (int c
= 0; c
< ncon
; c
++)
686 /* Use a hash table for the global to local index.
687 * The number of keys is a rough estimate, it will be optimized later.
689 dc
->ga2la
= gmx_hash_init(std::min(mtop
->natoms
/20,
690 mtop
->natoms
/(2*dd
->nnodes
)));
692 dc
->nthread
= gmx_omp_nthreads_get(emntDomdec
);
693 snew(dc
->ils
, dc
->nthread
);
695 dd
->constraint_comm
= specat_comm_init(dc
->nthread
);