2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, 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.
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.
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.
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.
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.
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.
37 /* This file is completely threadsafe - keep it that way! */
47 #include "types/commrec.h"
49 #include "gromacs/utility/smalloc.h"
50 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/math/vec.h"
55 #define block_bc(cr, d) gmx_bcast( sizeof(d), &(d), (cr))
56 /* Probably the test for (nr) > 0 in the next macro is only needed
57 * on BlueGene(/L), where IBM's MPI_Bcast will segfault after
58 * dereferencing a null pointer, even when no data is to be transferred. */
59 #define nblock_bc(cr, nr, d) { if ((nr) > 0) {gmx_bcast((nr)*sizeof((d)[0]), (d), (cr)); }}
60 #define snew_bc(cr, d, nr) { if (!MASTER(cr)) {snew((d), (nr)); }}
61 /* Dirty macro with bAlloc not as an argument */
62 #define nblock_abc(cr, nr, d) { if (bAlloc) {snew((d), (nr)); } nblock_bc(cr, (nr), (d)); }
64 static void bc_string(const t_commrec *cr, t_symtab *symtab, char ***s)
70 handle = lookup_symtab(symtab, *s);
75 *s = get_symtab_handle(symtab, handle);
79 static void bc_strings(const t_commrec *cr, t_symtab *symtab, int nr, char ****nm)
89 for (i = 0; (i < nr); i++)
91 handle[i] = lookup_symtab(symtab, NM[i]);
94 nblock_bc(cr, nr, handle);
100 for (i = 0; (i < nr); i++)
102 (*nm)[i] = get_symtab_handle(symtab, handle[i]);
108 static void bc_strings_resinfo(const t_commrec *cr, t_symtab *symtab,
109 int nr, t_resinfo *resinfo)
117 for (i = 0; (i < nr); i++)
119 handle[i] = lookup_symtab(symtab, resinfo[i].name);
122 nblock_bc(cr, nr, handle);
126 for (i = 0; (i < nr); i++)
128 resinfo[i].name = get_symtab_handle(symtab, handle[i]);
134 static void bc_symtab(const t_commrec *cr, t_symtab *symtab)
139 block_bc(cr, symtab->nr);
141 snew_bc(cr, symtab->symbuf, 1);
142 symbuf = symtab->symbuf;
143 symbuf->bufsize = nr;
144 snew_bc(cr, symbuf->buf, nr);
145 for (i = 0; i < nr; i++)
149 len = strlen(symbuf->buf[i]) + 1;
152 snew_bc(cr, symbuf->buf[i], len);
153 nblock_bc(cr, len, symbuf->buf[i]);
157 static void bc_block(const t_commrec *cr, t_block *block)
159 block_bc(cr, block->nr);
160 snew_bc(cr, block->index, block->nr+1);
161 nblock_bc(cr, block->nr+1, block->index);
164 static void bc_blocka(const t_commrec *cr, t_blocka *block)
166 block_bc(cr, block->nr);
167 snew_bc(cr, block->index, block->nr+1);
168 nblock_bc(cr, block->nr+1, block->index);
169 block_bc(cr, block->nra);
172 snew_bc(cr, block->a, block->nra);
173 nblock_bc(cr, block->nra, block->a);
177 static void bc_grps(const t_commrec *cr, t_grps grps[])
181 for (i = 0; (i < egcNR); i++)
183 block_bc(cr, grps[i].nr);
184 snew_bc(cr, grps[i].nm_ind, grps[i].nr);
185 nblock_bc(cr, grps[i].nr, grps[i].nm_ind);
189 static void bc_atoms(const t_commrec *cr, t_symtab *symtab, t_atoms *atoms)
193 block_bc(cr, atoms->nr);
194 snew_bc(cr, atoms->atom, atoms->nr);
195 nblock_bc(cr, atoms->nr, atoms->atom);
196 bc_strings(cr, symtab, atoms->nr, &atoms->atomname);
197 block_bc(cr, atoms->nres);
198 snew_bc(cr, atoms->resinfo, atoms->nres);
199 nblock_bc(cr, atoms->nres, atoms->resinfo);
200 bc_strings_resinfo(cr, symtab, atoms->nres, atoms->resinfo);
201 /* QMMM requires atomtypes to be known on all nodes as well */
202 bc_strings(cr, symtab, atoms->nr, &atoms->atomtype);
203 bc_strings(cr, symtab, atoms->nr, &atoms->atomtypeB);
206 static void bc_groups(const t_commrec *cr, t_symtab *symtab,
207 int natoms, gmx_groups_t *groups)
212 bc_grps(cr, groups->grps);
213 block_bc(cr, groups->ngrpname);
214 bc_strings(cr, symtab, groups->ngrpname, &groups->grpname);
215 for (g = 0; g < egcNR; g++)
219 if (groups->grpnr[g])
231 groups->grpnr[g] = NULL;
235 snew_bc(cr, groups->grpnr[g], n);
236 nblock_bc(cr, n, groups->grpnr[g]);
241 fprintf(debug, "after bc_groups\n");
245 void bcast_state(const t_commrec *cr, t_state *state)
255 /* Broadcasts the state sizes and flags from the master to all nodes
256 * in cr->mpi_comm_mygroup. The arrays are not broadcasted. */
257 block_bc(cr, state->natoms);
258 block_bc(cr, state->ngtc);
259 block_bc(cr, state->nnhpres);
260 block_bc(cr, state->nhchainlength);
261 block_bc(cr, state->flags);
262 if (state->lambda == NULL)
264 snew_bc(cr, state->lambda, efptNR)
269 /* We allocate dynamically in dd_partition_system. */
272 /* The code below is reachable only by TPI and NM, so it is not
273 tested by anything. */
275 nnht = (state->ngtc)*(state->nhchainlength);
276 nnhtp = (state->nnhpres)*(state->nhchainlength);
278 /* We still need to allocate the arrays in state for non-master
279 * ranks, which is done (implicitly via bAlloc) in the dirty,
280 * dirty nblock_abc macro. */
281 bAlloc = !MASTER(cr);
284 state->nalloc = state->natoms;
286 for (i = 0; i < estNR; i++)
288 if (state->flags & (1<<i))
292 case estLAMBDA: nblock_bc(cr, efptNR, state->lambda); break;
293 case estFEPSTATE: block_bc(cr, state->fep_state); break;
294 case estBOX: block_bc(cr, state->box); break;
295 case estBOX_REL: block_bc(cr, state->box_rel); break;
296 case estBOXV: block_bc(cr, state->boxv); break;
297 case estPRES_PREV: block_bc(cr, state->pres_prev); break;
298 case estSVIR_PREV: block_bc(cr, state->svir_prev); break;
299 case estFVIR_PREV: block_bc(cr, state->fvir_prev); break;
300 case estNH_XI: nblock_abc(cr, nnht, state->nosehoover_xi); break;
301 case estNH_VXI: nblock_abc(cr, nnht, state->nosehoover_vxi); break;
302 case estNHPRES_XI: nblock_abc(cr, nnhtp, state->nhpres_xi); break;
303 case estNHPRES_VXI: nblock_abc(cr, nnhtp, state->nhpres_vxi); break;
304 case estTC_INT: nblock_abc(cr, state->ngtc, state->therm_integral); break;
305 case estVETA: block_bc(cr, state->veta); break;
306 case estVOL0: block_bc(cr, state->vol0); break;
307 case estX: nblock_abc(cr, state->natoms, state->x); break;
308 case estV: nblock_abc(cr, state->natoms, state->v); break;
309 case estSDX: nblock_abc(cr, state->natoms, state->sd_X); break;
310 case estCGP: nblock_abc(cr, state->natoms, state->cg_p); break;
311 case estDISRE_INITF: block_bc(cr, state->hist.disre_initf); break;
312 case estDISRE_RM3TAV:
313 block_bc(cr, state->hist.ndisrepairs);
314 nblock_abc(cr, state->hist.ndisrepairs, state->hist.disre_rm3tav);
316 case estORIRE_INITF: block_bc(cr, state->hist.orire_initf); break;
318 block_bc(cr, state->hist.norire_Dtav);
319 nblock_abc(cr, state->hist.norire_Dtav, state->hist.orire_Dtav);
323 "Communication is not implemented for %s in bcast_state",
330 static void bc_ilists(const t_commrec *cr, t_ilist *ilist)
334 /* Here we only communicate the non-zero length ilists */
337 for (ftype = 0; ftype < F_NRE; ftype++)
339 if (ilist[ftype].nr > 0)
342 block_bc(cr, ilist[ftype].nr);
343 nblock_bc(cr, ilist[ftype].nr, ilist[ftype].iatoms);
351 for (ftype = 0; ftype < F_NRE; ftype++)
360 block_bc(cr, ilist[ftype].nr);
361 snew_bc(cr, ilist[ftype].iatoms, ilist[ftype].nr);
362 nblock_bc(cr, ilist[ftype].nr, ilist[ftype].iatoms);
370 fprintf(debug, "after bc_ilists\n");
374 static void bc_cmap(const t_commrec *cr, gmx_cmap_t *cmap_grid)
376 int i, j, nelem, ngrid;
378 block_bc(cr, cmap_grid->ngrid);
379 block_bc(cr, cmap_grid->grid_spacing);
381 ngrid = cmap_grid->ngrid;
382 nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
386 snew_bc(cr, cmap_grid->cmapdata, ngrid);
388 for (i = 0; i < ngrid; i++)
390 snew_bc(cr, cmap_grid->cmapdata[i].cmap, 4*nelem);
391 nblock_bc(cr, 4*nelem, cmap_grid->cmapdata[i].cmap);
396 static void bc_ffparams(const t_commrec *cr, gmx_ffparams_t *ffp)
400 block_bc(cr, ffp->ntypes);
401 block_bc(cr, ffp->atnr);
402 snew_bc(cr, ffp->functype, ffp->ntypes);
403 snew_bc(cr, ffp->iparams, ffp->ntypes);
404 nblock_bc(cr, ffp->ntypes, ffp->functype);
405 nblock_bc(cr, ffp->ntypes, ffp->iparams);
406 block_bc(cr, ffp->reppow);
407 block_bc(cr, ffp->fudgeQQ);
408 bc_cmap(cr, &ffp->cmap_grid);
411 static void bc_grpopts(const t_commrec *cr, t_grpopts *g)
415 block_bc(cr, g->ngtc);
416 block_bc(cr, g->ngacc);
417 block_bc(cr, g->ngfrz);
418 block_bc(cr, g->ngener);
419 snew_bc(cr, g->nrdf, g->ngtc);
420 snew_bc(cr, g->tau_t, g->ngtc);
421 snew_bc(cr, g->ref_t, g->ngtc);
422 snew_bc(cr, g->acc, g->ngacc);
423 snew_bc(cr, g->nFreeze, g->ngfrz);
424 snew_bc(cr, g->egp_flags, g->ngener*g->ngener);
426 nblock_bc(cr, g->ngtc, g->nrdf);
427 nblock_bc(cr, g->ngtc, g->tau_t);
428 nblock_bc(cr, g->ngtc, g->ref_t);
429 nblock_bc(cr, g->ngacc, g->acc);
430 nblock_bc(cr, g->ngfrz, g->nFreeze);
431 nblock_bc(cr, g->ngener*g->ngener, g->egp_flags);
432 snew_bc(cr, g->annealing, g->ngtc);
433 snew_bc(cr, g->anneal_npoints, g->ngtc);
434 snew_bc(cr, g->anneal_time, g->ngtc);
435 snew_bc(cr, g->anneal_temp, g->ngtc);
436 nblock_bc(cr, g->ngtc, g->annealing);
437 nblock_bc(cr, g->ngtc, g->anneal_npoints);
438 for (i = 0; (i < g->ngtc); i++)
440 n = g->anneal_npoints[i];
443 snew_bc(cr, g->anneal_time[i], n);
444 snew_bc(cr, g->anneal_temp[i], n);
445 nblock_bc(cr, n, g->anneal_time[i]);
446 nblock_bc(cr, n, g->anneal_temp[i]);
450 /* QMMM stuff, see inputrec */
451 block_bc(cr, g->ngQM);
452 snew_bc(cr, g->QMmethod, g->ngQM);
453 snew_bc(cr, g->QMbasis, g->ngQM);
454 snew_bc(cr, g->QMcharge, g->ngQM);
455 snew_bc(cr, g->QMmult, g->ngQM);
456 snew_bc(cr, g->bSH, g->ngQM);
457 snew_bc(cr, g->CASorbitals, g->ngQM);
458 snew_bc(cr, g->CASelectrons, g->ngQM);
459 snew_bc(cr, g->SAon, g->ngQM);
460 snew_bc(cr, g->SAoff, g->ngQM);
461 snew_bc(cr, g->SAsteps, g->ngQM);
465 nblock_bc(cr, g->ngQM, g->QMmethod);
466 nblock_bc(cr, g->ngQM, g->QMbasis);
467 nblock_bc(cr, g->ngQM, g->QMcharge);
468 nblock_bc(cr, g->ngQM, g->QMmult);
469 nblock_bc(cr, g->ngQM, g->bSH);
470 nblock_bc(cr, g->ngQM, g->CASorbitals);
471 nblock_bc(cr, g->ngQM, g->CASelectrons);
472 nblock_bc(cr, g->ngQM, g->SAon);
473 nblock_bc(cr, g->ngQM, g->SAoff);
474 nblock_bc(cr, g->ngQM, g->SAsteps);
475 /* end of QMMM stuff */
479 static void bc_cosines(const t_commrec *cr, t_cosines *cs)
482 snew_bc(cr, cs->a, cs->n);
483 snew_bc(cr, cs->phi, cs->n);
486 nblock_bc(cr, cs->n, cs->a);
487 nblock_bc(cr, cs->n, cs->phi);
491 static void bc_pull_group(const t_commrec *cr, t_pull_group *pgrp)
496 snew_bc(cr, pgrp->ind, pgrp->nat);
497 nblock_bc(cr, pgrp->nat, pgrp->ind);
499 if (pgrp->nweight > 0)
501 snew_bc(cr, pgrp->weight, pgrp->nweight);
502 nblock_bc(cr, pgrp->nweight, pgrp->weight);
506 static void bc_pull(const t_commrec *cr, t_pull *pull)
511 snew_bc(cr, pull->group, pull->ngroup);
512 for (g = 0; g < pull->ngroup; g++)
514 bc_pull_group(cr, &pull->group[g]);
516 snew_bc(cr, pull->coord, pull->ncoord);
517 nblock_bc(cr, pull->ncoord, pull->coord);
520 static void bc_rotgrp(const t_commrec *cr, t_rotgrp *rotg)
525 snew_bc(cr, rotg->ind, rotg->nat);
526 nblock_bc(cr, rotg->nat, rotg->ind);
527 snew_bc(cr, rotg->x_ref, rotg->nat);
528 nblock_bc(cr, rotg->nat, rotg->x_ref);
532 static void bc_rot(const t_commrec *cr, t_rot *rot)
537 snew_bc(cr, rot->grp, rot->ngrp);
538 for (g = 0; g < rot->ngrp; g++)
540 bc_rotgrp(cr, &rot->grp[g]);
544 static void bc_adress(const t_commrec *cr, t_adress *adress)
546 block_bc(cr, *adress);
547 if (adress->n_tf_grps > 0)
549 snew_bc(cr, adress->tf_table_index, adress->n_tf_grps);
550 nblock_bc(cr, adress->n_tf_grps, adress->tf_table_index);
552 if (adress->n_energy_grps > 0)
554 snew_bc(cr, adress->group_explicit, adress->n_energy_grps);
555 nblock_bc(cr, adress->n_energy_grps, adress->group_explicit);
559 static void bc_imd(const t_commrec *cr, t_IMD *imd)
564 snew_bc(cr, imd->ind, imd->nat);
565 nblock_bc(cr, imd->nat, imd->ind);
568 static void bc_fepvals(const t_commrec *cr, t_lambda *fep)
570 gmx_bool bAlloc = TRUE;
573 block_bc(cr, fep->nstdhdl);
574 block_bc(cr, fep->init_lambda);
575 block_bc(cr, fep->init_fep_state);
576 block_bc(cr, fep->delta_lambda);
577 block_bc(cr, fep->bPrintEnergy);
578 block_bc(cr, fep->n_lambda);
579 if (fep->n_lambda > 0)
581 snew_bc(cr, fep->all_lambda, efptNR);
582 nblock_bc(cr, efptNR, fep->all_lambda);
583 for (i = 0; i < efptNR; i++)
585 snew_bc(cr, fep->all_lambda[i], fep->n_lambda);
586 nblock_bc(cr, fep->n_lambda, fep->all_lambda[i]);
589 block_bc(cr, fep->sc_alpha);
590 block_bc(cr, fep->sc_power);
591 block_bc(cr, fep->sc_r_power);
592 block_bc(cr, fep->sc_sigma);
593 block_bc(cr, fep->sc_sigma_min);
594 block_bc(cr, fep->bScCoul);
595 nblock_bc(cr, efptNR, &(fep->separate_dvdl[0]));
596 block_bc(cr, fep->dhdl_derivatives);
597 block_bc(cr, fep->dh_hist_size);
598 block_bc(cr, fep->dh_hist_spacing);
601 fprintf(debug, "after bc_fepvals\n");
605 static void bc_expandedvals(const t_commrec *cr, t_expanded *expand, int n_lambda)
607 gmx_bool bAlloc = TRUE;
610 block_bc(cr, expand->nstexpanded);
611 block_bc(cr, expand->elamstats);
612 block_bc(cr, expand->elmcmove);
613 block_bc(cr, expand->elmceq);
614 block_bc(cr, expand->equil_n_at_lam);
615 block_bc(cr, expand->equil_wl_delta);
616 block_bc(cr, expand->equil_ratio);
617 block_bc(cr, expand->equil_steps);
618 block_bc(cr, expand->equil_samples);
619 block_bc(cr, expand->lmc_seed);
620 block_bc(cr, expand->minvar);
621 block_bc(cr, expand->minvar_const);
622 block_bc(cr, expand->c_range);
623 block_bc(cr, expand->bSymmetrizedTMatrix);
624 block_bc(cr, expand->nstTij);
625 block_bc(cr, expand->lmc_repeats);
626 block_bc(cr, expand->lmc_forced_nstart);
627 block_bc(cr, expand->gibbsdeltalam);
628 block_bc(cr, expand->wl_scale);
629 block_bc(cr, expand->wl_ratio);
630 block_bc(cr, expand->init_wl_delta);
631 block_bc(cr, expand->bInit_weights);
632 snew_bc(cr, expand->init_lambda_weights, n_lambda);
633 nblock_bc(cr, n_lambda, expand->init_lambda_weights);
634 block_bc(cr, expand->mc_temp);
637 fprintf(debug, "after bc_expandedvals\n");
641 static void bc_simtempvals(const t_commrec *cr, t_simtemp *simtemp, int n_lambda)
643 gmx_bool bAlloc = TRUE;
646 block_bc(cr, simtemp->simtemp_low);
647 block_bc(cr, simtemp->simtemp_high);
648 block_bc(cr, simtemp->eSimTempScale);
649 snew_bc(cr, simtemp->temperatures, n_lambda);
650 nblock_bc(cr, n_lambda, simtemp->temperatures);
653 fprintf(debug, "after bc_simtempvals\n");
658 static void bc_swapions(const t_commrec *cr, t_swapcoords *swap)
665 /* Broadcast ion group atom indices */
666 snew_bc(cr, swap->ind, swap->nat);
667 nblock_bc(cr, swap->nat, swap->ind);
669 /* Broadcast split groups atom indices */
670 for (i = 0; i < 2; i++)
672 snew_bc(cr, swap->ind_split[i], swap->nat_split[i]);
673 nblock_bc(cr, swap->nat_split[i], swap->ind_split[i]);
676 /* Broadcast solvent group atom indices */
677 snew_bc(cr, swap->ind_sol, swap->nat_sol);
678 nblock_bc(cr, swap->nat_sol, swap->ind_sol);
682 static void bc_inputrec(const t_commrec *cr, t_inputrec *inputrec)
684 gmx_bool bAlloc = TRUE;
687 block_bc(cr, *inputrec);
689 bc_grpopts(cr, &(inputrec->opts));
691 /* even if efep is efepNO, we need to initialize to make sure that
692 * n_lambda is set to zero */
694 snew_bc(cr, inputrec->fepvals, 1);
695 if (inputrec->efep != efepNO || inputrec->bSimTemp)
697 bc_fepvals(cr, inputrec->fepvals);
699 /* need to initialize this as well because of data checked for in the logic */
700 snew_bc(cr, inputrec->expandedvals, 1);
701 if (inputrec->bExpanded)
703 bc_expandedvals(cr, inputrec->expandedvals, inputrec->fepvals->n_lambda);
705 snew_bc(cr, inputrec->simtempvals, 1);
706 if (inputrec->bSimTemp)
708 bc_simtempvals(cr, inputrec->simtempvals, inputrec->fepvals->n_lambda);
710 if (inputrec->ePull != epullNO)
712 snew_bc(cr, inputrec->pull, 1);
713 bc_pull(cr, inputrec->pull);
717 snew_bc(cr, inputrec->rot, 1);
718 bc_rot(cr, inputrec->rot);
722 snew_bc(cr, inputrec->imd, 1);
723 bc_imd(cr, inputrec->imd);
725 for (i = 0; (i < DIM); i++)
727 bc_cosines(cr, &(inputrec->ex[i]));
728 bc_cosines(cr, &(inputrec->et[i]));
730 if (inputrec->eSwapCoords != eswapNO)
732 snew_bc(cr, inputrec->swap, 1);
733 bc_swapions(cr, inputrec->swap);
735 if (inputrec->bAdress)
737 snew_bc(cr, inputrec->adress, 1);
738 bc_adress(cr, inputrec->adress);
742 static void bc_moltype(const t_commrec *cr, t_symtab *symtab,
743 gmx_moltype_t *moltype)
745 bc_string(cr, symtab, &moltype->name);
746 bc_atoms(cr, symtab, &moltype->atoms);
749 fprintf(debug, "after bc_atoms\n");
752 bc_ilists(cr, moltype->ilist);
753 bc_block(cr, &moltype->cgs);
754 bc_blocka(cr, &moltype->excls);
757 static void bc_molblock(const t_commrec *cr, gmx_molblock_t *molb)
759 gmx_bool bAlloc = TRUE;
761 block_bc(cr, molb->type);
762 block_bc(cr, molb->nmol);
763 block_bc(cr, molb->natoms_mol);
764 block_bc(cr, molb->nposres_xA);
765 if (molb->nposres_xA > 0)
767 snew_bc(cr, molb->posres_xA, molb->nposres_xA);
768 nblock_bc(cr, molb->nposres_xA*DIM, molb->posres_xA[0]);
770 block_bc(cr, molb->nposres_xB);
771 if (molb->nposres_xB > 0)
773 snew_bc(cr, molb->posres_xB, molb->nposres_xB);
774 nblock_bc(cr, molb->nposres_xB*DIM, molb->posres_xB[0]);
778 fprintf(debug, "after bc_molblock\n");
782 static void bc_atomtypes(const t_commrec *cr, t_atomtypes *atomtypes)
786 block_bc(cr, atomtypes->nr);
790 snew_bc(cr, atomtypes->radius, nr);
791 snew_bc(cr, atomtypes->vol, nr);
792 snew_bc(cr, atomtypes->surftens, nr);
793 snew_bc(cr, atomtypes->gb_radius, nr);
794 snew_bc(cr, atomtypes->S_hct, nr);
796 nblock_bc(cr, nr, atomtypes->radius);
797 nblock_bc(cr, nr, atomtypes->vol);
798 nblock_bc(cr, nr, atomtypes->surftens);
799 nblock_bc(cr, nr, atomtypes->gb_radius);
800 nblock_bc(cr, nr, atomtypes->S_hct);
804 void bcast_ir_mtop(const t_commrec *cr, t_inputrec *inputrec, gmx_mtop_t *mtop)
809 fprintf(debug, "in bc_data\n");
811 bc_inputrec(cr, inputrec);
814 fprintf(debug, "after bc_inputrec\n");
816 bc_symtab(cr, &mtop->symtab);
819 fprintf(debug, "after bc_symtab\n");
821 bc_string(cr, &mtop->symtab, &mtop->name);
824 fprintf(debug, "after bc_name\n");
827 bc_ffparams(cr, &mtop->ffparams);
829 block_bc(cr, mtop->nmoltype);
830 snew_bc(cr, mtop->moltype, mtop->nmoltype);
831 for (i = 0; i < mtop->nmoltype; i++)
833 bc_moltype(cr, &mtop->symtab, &mtop->moltype[i]);
836 block_bc(cr, mtop->nmolblock);
837 snew_bc(cr, mtop->molblock, mtop->nmolblock);
838 for (i = 0; i < mtop->nmolblock; i++)
840 bc_molblock(cr, &mtop->molblock[i]);
843 block_bc(cr, mtop->natoms);
845 bc_atomtypes(cr, &mtop->atomtypes);
847 bc_block(cr, &mtop->mols);
848 bc_groups(cr, &mtop->symtab, mtop->natoms, &mtop->groups);