Convert gmx_mtop_t to C++
[alexxy/gromacs.git] / src / gromacs / mdlib / broadcaststructs.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
39
40 #include "broadcaststructs.h"
41
42 #include <string.h>
43
44 #include "gromacs/gmxlib/network.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/mdlib/mdrun.h"
47 #include "gromacs/mdlib/tgroup.h"
48 #include "gromacs/mdtypes/awh-params.h"
49 #include "gromacs/mdtypes/commrec.h"
50 #include "gromacs/mdtypes/inputrec.h"
51 #include "gromacs/mdtypes/md_enums.h"
52 #include "gromacs/mdtypes/pull-params.h"
53 #include "gromacs/mdtypes/state.h"
54 #include "gromacs/topology/mtop_util.h"
55 #include "gromacs/topology/symtab.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/inmemoryserializer.h"
59 #include "gromacs/utility/keyvaluetree.h"
60 #include "gromacs/utility/keyvaluetreeserializer.h"
61 #include "gromacs/utility/smalloc.h"
62
63 static void bc_cstring(const t_commrec *cr, char **s)
64 {
65     int size = 0;
66
67     if (MASTER(cr) && *s != nullptr)
68     {
69         /* Size of the char buffer is string length + 1 for '\0' */
70         size = strlen(*s) + 1;
71     }
72     block_bc(cr, size);
73     if (size > 0)
74     {
75         if (!MASTER(cr))
76         {
77             srenew(*s, size);
78         }
79         nblock_bc(cr, size, *s);
80     }
81     else if (!MASTER(cr) && *s != nullptr)
82     {
83         sfree(*s);
84         *s = nullptr;
85     }
86 }
87
88 static void bc_string(const t_commrec *cr, t_symtab *symtab, char ***s)
89 {
90     int handle;
91
92     if (MASTER(cr))
93     {
94         handle = lookup_symtab(symtab, *s);
95     }
96     block_bc(cr, handle);
97     if (!MASTER(cr))
98     {
99         *s = get_symtab_handle(symtab, handle);
100     }
101 }
102
103 static void bc_strings(const t_commrec *cr, t_symtab *symtab, int nr, char ****nm)
104 {
105     int     i;
106     int    *handle;
107
108     snew(handle, nr);
109     if (MASTER(cr))
110     {
111         for (i = 0; (i < nr); i++)
112         {
113             handle[i] = lookup_symtab(symtab, (*nm)[i]);
114         }
115     }
116     nblock_bc(cr, nr, handle);
117
118     if (!MASTER(cr))
119     {
120         snew_bc(cr, *nm, nr);
121         for (i = 0; (i < nr); i++)
122         {
123             (*nm)[i] = get_symtab_handle(symtab, handle[i]);
124         }
125     }
126     sfree(handle);
127 }
128
129 static void bc_strings_resinfo(const t_commrec *cr, t_symtab *symtab,
130                                int nr, t_resinfo *resinfo)
131 {
132     int   i;
133     int  *handle;
134
135     snew(handle, nr);
136     if (MASTER(cr))
137     {
138         for (i = 0; (i < nr); i++)
139         {
140             handle[i] = lookup_symtab(symtab, resinfo[i].name);
141         }
142     }
143     nblock_bc(cr, nr, handle);
144
145     if (!MASTER(cr))
146     {
147         for (i = 0; (i < nr); i++)
148         {
149             resinfo[i].name = get_symtab_handle(symtab, handle[i]);
150         }
151     }
152     sfree(handle);
153 }
154
155 static void bc_symtab(const t_commrec *cr, t_symtab *symtab)
156 {
157     int       i, nr, len;
158     t_symbuf *symbuf;
159
160     block_bc(cr, symtab->nr);
161     nr = symtab->nr;
162     snew_bc(cr, symtab->symbuf, 1);
163     symbuf          = symtab->symbuf;
164     symbuf->bufsize = nr;
165     snew_bc(cr, symbuf->buf, nr);
166     for (i = 0; i < nr; i++)
167     {
168         if (MASTER(cr))
169         {
170             len = strlen(symbuf->buf[i]) + 1;
171         }
172         block_bc(cr, len);
173         snew_bc(cr, symbuf->buf[i], len);
174         nblock_bc(cr, len, symbuf->buf[i]);
175     }
176 }
177
178 static void bc_block(const t_commrec *cr, t_block *block)
179 {
180     block_bc(cr, block->nr);
181     snew_bc(cr, block->index, block->nr+1);
182     nblock_bc(cr, block->nr+1, block->index);
183 }
184
185 static void bc_blocka(const t_commrec *cr, t_blocka *block)
186 {
187     block_bc(cr, block->nr);
188     snew_bc(cr, block->index, block->nr+1);
189     nblock_bc(cr, block->nr+1, block->index);
190     block_bc(cr, block->nra);
191     if (block->nra)
192     {
193         snew_bc(cr, block->a, block->nra);
194         nblock_bc(cr, block->nra, block->a);
195     }
196 }
197
198 static void bc_grps(const t_commrec *cr, t_grps grps[])
199 {
200     int i;
201
202     for (i = 0; (i < egcNR); i++)
203     {
204         block_bc(cr, grps[i].nr);
205         snew_bc(cr, grps[i].nm_ind, grps[i].nr);
206         nblock_bc(cr, grps[i].nr, grps[i].nm_ind);
207     }
208 }
209
210 static void bc_atoms(const t_commrec *cr, t_symtab *symtab, t_atoms *atoms)
211 {
212     block_bc(cr, atoms->nr);
213     snew_bc(cr, atoms->atom, atoms->nr);
214     nblock_bc(cr, atoms->nr, atoms->atom);
215     bc_strings(cr, symtab, atoms->nr, &atoms->atomname);
216     block_bc(cr, atoms->nres);
217     snew_bc(cr, atoms->resinfo, atoms->nres);
218     nblock_bc(cr, atoms->nres, atoms->resinfo);
219     bc_strings_resinfo(cr, symtab, atoms->nres, atoms->resinfo);
220     /* QMMM requires atomtypes to be known on all nodes as well */
221     bc_strings(cr, symtab, atoms->nr, &atoms->atomtype);
222     bc_strings(cr, symtab, atoms->nr, &atoms->atomtypeB);
223 }
224
225 static void bc_groups(const t_commrec *cr, t_symtab *symtab,
226                       int natoms, gmx_groups_t *groups)
227 {
228     int g, n;
229
230     bc_grps(cr, groups->grps);
231     block_bc(cr, groups->ngrpname);
232     bc_strings(cr, symtab, groups->ngrpname, &groups->grpname);
233     for (g = 0; g < egcNR; g++)
234     {
235         if (MASTER(cr))
236         {
237             if (groups->grpnr[g])
238             {
239                 n = natoms;
240             }
241             else
242             {
243                 n = 0;
244             }
245         }
246         block_bc(cr, n);
247         if (n == 0)
248         {
249             groups->grpnr[g] = nullptr;
250         }
251         else
252         {
253             snew_bc(cr, groups->grpnr[g], n);
254             nblock_bc(cr, n, groups->grpnr[g]);
255         }
256     }
257     if (debug)
258     {
259         fprintf(debug, "after bc_groups\n");
260     }
261 }
262
263 template <typename AllocatorType>
264 static void bcastPaddedRVecVector(const t_commrec *cr, std::vector<gmx::RVec, AllocatorType> *v, int numAtoms)
265 {
266     v->resize(gmx::paddedRVecVectorSize(numAtoms));
267     nblock_bc(cr, numAtoms, as_rvec_array(v->data()));
268 }
269
270 void broadcastStateWithoutDynamics(const t_commrec *cr, t_state *state)
271 {
272     GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr), "broadcastStateWithoutDynamics should only be used for special cases without domain decomposition");
273
274     if (!PAR(cr))
275     {
276         return;
277     }
278
279     /* Broadcasts the state sizes and flags from the master to all ranks
280      * in cr->mpi_comm_mygroup.
281      */
282     block_bc(cr, state->natoms);
283     block_bc(cr, state->flags);
284
285     for (int i = 0; i < estNR; i++)
286     {
287         if (state->flags & (1 << i))
288         {
289             switch (i)
290             {
291                 case estLAMBDA:
292                     nblock_bc(cr, efptNR, state->lambda.data());
293                     break;
294                 case estFEPSTATE:
295                     block_bc(cr, state->fep_state);
296                     break;
297                 case estBOX:
298                     block_bc(cr, state->box);
299                     break;
300                 case estX:
301                     bcastPaddedRVecVector(cr, &state->x, state->natoms);
302                     break;
303                 default:
304                     GMX_RELEASE_ASSERT(false, "The state has a dynamic entry, while no dynamic entries should be present");
305                     break;
306             }
307         }
308     }
309 }
310
311 static void bc_ilists(const t_commrec *cr, t_ilist *ilist)
312 {
313     int ftype;
314
315     /* Here we only communicate the non-zero length ilists */
316     if (MASTER(cr))
317     {
318         for (ftype = 0; ftype < F_NRE; ftype++)
319         {
320             if (ilist[ftype].nr > 0)
321             {
322                 block_bc(cr, ftype);
323                 block_bc(cr, ilist[ftype].nr);
324                 nblock_bc(cr, ilist[ftype].nr, ilist[ftype].iatoms);
325             }
326         }
327         ftype = -1;
328         block_bc(cr, ftype);
329     }
330     else
331     {
332         for (ftype = 0; ftype < F_NRE; ftype++)
333         {
334             ilist[ftype].nr = 0;
335         }
336         do
337         {
338             block_bc(cr, ftype);
339             if (ftype >= 0)
340             {
341                 block_bc(cr, ilist[ftype].nr);
342                 snew_bc(cr, ilist[ftype].iatoms, ilist[ftype].nr);
343                 nblock_bc(cr, ilist[ftype].nr, ilist[ftype].iatoms);
344             }
345         }
346         while (ftype >= 0);
347     }
348
349     if (debug)
350     {
351         fprintf(debug, "after bc_ilists\n");
352     }
353 }
354
355 static void bc_cmap(const t_commrec *cr, gmx_cmap_t *cmap_grid)
356 {
357     int i, nelem, ngrid;
358
359     block_bc(cr, cmap_grid->ngrid);
360     block_bc(cr, cmap_grid->grid_spacing);
361
362     ngrid = cmap_grid->ngrid;
363     nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
364
365     if (ngrid > 0)
366     {
367         snew_bc(cr, cmap_grid->cmapdata, ngrid);
368
369         for (i = 0; i < ngrid; i++)
370         {
371             snew_bc(cr, cmap_grid->cmapdata[i].cmap, 4*nelem);
372             nblock_bc(cr, 4*nelem, cmap_grid->cmapdata[i].cmap);
373         }
374     }
375 }
376
377 static void bc_ffparams(const t_commrec *cr, gmx_ffparams_t *ffp)
378 {
379     block_bc(cr, ffp->ntypes);
380     block_bc(cr, ffp->atnr);
381     snew_bc(cr, ffp->functype, ffp->ntypes);
382     snew_bc(cr, ffp->iparams, ffp->ntypes);
383     nblock_bc(cr, ffp->ntypes, ffp->functype);
384     nblock_bc(cr, ffp->ntypes, ffp->iparams);
385     block_bc(cr, ffp->reppow);
386     block_bc(cr, ffp->fudgeQQ);
387     bc_cmap(cr, &ffp->cmap_grid);
388 }
389
390 static void bc_grpopts(const t_commrec *cr, t_grpopts *g)
391 {
392     int i, n;
393
394     block_bc(cr, g->ngtc);
395     block_bc(cr, g->ngacc);
396     block_bc(cr, g->ngfrz);
397     block_bc(cr, g->ngener);
398     snew_bc(cr, g->nrdf, g->ngtc);
399     snew_bc(cr, g->tau_t, g->ngtc);
400     snew_bc(cr, g->ref_t, g->ngtc);
401     snew_bc(cr, g->acc, g->ngacc);
402     snew_bc(cr, g->nFreeze, g->ngfrz);
403     snew_bc(cr, g->egp_flags, g->ngener*g->ngener);
404
405     nblock_bc(cr, g->ngtc, g->nrdf);
406     nblock_bc(cr, g->ngtc, g->tau_t);
407     nblock_bc(cr, g->ngtc, g->ref_t);
408     nblock_bc(cr, g->ngacc, g->acc);
409     nblock_bc(cr, g->ngfrz, g->nFreeze);
410     nblock_bc(cr, g->ngener*g->ngener, g->egp_flags);
411     snew_bc(cr, g->annealing, g->ngtc);
412     snew_bc(cr, g->anneal_npoints, g->ngtc);
413     snew_bc(cr, g->anneal_time, g->ngtc);
414     snew_bc(cr, g->anneal_temp, g->ngtc);
415     nblock_bc(cr, g->ngtc, g->annealing);
416     nblock_bc(cr, g->ngtc, g->anneal_npoints);
417     for (i = 0; (i < g->ngtc); i++)
418     {
419         n = g->anneal_npoints[i];
420         if (n > 0)
421         {
422             snew_bc(cr, g->anneal_time[i], n);
423             snew_bc(cr, g->anneal_temp[i], n);
424             nblock_bc(cr, n, g->anneal_time[i]);
425             nblock_bc(cr, n, g->anneal_temp[i]);
426         }
427     }
428
429     /* QMMM stuff, see inputrec */
430     block_bc(cr, g->ngQM);
431     snew_bc(cr, g->QMmethod, g->ngQM);
432     snew_bc(cr, g->QMbasis, g->ngQM);
433     snew_bc(cr, g->QMcharge, g->ngQM);
434     snew_bc(cr, g->QMmult, g->ngQM);
435     snew_bc(cr, g->bSH, g->ngQM);
436     snew_bc(cr, g->CASorbitals, g->ngQM);
437     snew_bc(cr, g->CASelectrons, g->ngQM);
438     snew_bc(cr, g->SAon, g->ngQM);
439     snew_bc(cr, g->SAoff, g->ngQM);
440     snew_bc(cr, g->SAsteps, g->ngQM);
441
442     if (g->ngQM)
443     {
444         nblock_bc(cr, g->ngQM, g->QMmethod);
445         nblock_bc(cr, g->ngQM, g->QMbasis);
446         nblock_bc(cr, g->ngQM, g->QMcharge);
447         nblock_bc(cr, g->ngQM, g->QMmult);
448         nblock_bc(cr, g->ngQM, g->bSH);
449         nblock_bc(cr, g->ngQM, g->CASorbitals);
450         nblock_bc(cr, g->ngQM, g->CASelectrons);
451         nblock_bc(cr, g->ngQM, g->SAon);
452         nblock_bc(cr, g->ngQM, g->SAoff);
453         nblock_bc(cr, g->ngQM, g->SAsteps);
454         /* end of QMMM stuff */
455     }
456 }
457
458 static void bc_awhBias(const t_commrec *cr, gmx::AwhBiasParams *awhBiasParams)
459 {
460     block_bc(cr, *awhBiasParams);
461
462     snew_bc(cr, awhBiasParams->dimParams, awhBiasParams->ndim);
463     nblock_bc(cr, awhBiasParams->ndim, awhBiasParams->dimParams);
464 }
465
466 static void bc_awh(const t_commrec *cr, gmx::AwhParams *awhParams)
467 {
468     int k;
469
470     block_bc(cr, *awhParams);
471     snew_bc(cr, awhParams->awhBiasParams, awhParams->numBias);
472     for (k = 0; k < awhParams->numBias; k++)
473     {
474         bc_awhBias(cr, &awhParams->awhBiasParams[k]);
475     }
476 }
477
478 static void bc_pull_group(const t_commrec *cr, t_pull_group *pgrp)
479 {
480     block_bc(cr, *pgrp);
481     if (pgrp->nat > 0)
482     {
483         snew_bc(cr, pgrp->ind, pgrp->nat);
484         nblock_bc(cr, pgrp->nat, pgrp->ind);
485     }
486     if (pgrp->nweight > 0)
487     {
488         snew_bc(cr, pgrp->weight, pgrp->nweight);
489         nblock_bc(cr, pgrp->nweight, pgrp->weight);
490     }
491 }
492
493 static void bc_pull(const t_commrec *cr, pull_params_t *pull)
494 {
495     int g;
496
497     block_bc(cr, *pull);
498     snew_bc(cr, pull->group, pull->ngroup);
499     for (g = 0; g < pull->ngroup; g++)
500     {
501         bc_pull_group(cr, &pull->group[g]);
502     }
503     snew_bc(cr, pull->coord, pull->ncoord);
504     nblock_bc(cr, pull->ncoord, pull->coord);
505     for (int c = 0; c < pull->ncoord; c++)
506     {
507         if (!MASTER(cr))
508         {
509             pull->coord[c].externalPotentialProvider = nullptr;
510         }
511         if (pull->coord[c].eType == epullEXTERNAL)
512         {
513             bc_cstring(cr, &pull->coord[c].externalPotentialProvider);
514         }
515     }
516 }
517
518 static void bc_rotgrp(const t_commrec *cr, t_rotgrp *rotg)
519 {
520     block_bc(cr, *rotg);
521     if (rotg->nat > 0)
522     {
523         snew_bc(cr, rotg->ind, rotg->nat);
524         nblock_bc(cr, rotg->nat, rotg->ind);
525         snew_bc(cr, rotg->x_ref, rotg->nat);
526         nblock_bc(cr, rotg->nat, rotg->x_ref);
527     }
528 }
529
530 static void bc_rot(const t_commrec *cr, t_rot *rot)
531 {
532     int g;
533
534     block_bc(cr, *rot);
535     snew_bc(cr, rot->grp, rot->ngrp);
536     for (g = 0; g < rot->ngrp; g++)
537     {
538         bc_rotgrp(cr, &rot->grp[g]);
539     }
540 }
541
542 static void bc_imd(const t_commrec *cr, t_IMD *imd)
543 {
544     block_bc(cr, *imd);
545     snew_bc(cr, imd->ind, imd->nat);
546     nblock_bc(cr, imd->nat, imd->ind);
547 }
548
549 static void bc_fepvals(const t_commrec *cr, t_lambda *fep)
550 {
551     int      i;
552
553     block_bc(cr, fep->nstdhdl);
554     block_bc(cr, fep->init_lambda);
555     block_bc(cr, fep->init_fep_state);
556     block_bc(cr, fep->delta_lambda);
557     block_bc(cr, fep->edHdLPrintEnergy);
558     block_bc(cr, fep->n_lambda);
559     if (fep->n_lambda > 0)
560     {
561         snew_bc(cr, fep->all_lambda, efptNR);
562         nblock_bc(cr, efptNR, fep->all_lambda);
563         for (i = 0; i < efptNR; i++)
564         {
565             snew_bc(cr, fep->all_lambda[i], fep->n_lambda);
566             nblock_bc(cr, fep->n_lambda, fep->all_lambda[i]);
567         }
568     }
569     block_bc(cr, fep->sc_alpha);
570     block_bc(cr, fep->sc_power);
571     block_bc(cr, fep->sc_r_power);
572     block_bc(cr, fep->sc_sigma);
573     block_bc(cr, fep->sc_sigma_min);
574     block_bc(cr, fep->bScCoul);
575     nblock_bc(cr, efptNR, &(fep->separate_dvdl[0]));
576     block_bc(cr, fep->dhdl_derivatives);
577     block_bc(cr, fep->dh_hist_size);
578     block_bc(cr, fep->dh_hist_spacing);
579     if (debug)
580     {
581         fprintf(debug, "after bc_fepvals\n");
582     }
583 }
584
585 static void bc_expandedvals(const t_commrec *cr, t_expanded *expand, int n_lambda)
586 {
587     block_bc(cr, expand->nstexpanded);
588     block_bc(cr, expand->elamstats);
589     block_bc(cr, expand->elmcmove);
590     block_bc(cr, expand->elmceq);
591     block_bc(cr, expand->equil_n_at_lam);
592     block_bc(cr, expand->equil_wl_delta);
593     block_bc(cr, expand->equil_ratio);
594     block_bc(cr, expand->equil_steps);
595     block_bc(cr, expand->equil_samples);
596     block_bc(cr, expand->lmc_seed);
597     block_bc(cr, expand->minvar);
598     block_bc(cr, expand->minvar_const);
599     block_bc(cr, expand->c_range);
600     block_bc(cr, expand->bSymmetrizedTMatrix);
601     block_bc(cr, expand->nstTij);
602     block_bc(cr, expand->lmc_repeats);
603     block_bc(cr, expand->lmc_forced_nstart);
604     block_bc(cr, expand->gibbsdeltalam);
605     block_bc(cr, expand->wl_scale);
606     block_bc(cr, expand->wl_ratio);
607     block_bc(cr, expand->init_wl_delta);
608     block_bc(cr, expand->bInit_weights);
609     snew_bc(cr, expand->init_lambda_weights, n_lambda);
610     nblock_bc(cr, n_lambda, expand->init_lambda_weights);
611     block_bc(cr, expand->mc_temp);
612     if (debug)
613     {
614         fprintf(debug, "after bc_expandedvals\n");
615     }
616 }
617
618 static void bc_simtempvals(const t_commrec *cr, t_simtemp *simtemp, int n_lambda)
619 {
620     block_bc(cr, simtemp->simtemp_low);
621     block_bc(cr, simtemp->simtemp_high);
622     block_bc(cr, simtemp->eSimTempScale);
623     snew_bc(cr, simtemp->temperatures, n_lambda);
624     nblock_bc(cr, n_lambda, simtemp->temperatures);
625     if (debug)
626     {
627         fprintf(debug, "after bc_simtempvals\n");
628     }
629 }
630
631
632 static void bc_swapions(const t_commrec *cr, t_swapcoords *swap)
633 {
634     block_bc(cr, *swap);
635
636     /* Broadcast atom indices for split groups, solvent group, and for all user-defined swap groups */
637     snew_bc(cr, swap->grp, swap->ngrp);
638     for (int i = 0; i < swap->ngrp; i++)
639     {
640         t_swapGroup *g = &swap->grp[i];
641
642         block_bc(cr, *g);
643         snew_bc(cr, g->ind, g->nat);
644         nblock_bc(cr, g->nat, g->ind);
645
646         int len = 0;
647         if (MASTER(cr))
648         {
649             len = strlen(g->molname);
650         }
651         block_bc(cr, len);
652         snew_bc(cr, g->molname, len);
653         nblock_bc(cr, len, g->molname);
654     }
655 }
656
657
658 static void bc_inputrec(const t_commrec *cr, t_inputrec *inputrec)
659 {
660     // Note that this overwrites pointers in inputrec, so all pointer fields
661     // Must be initialized separately below.
662     block_bc(cr, *inputrec);
663     if (SIMMASTER(cr))
664     {
665         gmx::InMemorySerializer serializer;
666         gmx::serializeKeyValueTree(*inputrec->params, &serializer);
667         std::vector<char>       buffer = serializer.finishAndGetBuffer();
668         size_t                  size   = buffer.size();
669         block_bc(cr, size);
670         nblock_bc(cr, size, buffer.data());
671     }
672     else
673     {
674         // block_bc() above overwrites the old pointer, so set it to a
675         // reasonable value in case code below throws.
676         // cppcheck-suppress redundantAssignment
677         inputrec->params = nullptr;
678         std::vector<char> buffer;
679         size_t            size;
680         block_bc(cr, size);
681         nblock_abc(cr, size, &buffer);
682         gmx::InMemoryDeserializer serializer(buffer);
683         // cppcheck-suppress redundantAssignment
684         inputrec->params = new gmx::KeyValueTreeObject(
685                     gmx::deserializeKeyValueTree(&serializer));
686     }
687
688     bc_grpopts(cr, &(inputrec->opts));
689
690     /* even if efep is efepNO, we need to initialize to make sure that
691      * n_lambda is set to zero */
692
693     snew_bc(cr, inputrec->fepvals, 1);
694     if (inputrec->efep != efepNO || inputrec->bSimTemp)
695     {
696         bc_fepvals(cr, inputrec->fepvals);
697     }
698     /* need to initialize this as well because of data checked for in the logic */
699     snew_bc(cr, inputrec->expandedvals, 1);
700     if (inputrec->bExpanded)
701     {
702         bc_expandedvals(cr, inputrec->expandedvals, inputrec->fepvals->n_lambda);
703     }
704     snew_bc(cr, inputrec->simtempvals, 1);
705     if (inputrec->bSimTemp)
706     {
707         bc_simtempvals(cr, inputrec->simtempvals, inputrec->fepvals->n_lambda);
708     }
709     if (inputrec->bPull)
710     {
711         snew_bc(cr, inputrec->pull, 1);
712         bc_pull(cr, inputrec->pull);
713     }
714     if (inputrec->bDoAwh)
715     {
716         snew_bc(cr, inputrec->awhParams, 1);
717         bc_awh(cr, inputrec->awhParams);
718     }
719
720     if (inputrec->bRot)
721     {
722         snew_bc(cr, inputrec->rot, 1);
723         bc_rot(cr, inputrec->rot);
724     }
725     if (inputrec->bIMD)
726     {
727         snew_bc(cr, inputrec->imd, 1);
728         bc_imd(cr, inputrec->imd);
729     }
730     if (inputrec->eSwapCoords != eswapNO)
731     {
732         snew_bc(cr, inputrec->swap, 1);
733         bc_swapions(cr, inputrec->swap);
734     }
735 }
736
737 static void bc_moltype(const t_commrec *cr, t_symtab *symtab,
738                        gmx_moltype_t *moltype)
739 {
740     bc_string(cr, symtab, &moltype->name);
741     bc_atoms(cr, symtab, &moltype->atoms);
742     if (debug)
743     {
744         fprintf(debug, "after bc_atoms\n");
745     }
746
747     bc_ilists(cr, moltype->ilist);
748     bc_block(cr, &moltype->cgs);
749     bc_blocka(cr, &moltype->excls);
750 }
751
752 static void bc_molblock(const t_commrec *cr, gmx_molblock_t *molb)
753 {
754     block_bc(cr, *molb);
755     if (molb->nposres_xA > 0)
756     {
757         snew_bc(cr, molb->posres_xA, molb->nposres_xA);
758         nblock_bc(cr, molb->nposres_xA*DIM, molb->posres_xA[0]);
759     }
760     if (molb->nposres_xB > 0)
761     {
762         snew_bc(cr, molb->posres_xB, molb->nposres_xB);
763         nblock_bc(cr, molb->nposres_xB*DIM, molb->posres_xB[0]);
764     }
765     if (debug)
766     {
767         fprintf(debug, "after bc_molblock\n");
768     }
769 }
770
771 static void bc_atomtypes(const t_commrec *cr, t_atomtypes *atomtypes)
772 {
773     block_bc(cr, atomtypes->nr);
774 }
775
776 /*! \brief Broadcasts ir and mtop from the master to all nodes in
777  * cr->mpi_comm_mygroup. */
778 static
779 void bcast_ir_mtop(const t_commrec *cr, t_inputrec *inputrec, gmx_mtop_t *mtop)
780 {
781     if (debug)
782     {
783         fprintf(debug, "in bc_data\n");
784     }
785     bc_inputrec(cr, inputrec);
786     if (debug)
787     {
788         fprintf(debug, "after bc_inputrec\n");
789     }
790     bc_symtab(cr, &mtop->symtab);
791     if (debug)
792     {
793         fprintf(debug, "after bc_symtab\n");
794     }
795     bc_string(cr, &mtop->symtab, &mtop->name);
796     if (debug)
797     {
798         fprintf(debug, "after bc_name\n");
799     }
800
801     bc_ffparams(cr, &mtop->ffparams);
802
803     int nmoltype = mtop->moltype.size();
804     block_bc(cr, nmoltype);
805     mtop->moltype.resize(nmoltype);
806     for (gmx_moltype_t &moltype : mtop->moltype)
807     {
808         bc_moltype(cr, &mtop->symtab, &moltype);
809     }
810
811     block_bc(cr, mtop->bIntermolecularInteractions);
812     if (mtop->bIntermolecularInteractions)
813     {
814         snew_bc(cr, mtop->intermolecular_ilist, F_NRE);
815         bc_ilists(cr, mtop->intermolecular_ilist);
816     }
817
818     int nmolblock = mtop->molblock.size();
819     block_bc(cr, nmolblock);
820     mtop->molblock.resize(nmolblock);
821     for (gmx_molblock_t &molblock : mtop->molblock)
822     {
823         bc_molblock(cr, &molblock);
824     }
825
826     block_bc(cr, mtop->natoms);
827
828     bc_atomtypes(cr, &mtop->atomtypes);
829
830     bc_groups(cr, &mtop->symtab, mtop->natoms, &mtop->groups);
831
832     GMX_RELEASE_ASSERT(!MASTER(cr) || mtop->haveMoleculeIndices, "mtop should have valid molecule indices");
833     mtop->haveMoleculeIndices = true;
834 }
835
836 void init_parallel(t_commrec *cr, t_inputrec *inputrec,
837                    gmx_mtop_t *mtop)
838 {
839     bcast_ir_mtop(cr, inputrec, mtop);
840 }