1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * GROningen Mixture of Alchemy and Childrens' Stories
35 /* This file is completely threadsafe - keep it that way! */
47 #include "gmx_fatal.h"
52 #define block_bc(cr, d) gmx_bcast( sizeof(d), &(d), (cr))
53 /* Probably the test for (nr) > 0 in the next macro is only needed
54 * on BlueGene(/L), where IBM's MPI_Bcast will segfault after
55 * dereferencing a null pointer, even when no data is to be transferred. */
56 #define nblock_bc(cr, nr, d) { if ((nr) > 0) {gmx_bcast((nr)*sizeof((d)[0]), (d), (cr)); }}
57 #define snew_bc(cr, d, nr) { if (!MASTER(cr)) {snew((d), (nr)); }}
58 /* Dirty macro with bAlloc not as an argument */
59 #define nblock_abc(cr, nr, d) { if (bAlloc) {snew((d), (nr)); } nblock_bc(cr, (nr), (d)); }
61 static void bc_string(const t_commrec *cr, t_symtab *symtab, char ***s)
67 handle = lookup_symtab(symtab, *s);
72 *s = get_symtab_handle(symtab, handle);
76 static void bc_strings(const t_commrec *cr, t_symtab *symtab, int nr, char ****nm)
86 for (i = 0; (i < nr); i++)
88 handle[i] = lookup_symtab(symtab, NM[i]);
91 nblock_bc(cr, nr, handle);
97 for (i = 0; (i < nr); i++)
99 (*nm)[i] = get_symtab_handle(symtab, handle[i]);
105 static void bc_strings_resinfo(const t_commrec *cr, t_symtab *symtab,
106 int nr, t_resinfo *resinfo)
114 for (i = 0; (i < nr); i++)
116 handle[i] = lookup_symtab(symtab, resinfo[i].name);
119 nblock_bc(cr, nr, handle);
123 for (i = 0; (i < nr); i++)
125 resinfo[i].name = get_symtab_handle(symtab, handle[i]);
131 static void bc_symtab(const t_commrec *cr, t_symtab *symtab)
136 block_bc(cr, symtab->nr);
138 snew_bc(cr, symtab->symbuf, 1);
139 symbuf = symtab->symbuf;
140 symbuf->bufsize = nr;
141 snew_bc(cr, symbuf->buf, nr);
142 for (i = 0; i < nr; i++)
146 len = strlen(symbuf->buf[i]) + 1;
149 snew_bc(cr, symbuf->buf[i], len);
150 nblock_bc(cr, len, symbuf->buf[i]);
154 static void bc_block(const t_commrec *cr, t_block *block)
156 block_bc(cr, block->nr);
157 snew_bc(cr, block->index, block->nr+1);
158 nblock_bc(cr, block->nr+1, block->index);
161 static void bc_blocka(const t_commrec *cr, t_blocka *block)
163 block_bc(cr, block->nr);
164 snew_bc(cr, block->index, block->nr+1);
165 nblock_bc(cr, block->nr+1, block->index);
166 block_bc(cr, block->nra);
169 snew_bc(cr, block->a, block->nra);
170 nblock_bc(cr, block->nra, block->a);
174 static void bc_grps(const t_commrec *cr, t_grps grps[])
178 for (i = 0; (i < egcNR); i++)
180 block_bc(cr, grps[i].nr);
181 snew_bc(cr, grps[i].nm_ind, grps[i].nr);
182 nblock_bc(cr, grps[i].nr, grps[i].nm_ind);
186 static void bc_atoms(const t_commrec *cr, t_symtab *symtab, t_atoms *atoms)
190 block_bc(cr, atoms->nr);
191 snew_bc(cr, atoms->atom, atoms->nr);
192 nblock_bc(cr, atoms->nr, atoms->atom);
193 bc_strings(cr, symtab, atoms->nr, &atoms->atomname);
194 block_bc(cr, atoms->nres);
195 snew_bc(cr, atoms->resinfo, atoms->nres);
196 nblock_bc(cr, atoms->nres, atoms->resinfo);
197 bc_strings_resinfo(cr, symtab, atoms->nres, atoms->resinfo);
198 /* QMMM requires atomtypes to be known on all nodes as well */
199 bc_strings(cr, symtab, atoms->nr, &atoms->atomtype);
200 bc_strings(cr, symtab, atoms->nr, &atoms->atomtypeB);
203 static void bc_groups(const t_commrec *cr, t_symtab *symtab,
204 int natoms, gmx_groups_t *groups)
209 bc_grps(cr, groups->grps);
210 block_bc(cr, groups->ngrpname);
211 bc_strings(cr, symtab, groups->ngrpname, &groups->grpname);
212 for (g = 0; g < egcNR; g++)
216 if (groups->grpnr[g])
228 groups->grpnr[g] = NULL;
232 snew_bc(cr, groups->grpnr[g], n);
233 nblock_bc(cr, n, groups->grpnr[g]);
238 fprintf(debug, "after bc_groups\n");
242 void bcast_state_setup(const t_commrec *cr, t_state *state)
244 block_bc(cr, state->natoms);
245 block_bc(cr, state->ngtc);
246 block_bc(cr, state->nnhpres);
247 block_bc(cr, state->nhchainlength);
248 block_bc(cr, state->nrng);
249 block_bc(cr, state->nrngi);
250 block_bc(cr, state->flags);
251 if (state->lambda == NULL)
253 snew_bc(cr, state->lambda, efptNR)
257 void bcast_state(const t_commrec *cr, t_state *state, gmx_bool bAlloc)
261 bcast_state_setup(cr, state);
263 nnht = (state->ngtc)*(state->nhchainlength);
264 nnhtp = (state->nnhpres)*(state->nhchainlength);
272 state->nalloc = state->natoms;
274 for (i = 0; i < estNR; i++)
276 if (state->flags & (1<<i))
280 case estLAMBDA: nblock_bc(cr, efptNR, state->lambda); break;
281 case estFEPSTATE: block_bc(cr, state->fep_state); break;
282 case estBOX: block_bc(cr, state->box); break;
283 case estBOX_REL: block_bc(cr, state->box_rel); break;
284 case estBOXV: block_bc(cr, state->boxv); break;
285 case estPRES_PREV: block_bc(cr, state->pres_prev); break;
286 case estSVIR_PREV: block_bc(cr, state->svir_prev); break;
287 case estFVIR_PREV: block_bc(cr, state->fvir_prev); break;
288 case estNH_XI: nblock_abc(cr, nnht, state->nosehoover_xi); break;
289 case estNH_VXI: nblock_abc(cr, nnht, state->nosehoover_vxi); break;
290 case estNHPRES_XI: nblock_abc(cr, nnhtp, state->nhpres_xi); break;
291 case estNHPRES_VXI: nblock_abc(cr, nnhtp, state->nhpres_vxi); break;
292 case estTC_INT: nblock_abc(cr, state->ngtc, state->therm_integral); break;
293 case estVETA: block_bc(cr, state->veta); break;
294 case estVOL0: block_bc(cr, state->vol0); break;
295 case estX: nblock_abc(cr, state->natoms, state->x); break;
296 case estV: nblock_abc(cr, state->natoms, state->v); break;
297 case estSDX: nblock_abc(cr, state->natoms, state->sd_X); break;
298 case estCGP: nblock_abc(cr, state->natoms, state->cg_p); break;
299 case estLD_RNG: if (state->nrngi == 1)
301 nblock_abc(cr, state->nrng, state->ld_rng);
304 case estLD_RNGI: if (state->nrngi == 1)
306 nblock_abc(cr, state->nrngi, state->ld_rngi);
309 case estDISRE_INITF: block_bc(cr, state->hist.disre_initf); break;
310 case estDISRE_RM3TAV:
311 block_bc(cr, state->hist.ndisrepairs);
312 nblock_abc(cr, state->hist.ndisrepairs, state->hist.disre_rm3tav);
314 case estORIRE_INITF: block_bc(cr, state->hist.orire_initf); break;
316 block_bc(cr, state->hist.norire_Dtav);
317 nblock_abc(cr, state->hist.norire_Dtav, state->hist.orire_Dtav);
321 "Communication is not implemented for %s in bcast_state",
328 static void bc_ilists(const t_commrec *cr, t_ilist *ilist)
332 /* Here we only communicate the non-zero length ilists */
335 for (ftype = 0; ftype < F_NRE; ftype++)
337 if (ilist[ftype].nr > 0)
340 block_bc(cr, ilist[ftype].nr);
341 nblock_bc(cr, ilist[ftype].nr, ilist[ftype].iatoms);
349 for (ftype = 0; ftype < F_NRE; ftype++)
358 block_bc(cr, ilist[ftype].nr);
359 snew_bc(cr, ilist[ftype].iatoms, ilist[ftype].nr);
360 nblock_bc(cr, ilist[ftype].nr, ilist[ftype].iatoms);
368 fprintf(debug, "after bc_ilists\n");
372 static void bc_cmap(const t_commrec *cr, gmx_cmap_t *cmap_grid)
374 int i, j, nelem, ngrid;
376 block_bc(cr, cmap_grid->ngrid);
377 block_bc(cr, cmap_grid->grid_spacing);
379 ngrid = cmap_grid->ngrid;
380 nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
384 snew_bc(cr, cmap_grid->cmapdata, ngrid);
386 for (i = 0; i < ngrid; i++)
388 snew_bc(cr, cmap_grid->cmapdata[i].cmap, 4*nelem);
389 nblock_bc(cr, 4*nelem, cmap_grid->cmapdata[i].cmap);
394 static void bc_ffparams(const t_commrec *cr, gmx_ffparams_t *ffp)
398 block_bc(cr, ffp->ntypes);
399 block_bc(cr, ffp->atnr);
400 snew_bc(cr, ffp->functype, ffp->ntypes);
401 snew_bc(cr, ffp->iparams, ffp->ntypes);
402 nblock_bc(cr, ffp->ntypes, ffp->functype);
403 nblock_bc(cr, ffp->ntypes, ffp->iparams);
404 block_bc(cr, ffp->reppow);
405 block_bc(cr, ffp->fudgeQQ);
406 bc_cmap(cr, &ffp->cmap_grid);
409 static void bc_grpopts(const t_commrec *cr, t_grpopts *g)
413 block_bc(cr, g->ngtc);
414 block_bc(cr, g->ngacc);
415 block_bc(cr, g->ngfrz);
416 block_bc(cr, g->ngener);
417 snew_bc(cr, g->nrdf, g->ngtc);
418 snew_bc(cr, g->tau_t, g->ngtc);
419 snew_bc(cr, g->ref_t, g->ngtc);
420 snew_bc(cr, g->acc, g->ngacc);
421 snew_bc(cr, g->nFreeze, g->ngfrz);
422 snew_bc(cr, g->egp_flags, g->ngener*g->ngener);
424 nblock_bc(cr, g->ngtc, g->nrdf);
425 nblock_bc(cr, g->ngtc, g->tau_t);
426 nblock_bc(cr, g->ngtc, g->ref_t);
427 nblock_bc(cr, g->ngacc, g->acc);
428 nblock_bc(cr, g->ngfrz, g->nFreeze);
429 nblock_bc(cr, g->ngener*g->ngener, g->egp_flags);
430 snew_bc(cr, g->annealing, g->ngtc);
431 snew_bc(cr, g->anneal_npoints, g->ngtc);
432 snew_bc(cr, g->anneal_time, g->ngtc);
433 snew_bc(cr, g->anneal_temp, g->ngtc);
434 nblock_bc(cr, g->ngtc, g->annealing);
435 nblock_bc(cr, g->ngtc, g->anneal_npoints);
436 for (i = 0; (i < g->ngtc); i++)
438 n = g->anneal_npoints[i];
441 snew_bc(cr, g->anneal_time[i], n);
442 snew_bc(cr, g->anneal_temp[i], n);
443 nblock_bc(cr, n, g->anneal_time[i]);
444 nblock_bc(cr, n, g->anneal_temp[i]);
448 /* QMMM stuff, see inputrec */
449 block_bc(cr, g->ngQM);
450 snew_bc(cr, g->QMmethod, g->ngQM);
451 snew_bc(cr, g->QMbasis, g->ngQM);
452 snew_bc(cr, g->QMcharge, g->ngQM);
453 snew_bc(cr, g->QMmult, g->ngQM);
454 snew_bc(cr, g->bSH, g->ngQM);
455 snew_bc(cr, g->CASorbitals, g->ngQM);
456 snew_bc(cr, g->CASelectrons, g->ngQM);
457 snew_bc(cr, g->SAon, g->ngQM);
458 snew_bc(cr, g->SAoff, g->ngQM);
459 snew_bc(cr, g->SAsteps, g->ngQM);
463 nblock_bc(cr, g->ngQM, g->QMmethod);
464 nblock_bc(cr, g->ngQM, g->QMbasis);
465 nblock_bc(cr, g->ngQM, g->QMcharge);
466 nblock_bc(cr, g->ngQM, g->QMmult);
467 nblock_bc(cr, g->ngQM, g->bSH);
468 nblock_bc(cr, g->ngQM, g->CASorbitals);
469 nblock_bc(cr, g->ngQM, g->CASelectrons);
470 nblock_bc(cr, g->ngQM, g->SAon);
471 nblock_bc(cr, g->ngQM, g->SAoff);
472 nblock_bc(cr, g->ngQM, g->SAsteps);
473 /* end of QMMM stuff */
477 static void bc_cosines(const t_commrec *cr, t_cosines *cs)
480 snew_bc(cr, cs->a, cs->n);
481 snew_bc(cr, cs->phi, cs->n);
484 nblock_bc(cr, cs->n, cs->a);
485 nblock_bc(cr, cs->n, cs->phi);
489 static void bc_pull_group(const t_commrec *cr, t_pull_group *pgrp)
494 snew_bc(cr, pgrp->ind, pgrp->nat);
495 nblock_bc(cr, pgrp->nat, pgrp->ind);
497 if (pgrp->nweight > 0)
499 snew_bc(cr, pgrp->weight, pgrp->nweight);
500 nblock_bc(cr, pgrp->nweight, pgrp->weight);
504 static void bc_pull(const t_commrec *cr, t_pull *pull)
509 snew_bc(cr, pull->group, pull->ngroup);
510 for (g = 0; g < pull->ngroup; g++)
512 bc_pull_group(cr, &pull->group[g]);
514 snew_bc(cr, pull->coord, pull->ncoord);
515 nblock_bc(cr, pull->ncoord, pull->coord);
518 static void bc_rotgrp(const t_commrec *cr, t_rotgrp *rotg)
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);
530 static void bc_rot(const t_commrec *cr, t_rot *rot)
535 snew_bc(cr, rot->grp, rot->ngrp);
536 for (g = 0; g < rot->ngrp; g++)
538 bc_rotgrp(cr, &rot->grp[g]);
542 static void bc_adress(const t_commrec *cr, t_adress *adress)
544 block_bc(cr, *adress);
545 if (adress->n_tf_grps > 0)
547 snew_bc(cr, adress->tf_table_index, adress->n_tf_grps);
548 nblock_bc(cr, adress->n_tf_grps, adress->tf_table_index);
550 if (adress->n_energy_grps > 0)
552 snew_bc(cr, adress->group_explicit, adress->n_energy_grps);
553 nblock_bc(cr, adress->n_energy_grps, adress->group_explicit);
556 static void bc_fepvals(const t_commrec *cr, t_lambda *fep)
558 gmx_bool bAlloc = TRUE;
561 block_bc(cr, fep->nstdhdl);
562 block_bc(cr, fep->init_lambda);
563 block_bc(cr, fep->init_fep_state);
564 block_bc(cr, fep->delta_lambda);
565 block_bc(cr, fep->bPrintEnergy);
566 block_bc(cr, fep->n_lambda);
567 if (fep->n_lambda > 0)
569 snew_bc(cr, fep->all_lambda, efptNR);
570 nblock_bc(cr, efptNR, fep->all_lambda);
571 for (i = 0; i < efptNR; i++)
573 snew_bc(cr, fep->all_lambda[i], fep->n_lambda);
574 nblock_bc(cr, fep->n_lambda, fep->all_lambda[i]);
577 block_bc(cr, fep->sc_alpha);
578 block_bc(cr, fep->sc_power);
579 block_bc(cr, fep->sc_r_power);
580 block_bc(cr, fep->sc_sigma);
581 block_bc(cr, fep->sc_sigma_min);
582 block_bc(cr, fep->bScCoul);
583 nblock_bc(cr, efptNR, &(fep->separate_dvdl[0]));
584 block_bc(cr, fep->dhdl_derivatives);
585 block_bc(cr, fep->dh_hist_size);
586 block_bc(cr, fep->dh_hist_spacing);
589 fprintf(debug, "after bc_fepvals\n");
593 static void bc_expandedvals(const t_commrec *cr, t_expanded *expand, int n_lambda)
595 gmx_bool bAlloc = TRUE;
598 block_bc(cr, expand->nstexpanded);
599 block_bc(cr, expand->elamstats);
600 block_bc(cr, expand->elmcmove);
601 block_bc(cr, expand->elmceq);
602 block_bc(cr, expand->equil_n_at_lam);
603 block_bc(cr, expand->equil_wl_delta);
604 block_bc(cr, expand->equil_ratio);
605 block_bc(cr, expand->equil_steps);
606 block_bc(cr, expand->equil_samples);
607 block_bc(cr, expand->lmc_seed);
608 block_bc(cr, expand->minvar);
609 block_bc(cr, expand->minvar_const);
610 block_bc(cr, expand->c_range);
611 block_bc(cr, expand->bSymmetrizedTMatrix);
612 block_bc(cr, expand->nstTij);
613 block_bc(cr, expand->lmc_repeats);
614 block_bc(cr, expand->lmc_forced_nstart);
615 block_bc(cr, expand->gibbsdeltalam);
616 block_bc(cr, expand->wl_scale);
617 block_bc(cr, expand->wl_ratio);
618 block_bc(cr, expand->init_wl_delta);
619 block_bc(cr, expand->bInit_weights);
620 snew_bc(cr, expand->init_lambda_weights, n_lambda);
621 nblock_bc(cr, n_lambda, expand->init_lambda_weights);
622 block_bc(cr, expand->mc_temp);
625 fprintf(debug, "after bc_expandedvals\n");
629 static void bc_simtempvals(const t_commrec *cr, t_simtemp *simtemp, int n_lambda)
631 gmx_bool bAlloc = TRUE;
634 block_bc(cr, simtemp->simtemp_low);
635 block_bc(cr, simtemp->simtemp_high);
636 block_bc(cr, simtemp->eSimTempScale);
637 snew_bc(cr, simtemp->temperatures, n_lambda);
638 nblock_bc(cr, n_lambda, simtemp->temperatures);
641 fprintf(debug, "after bc_simtempvals\n");
645 static void bc_inputrec(const t_commrec *cr, t_inputrec *inputrec)
647 gmx_bool bAlloc = TRUE;
650 block_bc(cr, *inputrec);
652 bc_grpopts(cr, &(inputrec->opts));
654 /* even if efep is efepNO, we need to initialize to make sure that
655 * n_lambda is set to zero */
657 snew_bc(cr, inputrec->fepvals, 1);
658 if (inputrec->efep != efepNO || inputrec->bSimTemp)
660 bc_fepvals(cr, inputrec->fepvals);
662 /* need to initialize this as well because of data checked for in the logic */
663 snew_bc(cr, inputrec->expandedvals, 1);
664 if (inputrec->bExpanded)
666 bc_expandedvals(cr, inputrec->expandedvals, inputrec->fepvals->n_lambda);
668 snew_bc(cr, inputrec->simtempvals, 1);
669 if (inputrec->bSimTemp)
671 bc_simtempvals(cr, inputrec->simtempvals, inputrec->fepvals->n_lambda);
673 if (inputrec->ePull != epullNO)
675 snew_bc(cr, inputrec->pull, 1);
676 bc_pull(cr, inputrec->pull);
680 snew_bc(cr, inputrec->rot, 1);
681 bc_rot(cr, inputrec->rot);
683 for (i = 0; (i < DIM); i++)
685 bc_cosines(cr, &(inputrec->ex[i]));
686 bc_cosines(cr, &(inputrec->et[i]));
688 if (inputrec->bAdress)
690 snew_bc(cr, inputrec->adress, 1);
691 bc_adress(cr, inputrec->adress);
695 static void bc_moltype(const t_commrec *cr, t_symtab *symtab,
696 gmx_moltype_t *moltype)
698 bc_string(cr, symtab, &moltype->name);
699 bc_atoms(cr, symtab, &moltype->atoms);
702 fprintf(debug, "after bc_atoms\n");
705 bc_ilists(cr, moltype->ilist);
706 bc_block(cr, &moltype->cgs);
707 bc_blocka(cr, &moltype->excls);
710 static void bc_molblock(const t_commrec *cr, gmx_molblock_t *molb)
712 gmx_bool bAlloc = TRUE;
714 block_bc(cr, molb->type);
715 block_bc(cr, molb->nmol);
716 block_bc(cr, molb->natoms_mol);
717 block_bc(cr, molb->nposres_xA);
718 if (molb->nposres_xA > 0)
720 snew_bc(cr, molb->posres_xA, molb->nposres_xA);
721 nblock_bc(cr, molb->nposres_xA*DIM, molb->posres_xA[0]);
723 block_bc(cr, molb->nposres_xB);
724 if (molb->nposres_xB > 0)
726 snew_bc(cr, molb->posres_xB, molb->nposres_xB);
727 nblock_bc(cr, molb->nposres_xB*DIM, molb->posres_xB[0]);
731 fprintf(debug, "after bc_molblock\n");
735 static void bc_atomtypes(const t_commrec *cr, t_atomtypes *atomtypes)
739 block_bc(cr, atomtypes->nr);
743 snew_bc(cr, atomtypes->radius, nr);
744 snew_bc(cr, atomtypes->vol, nr);
745 snew_bc(cr, atomtypes->surftens, nr);
746 snew_bc(cr, atomtypes->gb_radius, nr);
747 snew_bc(cr, atomtypes->S_hct, nr);
749 nblock_bc(cr, nr, atomtypes->radius);
750 nblock_bc(cr, nr, atomtypes->vol);
751 nblock_bc(cr, nr, atomtypes->surftens);
752 nblock_bc(cr, nr, atomtypes->gb_radius);
753 nblock_bc(cr, nr, atomtypes->S_hct);
757 void bcast_ir_mtop(const t_commrec *cr, t_inputrec *inputrec, gmx_mtop_t *mtop)
762 fprintf(debug, "in bc_data\n");
764 bc_inputrec(cr, inputrec);
767 fprintf(debug, "after bc_inputrec\n");
769 bc_symtab(cr, &mtop->symtab);
772 fprintf(debug, "after bc_symtab\n");
774 bc_string(cr, &mtop->symtab, &mtop->name);
777 fprintf(debug, "after bc_name\n");
780 bc_ffparams(cr, &mtop->ffparams);
782 block_bc(cr, mtop->nmoltype);
783 snew_bc(cr, mtop->moltype, mtop->nmoltype);
784 for (i = 0; i < mtop->nmoltype; i++)
786 bc_moltype(cr, &mtop->symtab, &mtop->moltype[i]);
789 block_bc(cr, mtop->nmolblock);
790 snew_bc(cr, mtop->molblock, mtop->nmolblock);
791 for (i = 0; i < mtop->nmolblock; i++)
793 bc_molblock(cr, &mtop->molblock[i]);
796 block_bc(cr, mtop->natoms);
798 bc_atomtypes(cr, &mtop->atomtypes);
800 bc_block(cr, &mtop->mols);
801 bc_groups(cr, &mtop->symtab, mtop->natoms, &mtop->groups);