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