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