91565d67117875cf2991b7fc20432ed0c4826863
[alexxy/gromacs.git] / src / gromacs / domdec / domdec_constraints.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2006,2007,2008,2009,2010,2012,2013,2014,2015,2016, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  *
38  * \brief This file implements functions for domdec to use
39  * while managing inter-atomic constraints.
40  *
41  * \author Berk Hess <hess@kth.se>
42  * \ingroup module_domdec
43  */
44
45 #include "gmxpre.h"
46
47 #include "domdec_constraints.h"
48
49 #include <assert.h>
50
51 #include <algorithm>
52
53 #include "gromacs/domdec/domdec.h"
54 #include "gromacs/domdec/domdec_struct.h"
55 #include "gromacs/domdec/ga2la.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/mdlib/constr.h"
58 #include "gromacs/mdlib/gmx_omp_nthreads.h"
59 #include "gromacs/mdtypes/commrec.h"
60 #include "gromacs/pbcutil/ishift.h"
61 #include "gromacs/topology/mtop_lookup.h"
62 #include "gromacs/utility/exceptions.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/gmxassert.h"
65 #include "gromacs/utility/smalloc.h"
66
67 #include "domdec_specatomcomm.h"
68 #include "hash.h"
69
70 /*! \brief Struct used during constraint setup with domain decomposition */
71 struct gmx_domdec_constraints_t {
72     //! @cond Doxygen_Suppress
73     int         *molb_con_offset; /**< Offset in the constraint array for each molblock */
74     int         *molb_ncon_mol;   /**< The number of constraints per molecule for each molblock */
75
76     int          ncon;            /**< The fully local and conneced constraints */
77     /* The global constraint number, only required for clearing gc_req */
78     int         *con_gl;          /**< Global constraint indices for local constraints */
79     int         *con_nlocat;      /**< Number of local atoms (2/1/0) for each constraint */
80     int          con_nalloc;      /**< Allocation size for \p con_gl and \p con_nlocat */
81
82     char        *gc_req;          /**< Boolean that tells if a global constraint index has been requested; note: size global #constraints */
83     gmx_hash_t  *ga2la;           /**< Global to local communicated constraint atom only index */
84
85     /* Multi-threading stuff */
86     int      nthread;           /**< Number of threads used for DD constraint setup */
87     t_ilist *ils;               /**< Constraint ilist working arrays, size \p nthread */
88     //! @endcond
89 };
90
91 void dd_move_x_constraints(gmx_domdec_t *dd, matrix box,
92                            rvec *x0, rvec *x1, gmx_bool bX1IsCoord)
93 {
94     if (dd->constraint_comm)
95     {
96         dd_move_x_specat(dd, dd->constraint_comm, box, x0, x1, bX1IsCoord);
97     }
98 }
99
100 int *dd_constraints_nlocalatoms(gmx_domdec_t *dd)
101 {
102     if (dd->constraints)
103     {
104         return dd->constraints->con_nlocat;
105     }
106     else
107     {
108         return NULL;
109     }
110 }
111
112 void dd_clear_local_constraint_indices(gmx_domdec_t *dd)
113 {
114     gmx_domdec_constraints_t *dc;
115     int i;
116
117     dc = dd->constraints;
118
119     for (i = 0; i < dc->ncon; i++)
120     {
121         dc->gc_req[dc->con_gl[i]] = 0;
122     }
123
124     if (dd->constraint_comm)
125     {
126         gmx_hash_clear_and_optimize(dc->ga2la);
127     }
128 }
129
130 void dd_clear_local_vsite_indices(gmx_domdec_t *dd)
131 {
132     if (dd->vsite_comm)
133     {
134         gmx_hash_clear_and_optimize(dd->ga2la_vsite);
135     }
136 }
137
138 /*! \brief Walks over the constraints out from the local atoms into the non-local atoms and adds them to a list */
139 static void walk_out(int con, int con_offset, int a, int offset, int nrec,
140                      int ncon1, const t_iatom *ia1, const t_iatom *ia2,
141                      const t_blocka *at2con,
142                      const gmx_ga2la_t *ga2la, gmx_bool bHomeConnect,
143                      gmx_domdec_constraints_t *dc,
144                      gmx_domdec_specat_comm_t *dcc,
145                      t_ilist *il_local,
146                      ind_req_t *ireq)
147 {
148     int            a1_gl, a2_gl, a_loc, i, coni, b;
149     const t_iatom *iap;
150
151     if (dc->gc_req[con_offset+con] == 0)
152     {
153         /* Add this non-home constraint to the list */
154         if (dc->ncon+1 > dc->con_nalloc)
155         {
156             dc->con_nalloc = over_alloc_large(dc->ncon+1);
157             srenew(dc->con_gl, dc->con_nalloc);
158             srenew(dc->con_nlocat, dc->con_nalloc);
159         }
160         dc->con_gl[dc->ncon]       = con_offset + con;
161         dc->con_nlocat[dc->ncon]   = (bHomeConnect ? 1 : 0);
162         dc->gc_req[con_offset+con] = 1;
163         if (il_local->nr + 3 > il_local->nalloc)
164         {
165             il_local->nalloc = over_alloc_dd(il_local->nr+3);
166             srenew(il_local->iatoms, il_local->nalloc);
167         }
168         iap = constr_iatomptr(ncon1, ia1, ia2, con);
169         il_local->iatoms[il_local->nr++] = iap[0];
170         a1_gl = offset + iap[1];
171         a2_gl = offset + iap[2];
172         /* The following indexing code can probably be optizimed */
173         if (ga2la_get_home(ga2la, a1_gl, &a_loc))
174         {
175             il_local->iatoms[il_local->nr++] = a_loc;
176         }
177         else
178         {
179             /* We set this index later */
180             il_local->iatoms[il_local->nr++] = -a1_gl - 1;
181         }
182         if (ga2la_get_home(ga2la, a2_gl, &a_loc))
183         {
184             il_local->iatoms[il_local->nr++] = a_loc;
185         }
186         else
187         {
188             /* We set this index later */
189             il_local->iatoms[il_local->nr++] = -a2_gl - 1;
190         }
191         dc->ncon++;
192     }
193     /* Check to not ask for the same atom more than once */
194     if (gmx_hash_get_minone(dc->ga2la, offset+a) == -1)
195     {
196         assert(dcc);
197         /* Add this non-home atom to the list */
198         if (ireq->n+1 > ireq->nalloc)
199         {
200             ireq->nalloc = over_alloc_large(ireq->n+1);
201             srenew(ireq->ind, ireq->nalloc);
202         }
203         ireq->ind[ireq->n++] = offset + a;
204         /* Temporarily mark with -2, we get the index later */
205         gmx_hash_set(dc->ga2la, offset+a, -2);
206     }
207
208     if (nrec > 0)
209     {
210         for (i = at2con->index[a]; i < at2con->index[a+1]; i++)
211         {
212             coni = at2con->a[i];
213             if (coni != con)
214             {
215                 /* Walk further */
216                 iap = constr_iatomptr(ncon1, ia1, ia2, coni);
217                 if (a == iap[1])
218                 {
219                     b = iap[2];
220                 }
221                 else
222                 {
223                     b = iap[1];
224                 }
225                 if (!ga2la_get_home(ga2la, offset+b, &a_loc))
226                 {
227                     walk_out(coni, con_offset, b, offset, nrec-1,
228                              ncon1, ia1, ia2, at2con,
229                              ga2la, FALSE, dc, dcc, il_local, ireq);
230                 }
231             }
232         }
233     }
234 }
235
236 /*! \brief Looks up SETTLE constraints for a range of charge-groups */
237 static void atoms_to_settles(gmx_domdec_t *dd,
238                              const gmx_mtop_t *mtop,
239                              const int *cginfo,
240                              const int **at2settle_mt,
241                              int cg_start, int cg_end,
242                              t_ilist *ils_local,
243                              ind_req_t *ireq)
244 {
245     gmx_ga2la_t *ga2la = dd->ga2la;
246     int          nral  = NRAL(F_SETTLE);
247
248     int          mb    = 0;
249     for (int cg = cg_start; cg < cg_end; cg++)
250     {
251         if (GET_CGINFO_SETTLE(cginfo[cg]))
252         {
253             for (int a = dd->cgindex[cg]; a < dd->cgindex[cg+1]; a++)
254             {
255                 int a_gl = dd->gatindex[a];
256                 int a_mol;
257                 mtopGetMolblockIndex(mtop, a_gl, &mb, NULL, &a_mol);
258
259                 const gmx_molblock_t *molb   = &mtop->molblock[mb];
260                 int                   settle = at2settle_mt[molb->type][a_mol];
261
262                 if (settle >= 0)
263                 {
264                     int      offset  = a_gl - a_mol;
265
266                     t_iatom *ia1     = mtop->moltype[molb->type].ilist[F_SETTLE].iatoms;
267
268                     int      a_gls[3], a_locs[3];
269                     gmx_bool bAssign = FALSE;
270                     int      nlocal  = 0;
271                     for (int sa = 0; sa < nral; sa++)
272                     {
273                         int a_glsa = offset + ia1[settle*(1+nral)+1+sa];
274                         a_gls[sa]  = a_glsa;
275                         if (ga2la_get_home(ga2la, a_glsa, &a_locs[sa]))
276                         {
277                             if (nlocal == 0 && a_gl == a_glsa)
278                             {
279                                 bAssign = TRUE;
280                             }
281                             nlocal++;
282                         }
283                     }
284
285                     if (bAssign)
286                     {
287                         if (ils_local->nr+1+nral > ils_local->nalloc)
288                         {
289                             ils_local->nalloc = over_alloc_dd(ils_local->nr+1+nral);
290                             srenew(ils_local->iatoms, ils_local->nalloc);
291                         }
292
293                         ils_local->iatoms[ils_local->nr++] = ia1[settle*4];
294
295                         for (int sa = 0; sa < nral; sa++)
296                         {
297                             if (ga2la_get_home(ga2la, a_gls[sa], &a_locs[sa]))
298                             {
299                                 ils_local->iatoms[ils_local->nr++] = a_locs[sa];
300                             }
301                             else
302                             {
303                                 ils_local->iatoms[ils_local->nr++] = -a_gls[sa] - 1;
304                                 /* Add this non-home atom to the list */
305                                 if (ireq->n+1 > ireq->nalloc)
306                                 {
307                                     ireq->nalloc = over_alloc_large(ireq->n+1);
308                                     srenew(ireq->ind, ireq->nalloc);
309                                 }
310                                 ireq->ind[ireq->n++] = a_gls[sa];
311                                 /* A check on double atom requests is
312                                  * not required for settle.
313                                  */
314                             }
315                         }
316                     }
317                 }
318             }
319         }
320     }
321 }
322
323 /*! \brief Looks up constraint for the local atoms */
324 static void atoms_to_constraints(gmx_domdec_t *dd,
325                                  const gmx_mtop_t *mtop,
326                                  const int *cginfo,
327                                  const t_blocka *at2con_mt, int nrec,
328                                  t_ilist *ilc_local,
329                                  ind_req_t *ireq)
330 {
331     const t_blocka             *at2con;
332     int                         ncon1;
333     t_iatom                    *ia1, *ia2, *iap;
334     int                         a_loc, b_lo, offset, b_mol, i, con, con_offset;
335
336     gmx_domdec_constraints_t   *dc     = dd->constraints;
337     gmx_domdec_specat_comm_t   *dcc    = dd->constraint_comm;
338
339     gmx_ga2la_t                *ga2la  = dd->ga2la;
340
341     int mb    = 0;
342     int nhome = 0;
343     for (int cg = 0; cg < dd->ncg_home; cg++)
344     {
345         if (GET_CGINFO_CONSTR(cginfo[cg]))
346         {
347             for (int a = dd->cgindex[cg]; a < dd->cgindex[cg+1]; a++)
348             {
349                 int a_gl = dd->gatindex[a];
350                 int molnr, a_mol;
351                 mtopGetMolblockIndex(mtop, a_gl, &mb, &molnr, &a_mol);
352
353                 const gmx_molblock_t *molb = &mtop->molblock[mb];
354
355                 ncon1 = mtop->moltype[molb->type].ilist[F_CONSTR].nr/NRAL(F_SETTLE);
356
357                 ia1 = mtop->moltype[molb->type].ilist[F_CONSTR].iatoms;
358                 ia2 = mtop->moltype[molb->type].ilist[F_CONSTRNC].iatoms;
359
360                 /* Calculate the global constraint number offset for the molecule.
361                  * This is only required for the global index to make sure
362                  * that we use each constraint only once.
363                  */
364                 con_offset =
365                     dc->molb_con_offset[mb] + molnr*dc->molb_ncon_mol[mb];
366
367                 /* The global atom number offset for this molecule */
368                 offset = a_gl - a_mol;
369                 at2con = &at2con_mt[molb->type];
370                 for (i = at2con->index[a_mol]; i < at2con->index[a_mol+1]; i++)
371                 {
372                     con = at2con->a[i];
373                     iap = constr_iatomptr(ncon1, ia1, ia2, con);
374                     if (a_mol == iap[1])
375                     {
376                         b_mol = iap[2];
377                     }
378                     else
379                     {
380                         b_mol = iap[1];
381                     }
382                     if (ga2la_get_home(ga2la, offset+b_mol, &a_loc))
383                     {
384                         /* Add this fully home constraint at the first atom */
385                         if (a_mol < b_mol)
386                         {
387                             if (dc->ncon+1 > dc->con_nalloc)
388                             {
389                                 dc->con_nalloc = over_alloc_large(dc->ncon+1);
390                                 srenew(dc->con_gl, dc->con_nalloc);
391                                 srenew(dc->con_nlocat, dc->con_nalloc);
392                             }
393                             dc->con_gl[dc->ncon]     = con_offset + con;
394                             dc->con_nlocat[dc->ncon] = 2;
395                             if (ilc_local->nr + 3 > ilc_local->nalloc)
396                             {
397                                 ilc_local->nalloc = over_alloc_dd(ilc_local->nr + 3);
398                                 srenew(ilc_local->iatoms, ilc_local->nalloc);
399                             }
400                             b_lo = a_loc;
401                             ilc_local->iatoms[ilc_local->nr++] = iap[0];
402                             ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? a    : b_lo);
403                             ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? b_lo : a   );
404                             dc->ncon++;
405                             nhome++;
406                         }
407                     }
408                     else
409                     {
410                         /* We need the nrec constraints coupled to this constraint,
411                          * so we need to walk out of the home cell by nrec+1 atoms,
412                          * since already atom bg is not locally present.
413                          * Therefore we call walk_out with nrec recursions to go
414                          * after this first call.
415                          */
416                         walk_out(con, con_offset, b_mol, offset, nrec,
417                                  ncon1, ia1, ia2, at2con,
418                                  dd->ga2la, TRUE, dc, dcc, ilc_local, ireq);
419                     }
420                 }
421             }
422         }
423     }
424
425     if (debug)
426     {
427         fprintf(debug,
428                 "Constraints: home %3d border %3d atoms: %3d\n",
429                 nhome, dc->ncon-nhome,
430                 dd->constraint_comm ? ireq->n : 0);
431     }
432 }
433
434 int dd_make_local_constraints(gmx_domdec_t *dd, int at_start,
435                               const struct gmx_mtop_t *mtop,
436                               const int *cginfo,
437                               gmx_constr_t constr, int nrec,
438                               t_ilist *il_local)
439 {
440     gmx_domdec_constraints_t   *dc;
441     t_ilist                    *ilc_local, *ils_local;
442     ind_req_t                  *ireq;
443     const t_blocka             *at2con_mt;
444     const int                 **at2settle_mt;
445     gmx_hash_t                 *ga2la_specat;
446     int at_end, i, j;
447     t_iatom                    *iap;
448
449     // This code should not be called unless this condition is true,
450     // because that's the only time init_domdec_constraints is
451     // called...
452     GMX_RELEASE_ASSERT(dd->bInterCGcons || dd->bInterCGsettles, "dd_make_local_constraints called when there are no local constraints");
453     // ... and init_domdec_constraints always sets
454     // dd->constraint_comm...
455     GMX_RELEASE_ASSERT(dd->constraint_comm, "Invalid use of dd_make_local_constraints before construction of constraint_comm");
456     // ... which static analysis needs to be reassured about, because
457     // otherwise, when dd->bInterCGsettles is
458     // true. dd->constraint_comm is unilaterally dereferenced before
459     // the call to atoms_to_settles.
460
461     dc = dd->constraints;
462
463     ilc_local = &il_local[F_CONSTR];
464     ils_local = &il_local[F_SETTLE];
465
466     dc->ncon      = 0;
467     ilc_local->nr = 0;
468     if (dd->constraint_comm)
469     {
470         at2con_mt = atom2constraints_moltype(constr);
471         ireq      = &dd->constraint_comm->ireq[0];
472         ireq->n   = 0;
473     }
474     else
475     {
476         // Currently unreachable
477         at2con_mt = NULL;
478         ireq      = NULL;
479     }
480
481     if (dd->bInterCGsettles)
482     {
483         at2settle_mt  = atom2settle_moltype(constr);
484         ils_local->nr = 0;
485     }
486     else
487     {
488         /* Settle works inside charge groups, we assigned them already */
489         at2settle_mt = NULL;
490     }
491
492     if (at2settle_mt == NULL)
493     {
494         atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec,
495                              ilc_local, ireq);
496     }
497     else
498     {
499         int t0_set;
500         int thread;
501
502         /* Do the constraints, if present, on the first thread.
503          * Do the settles on all other threads.
504          */
505         t0_set = ((at2con_mt != NULL && dc->nthread > 1) ? 1 : 0);
506
507 #pragma omp parallel for num_threads(dc->nthread) schedule(static)
508         for (thread = 0; thread < dc->nthread; thread++)
509         {
510             try
511             {
512                 if (at2con_mt && thread == 0)
513                 {
514                     atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec,
515                                          ilc_local, ireq);
516                 }
517
518                 if (thread >= t0_set)
519                 {
520                     int        cg0, cg1;
521                     t_ilist   *ilst;
522                     ind_req_t *ireqt;
523
524                     /* Distribute the settle check+assignments over
525                      * dc->nthread or dc->nthread-1 threads.
526                      */
527                     cg0 = (dd->ncg_home*(thread-t0_set  ))/(dc->nthread-t0_set);
528                     cg1 = (dd->ncg_home*(thread-t0_set+1))/(dc->nthread-t0_set);
529
530                     if (thread == t0_set)
531                     {
532                         ilst = ils_local;
533                     }
534                     else
535                     {
536                         ilst = &dc->ils[thread];
537                     }
538                     ilst->nr = 0;
539
540                     ireqt = &dd->constraint_comm->ireq[thread];
541                     if (thread > 0)
542                     {
543                         ireqt->n = 0;
544                     }
545
546                     atoms_to_settles(dd, mtop, cginfo, at2settle_mt,
547                                      cg0, cg1,
548                                      ilst, ireqt);
549                 }
550             }
551             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
552         }
553
554         /* Combine the generate settles and requested indices */
555         for (thread = 1; thread < dc->nthread; thread++)
556         {
557             t_ilist   *ilst;
558             ind_req_t *ireqt;
559             int        ia;
560
561             if (thread > t0_set)
562             {
563                 ilst = &dc->ils[thread];
564                 if (ils_local->nr + ilst->nr > ils_local->nalloc)
565                 {
566                     ils_local->nalloc = over_alloc_large(ils_local->nr + ilst->nr);
567                     srenew(ils_local->iatoms, ils_local->nalloc);
568                 }
569                 for (ia = 0; ia < ilst->nr; ia++)
570                 {
571                     ils_local->iatoms[ils_local->nr+ia] = ilst->iatoms[ia];
572                 }
573                 ils_local->nr += ilst->nr;
574             }
575
576             ireqt = &dd->constraint_comm->ireq[thread];
577             if (ireq->n+ireqt->n > ireq->nalloc)
578             {
579                 ireq->nalloc = over_alloc_large(ireq->n+ireqt->n);
580                 srenew(ireq->ind, ireq->nalloc);
581             }
582             for (ia = 0; ia < ireqt->n; ia++)
583             {
584                 ireq->ind[ireq->n+ia] = ireqt->ind[ia];
585             }
586             ireq->n += ireqt->n;
587         }
588
589         if (debug)
590         {
591             fprintf(debug, "Settles: total %3d\n", ils_local->nr/4);
592         }
593     }
594
595     if (dd->constraint_comm)
596     {
597         int nral1;
598
599         at_end =
600             setup_specat_communication(dd, ireq, dd->constraint_comm,
601                                        dd->constraints->ga2la,
602                                        at_start, 2,
603                                        "constraint", " or lincs-order");
604
605         /* Fill in the missing indices */
606         ga2la_specat = dd->constraints->ga2la;
607
608         nral1 = 1 + NRAL(F_CONSTR);
609         for (i = 0; i < ilc_local->nr; i += nral1)
610         {
611             iap = ilc_local->iatoms + i;
612             for (j = 1; j < nral1; j++)
613             {
614                 if (iap[j] < 0)
615                 {
616                     iap[j] = gmx_hash_get_minone(ga2la_specat, -iap[j]-1);
617                 }
618             }
619         }
620
621         nral1 = 1 + NRAL(F_SETTLE);
622         for (i = 0; i < ils_local->nr; i += nral1)
623         {
624             iap = ils_local->iatoms + i;
625             for (j = 1; j < nral1; j++)
626             {
627                 if (iap[j] < 0)
628                 {
629                     iap[j] = gmx_hash_get_minone(ga2la_specat, -iap[j]-1);
630                 }
631             }
632         }
633     }
634     else
635     {
636         // Currently unreachable
637         at_end = at_start;
638     }
639
640     return at_end;
641 }
642
643 void init_domdec_constraints(gmx_domdec_t     *dd,
644                              const gmx_mtop_t *mtop)
645 {
646     gmx_domdec_constraints_t *dc;
647     const gmx_molblock_t     *molb;
648     int mb, ncon, c;
649
650     if (debug)
651     {
652         fprintf(debug, "Begin init_domdec_constraints\n");
653     }
654
655     snew(dd->constraints, 1);
656     dc = dd->constraints;
657
658     snew(dc->molb_con_offset, mtop->nmolblock);
659     snew(dc->molb_ncon_mol, mtop->nmolblock);
660
661     ncon = 0;
662     for (mb = 0; mb < mtop->nmolblock; mb++)
663     {
664         molb                    = &mtop->molblock[mb];
665         dc->molb_con_offset[mb] = ncon;
666         dc->molb_ncon_mol[mb]   =
667             mtop->moltype[molb->type].ilist[F_CONSTR].nr/3 +
668             mtop->moltype[molb->type].ilist[F_CONSTRNC].nr/3;
669         ncon += molb->nmol*dc->molb_ncon_mol[mb];
670     }
671
672     if (ncon > 0)
673     {
674         snew(dc->gc_req, ncon);
675         for (c = 0; c < ncon; c++)
676         {
677             dc->gc_req[c] = 0;
678         }
679     }
680
681     /* Use a hash table for the global to local index.
682      * The number of keys is a rough estimate, it will be optimized later.
683      */
684     dc->ga2la = gmx_hash_init(std::min(mtop->natoms/20,
685                                        mtop->natoms/(2*dd->nnodes)));
686
687     dc->nthread = gmx_omp_nthreads_get(emntDomdec);
688     snew(dc->ils, dc->nthread);
689
690     dd->constraint_comm = specat_comm_init(dc->nthread);
691 }