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,2015, 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! */
42 #include "gromacs/legacyheaders/main.h"
43 #include "gromacs/legacyheaders/mdrun.h"
44 #include "gromacs/legacyheaders/network.h"
45 #include "gromacs/legacyheaders/tgroup.h"
46 #include "gromacs/legacyheaders/typedefs.h"
47 #include "gromacs/legacyheaders/types/commrec.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/topology/symtab.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/smalloc.h"
53 #define block_bc(cr, d) gmx_bcast( sizeof(d), &(d), (cr))
54 /* Probably the test for (nr) > 0 in the next macro is only needed
55 * on BlueGene(/L), where IBM's MPI_Bcast will segfault after
56 * dereferencing a null pointer, even when no data is to be transferred. */
57 #define nblock_bc(cr, nr, d) { if ((nr) > 0) {gmx_bcast((nr)*sizeof((d)[0]), (d), (cr)); }}
58 #define snew_bc(cr, d, nr) { if (!MASTER(cr)) {snew((d), (nr)); }}
59 /* Dirty macro with bAlloc not as an argument */
60 #define nblock_abc(cr, nr, d) { if (bAlloc) {snew((d), (nr)); } nblock_bc(cr, (nr), (d)); }
62 static void bc_string(const t_commrec *cr, t_symtab *symtab, char ***s)
68 handle = lookup_symtab(symtab, *s);
73 *s = get_symtab_handle(symtab, handle);
77 static void bc_strings(const t_commrec *cr, t_symtab *symtab, int nr, char ****nm)
85 for (i = 0; (i < nr); i++)
87 handle[i] = lookup_symtab(symtab, (*nm)[i]);
90 nblock_bc(cr, nr, handle);
95 for (i = 0; (i < nr); i++)
97 (*nm)[i] = get_symtab_handle(symtab, handle[i]);
103 static void bc_strings_resinfo(const t_commrec *cr, t_symtab *symtab,
104 int nr, t_resinfo *resinfo)
112 for (i = 0; (i < nr); i++)
114 handle[i] = lookup_symtab(symtab, resinfo[i].name);
117 nblock_bc(cr, nr, handle);
121 for (i = 0; (i < nr); i++)
123 resinfo[i].name = get_symtab_handle(symtab, handle[i]);
129 static void bc_symtab(const t_commrec *cr, t_symtab *symtab)
134 block_bc(cr, symtab->nr);
136 snew_bc(cr, symtab->symbuf, 1);
137 symbuf = symtab->symbuf;
138 symbuf->bufsize = nr;
139 snew_bc(cr, symbuf->buf, nr);
140 for (i = 0; i < nr; i++)
144 len = strlen(symbuf->buf[i]) + 1;
147 snew_bc(cr, symbuf->buf[i], len);
148 nblock_bc(cr, len, symbuf->buf[i]);
152 static void bc_block(const t_commrec *cr, t_block *block)
154 block_bc(cr, block->nr);
155 snew_bc(cr, block->index, block->nr+1);
156 nblock_bc(cr, block->nr+1, block->index);
159 static void bc_blocka(const t_commrec *cr, t_blocka *block)
161 block_bc(cr, block->nr);
162 snew_bc(cr, block->index, block->nr+1);
163 nblock_bc(cr, block->nr+1, block->index);
164 block_bc(cr, block->nra);
167 snew_bc(cr, block->a, block->nra);
168 nblock_bc(cr, block->nra, block->a);
172 static void bc_grps(const t_commrec *cr, t_grps grps[])
176 for (i = 0; (i < egcNR); i++)
178 block_bc(cr, grps[i].nr);
179 snew_bc(cr, grps[i].nm_ind, grps[i].nr);
180 nblock_bc(cr, grps[i].nr, grps[i].nm_ind);
184 static void bc_atoms(const t_commrec *cr, t_symtab *symtab, t_atoms *atoms)
186 block_bc(cr, atoms->nr);
187 snew_bc(cr, atoms->atom, atoms->nr);
188 nblock_bc(cr, atoms->nr, atoms->atom);
189 bc_strings(cr, symtab, atoms->nr, &atoms->atomname);
190 block_bc(cr, atoms->nres);
191 snew_bc(cr, atoms->resinfo, atoms->nres);
192 nblock_bc(cr, atoms->nres, atoms->resinfo);
193 bc_strings_resinfo(cr, symtab, atoms->nres, atoms->resinfo);
194 /* QMMM requires atomtypes to be known on all nodes as well */
195 bc_strings(cr, symtab, atoms->nr, &atoms->atomtype);
196 bc_strings(cr, symtab, atoms->nr, &atoms->atomtypeB);
199 static void bc_groups(const t_commrec *cr, t_symtab *symtab,
200 int natoms, gmx_groups_t *groups)
204 bc_grps(cr, groups->grps);
205 block_bc(cr, groups->ngrpname);
206 bc_strings(cr, symtab, groups->ngrpname, &groups->grpname);
207 for (g = 0; g < egcNR; g++)
211 if (groups->grpnr[g])
223 groups->grpnr[g] = NULL;
227 snew_bc(cr, groups->grpnr[g], n);
228 nblock_bc(cr, n, groups->grpnr[g]);
233 fprintf(debug, "after bc_groups\n");
237 void bcast_state(const t_commrec *cr, t_state *state)
247 /* Broadcasts the state sizes and flags from the master to all nodes
248 * in cr->mpi_comm_mygroup. The arrays are not broadcasted. */
249 block_bc(cr, state->natoms);
250 block_bc(cr, state->ngtc);
251 block_bc(cr, state->nnhpres);
252 block_bc(cr, state->nhchainlength);
253 block_bc(cr, state->flags);
254 if (state->lambda == NULL)
256 snew_bc(cr, state->lambda, efptNR)
261 /* We allocate dynamically in dd_partition_system. */
264 /* The code below is reachable only by TPI and NM, so it is not
265 tested by anything. */
267 nnht = (state->ngtc)*(state->nhchainlength);
268 nnhtp = (state->nnhpres)*(state->nhchainlength);
270 /* We still need to allocate the arrays in state for non-master
271 * ranks, which is done (implicitly via bAlloc) in the dirty,
272 * dirty nblock_abc macro. */
273 bAlloc = !MASTER(cr);
276 state->nalloc = state->natoms;
278 for (i = 0; i < estNR; i++)
280 if (state->flags & (1<<i))
284 case estLAMBDA: nblock_bc(cr, efptNR, state->lambda); break;
285 case estFEPSTATE: block_bc(cr, state->fep_state); break;
286 case estBOX: block_bc(cr, state->box); break;
287 case estBOX_REL: block_bc(cr, state->box_rel); break;
288 case estBOXV: block_bc(cr, state->boxv); break;
289 case estPRES_PREV: block_bc(cr, state->pres_prev); break;
290 case estSVIR_PREV: block_bc(cr, state->svir_prev); break;
291 case estFVIR_PREV: block_bc(cr, state->fvir_prev); break;
292 case estNH_XI: nblock_abc(cr, nnht, state->nosehoover_xi); break;
293 case estNH_VXI: nblock_abc(cr, nnht, state->nosehoover_vxi); break;
294 case estNHPRES_XI: nblock_abc(cr, nnhtp, state->nhpres_xi); break;
295 case estNHPRES_VXI: nblock_abc(cr, nnhtp, state->nhpres_vxi); break;
296 case estTC_INT: nblock_abc(cr, state->ngtc, state->therm_integral); break;
297 case estVETA: block_bc(cr, state->veta); break;
298 case estVOL0: block_bc(cr, state->vol0); break;
299 case estX: nblock_abc(cr, state->natoms, state->x); break;
300 case estV: nblock_abc(cr, state->natoms, state->v); break;
301 case estSDX: nblock_abc(cr, state->natoms, state->sd_X); break;
302 case estCGP: nblock_abc(cr, state->natoms, state->cg_p); break;
303 case estDISRE_INITF: block_bc(cr, state->hist.disre_initf); break;
304 case estDISRE_RM3TAV:
305 block_bc(cr, state->hist.ndisrepairs);
306 nblock_abc(cr, state->hist.ndisrepairs, state->hist.disre_rm3tav);
308 case estORIRE_INITF: block_bc(cr, state->hist.orire_initf); break;
310 block_bc(cr, state->hist.norire_Dtav);
311 nblock_abc(cr, state->hist.norire_Dtav, state->hist.orire_Dtav);
315 "Communication is not implemented for %s in bcast_state",
322 static void bc_ilists(const t_commrec *cr, t_ilist *ilist)
326 /* Here we only communicate the non-zero length ilists */
329 for (ftype = 0; ftype < F_NRE; ftype++)
331 if (ilist[ftype].nr > 0)
334 block_bc(cr, ilist[ftype].nr);
335 nblock_bc(cr, ilist[ftype].nr, ilist[ftype].iatoms);
343 for (ftype = 0; ftype < F_NRE; ftype++)
352 block_bc(cr, ilist[ftype].nr);
353 snew_bc(cr, ilist[ftype].iatoms, ilist[ftype].nr);
354 nblock_bc(cr, ilist[ftype].nr, ilist[ftype].iatoms);
362 fprintf(debug, "after bc_ilists\n");
366 static void bc_cmap(const t_commrec *cr, gmx_cmap_t *cmap_grid)
370 block_bc(cr, cmap_grid->ngrid);
371 block_bc(cr, cmap_grid->grid_spacing);
373 ngrid = cmap_grid->ngrid;
374 nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
378 snew_bc(cr, cmap_grid->cmapdata, ngrid);
380 for (i = 0; i < ngrid; i++)
382 snew_bc(cr, cmap_grid->cmapdata[i].cmap, 4*nelem);
383 nblock_bc(cr, 4*nelem, cmap_grid->cmapdata[i].cmap);
388 static void bc_ffparams(const t_commrec *cr, gmx_ffparams_t *ffp)
390 block_bc(cr, ffp->ntypes);
391 block_bc(cr, ffp->atnr);
392 snew_bc(cr, ffp->functype, ffp->ntypes);
393 snew_bc(cr, ffp->iparams, ffp->ntypes);
394 nblock_bc(cr, ffp->ntypes, ffp->functype);
395 nblock_bc(cr, ffp->ntypes, ffp->iparams);
396 block_bc(cr, ffp->reppow);
397 block_bc(cr, ffp->fudgeQQ);
398 bc_cmap(cr, &ffp->cmap_grid);
401 static void bc_grpopts(const t_commrec *cr, t_grpopts *g)
405 block_bc(cr, g->ngtc);
406 block_bc(cr, g->ngacc);
407 block_bc(cr, g->ngfrz);
408 block_bc(cr, g->ngener);
409 snew_bc(cr, g->nrdf, g->ngtc);
410 snew_bc(cr, g->tau_t, g->ngtc);
411 snew_bc(cr, g->ref_t, g->ngtc);
412 snew_bc(cr, g->acc, g->ngacc);
413 snew_bc(cr, g->nFreeze, g->ngfrz);
414 snew_bc(cr, g->egp_flags, g->ngener*g->ngener);
416 nblock_bc(cr, g->ngtc, g->nrdf);
417 nblock_bc(cr, g->ngtc, g->tau_t);
418 nblock_bc(cr, g->ngtc, g->ref_t);
419 nblock_bc(cr, g->ngacc, g->acc);
420 nblock_bc(cr, g->ngfrz, g->nFreeze);
421 nblock_bc(cr, g->ngener*g->ngener, g->egp_flags);
422 snew_bc(cr, g->annealing, g->ngtc);
423 snew_bc(cr, g->anneal_npoints, g->ngtc);
424 snew_bc(cr, g->anneal_time, g->ngtc);
425 snew_bc(cr, g->anneal_temp, g->ngtc);
426 nblock_bc(cr, g->ngtc, g->annealing);
427 nblock_bc(cr, g->ngtc, g->anneal_npoints);
428 for (i = 0; (i < g->ngtc); i++)
430 n = g->anneal_npoints[i];
433 snew_bc(cr, g->anneal_time[i], n);
434 snew_bc(cr, g->anneal_temp[i], n);
435 nblock_bc(cr, n, g->anneal_time[i]);
436 nblock_bc(cr, n, g->anneal_temp[i]);
440 /* QMMM stuff, see inputrec */
441 block_bc(cr, g->ngQM);
442 snew_bc(cr, g->QMmethod, g->ngQM);
443 snew_bc(cr, g->QMbasis, g->ngQM);
444 snew_bc(cr, g->QMcharge, g->ngQM);
445 snew_bc(cr, g->QMmult, g->ngQM);
446 snew_bc(cr, g->bSH, g->ngQM);
447 snew_bc(cr, g->CASorbitals, g->ngQM);
448 snew_bc(cr, g->CASelectrons, g->ngQM);
449 snew_bc(cr, g->SAon, g->ngQM);
450 snew_bc(cr, g->SAoff, g->ngQM);
451 snew_bc(cr, g->SAsteps, g->ngQM);
455 nblock_bc(cr, g->ngQM, g->QMmethod);
456 nblock_bc(cr, g->ngQM, g->QMbasis);
457 nblock_bc(cr, g->ngQM, g->QMcharge);
458 nblock_bc(cr, g->ngQM, g->QMmult);
459 nblock_bc(cr, g->ngQM, g->bSH);
460 nblock_bc(cr, g->ngQM, g->CASorbitals);
461 nblock_bc(cr, g->ngQM, g->CASelectrons);
462 nblock_bc(cr, g->ngQM, g->SAon);
463 nblock_bc(cr, g->ngQM, g->SAoff);
464 nblock_bc(cr, g->ngQM, g->SAsteps);
465 /* end of QMMM stuff */
469 static void bc_cosines(const t_commrec *cr, t_cosines *cs)
472 snew_bc(cr, cs->a, cs->n);
473 snew_bc(cr, cs->phi, cs->n);
476 nblock_bc(cr, cs->n, cs->a);
477 nblock_bc(cr, cs->n, cs->phi);
481 static void bc_pull_group(const t_commrec *cr, t_pull_group *pgrp)
486 snew_bc(cr, pgrp->ind, pgrp->nat);
487 nblock_bc(cr, pgrp->nat, pgrp->ind);
489 if (pgrp->nweight > 0)
491 snew_bc(cr, pgrp->weight, pgrp->nweight);
492 nblock_bc(cr, pgrp->nweight, pgrp->weight);
496 static void bc_pull(const t_commrec *cr, t_pull *pull)
501 snew_bc(cr, pull->group, pull->ngroup);
502 for (g = 0; g < pull->ngroup; g++)
504 bc_pull_group(cr, &pull->group[g]);
506 snew_bc(cr, pull->coord, pull->ncoord);
507 nblock_bc(cr, pull->ncoord, pull->coord);
510 static void bc_rotgrp(const t_commrec *cr, t_rotgrp *rotg)
515 snew_bc(cr, rotg->ind, rotg->nat);
516 nblock_bc(cr, rotg->nat, rotg->ind);
517 snew_bc(cr, rotg->x_ref, rotg->nat);
518 nblock_bc(cr, rotg->nat, rotg->x_ref);
522 static void bc_rot(const t_commrec *cr, t_rot *rot)
527 snew_bc(cr, rot->grp, rot->ngrp);
528 for (g = 0; g < rot->ngrp; g++)
530 bc_rotgrp(cr, &rot->grp[g]);
534 static void bc_adress(const t_commrec *cr, t_adress *adress)
536 block_bc(cr, *adress);
537 if (adress->n_tf_grps > 0)
539 snew_bc(cr, adress->tf_table_index, adress->n_tf_grps);
540 nblock_bc(cr, adress->n_tf_grps, adress->tf_table_index);
542 if (adress->n_energy_grps > 0)
544 snew_bc(cr, adress->group_explicit, adress->n_energy_grps);
545 nblock_bc(cr, adress->n_energy_grps, adress->group_explicit);
549 static void bc_imd(const t_commrec *cr, t_IMD *imd)
552 snew_bc(cr, imd->ind, imd->nat);
553 nblock_bc(cr, imd->nat, imd->ind);
556 static void bc_fepvals(const t_commrec *cr, t_lambda *fep)
560 block_bc(cr, fep->nstdhdl);
561 block_bc(cr, fep->init_lambda);
562 block_bc(cr, fep->init_fep_state);
563 block_bc(cr, fep->delta_lambda);
564 block_bc(cr, fep->edHdLPrintEnergy);
565 block_bc(cr, fep->n_lambda);
566 if (fep->n_lambda > 0)
568 snew_bc(cr, fep->all_lambda, efptNR);
569 nblock_bc(cr, efptNR, fep->all_lambda);
570 for (i = 0; i < efptNR; i++)
572 snew_bc(cr, fep->all_lambda[i], fep->n_lambda);
573 nblock_bc(cr, fep->n_lambda, fep->all_lambda[i]);
576 block_bc(cr, fep->sc_alpha);
577 block_bc(cr, fep->sc_power);
578 block_bc(cr, fep->sc_r_power);
579 block_bc(cr, fep->sc_sigma);
580 block_bc(cr, fep->sc_sigma_min);
581 block_bc(cr, fep->bScCoul);
582 nblock_bc(cr, efptNR, &(fep->separate_dvdl[0]));
583 block_bc(cr, fep->dhdl_derivatives);
584 block_bc(cr, fep->dh_hist_size);
585 block_bc(cr, fep->dh_hist_spacing);
588 fprintf(debug, "after bc_fepvals\n");
592 static void bc_expandedvals(const t_commrec *cr, t_expanded *expand, int n_lambda)
594 block_bc(cr, expand->nstexpanded);
595 block_bc(cr, expand->elamstats);
596 block_bc(cr, expand->elmcmove);
597 block_bc(cr, expand->elmceq);
598 block_bc(cr, expand->equil_n_at_lam);
599 block_bc(cr, expand->equil_wl_delta);
600 block_bc(cr, expand->equil_ratio);
601 block_bc(cr, expand->equil_steps);
602 block_bc(cr, expand->equil_samples);
603 block_bc(cr, expand->lmc_seed);
604 block_bc(cr, expand->minvar);
605 block_bc(cr, expand->minvar_const);
606 block_bc(cr, expand->c_range);
607 block_bc(cr, expand->bSymmetrizedTMatrix);
608 block_bc(cr, expand->nstTij);
609 block_bc(cr, expand->lmc_repeats);
610 block_bc(cr, expand->lmc_forced_nstart);
611 block_bc(cr, expand->gibbsdeltalam);
612 block_bc(cr, expand->wl_scale);
613 block_bc(cr, expand->wl_ratio);
614 block_bc(cr, expand->init_wl_delta);
615 block_bc(cr, expand->bInit_weights);
616 snew_bc(cr, expand->init_lambda_weights, n_lambda);
617 nblock_bc(cr, n_lambda, expand->init_lambda_weights);
618 block_bc(cr, expand->mc_temp);
621 fprintf(debug, "after bc_expandedvals\n");
625 static void bc_simtempvals(const t_commrec *cr, t_simtemp *simtemp, int n_lambda)
627 block_bc(cr, simtemp->simtemp_low);
628 block_bc(cr, simtemp->simtemp_high);
629 block_bc(cr, simtemp->eSimTempScale);
630 snew_bc(cr, simtemp->temperatures, n_lambda);
631 nblock_bc(cr, n_lambda, simtemp->temperatures);
634 fprintf(debug, "after bc_simtempvals\n");
639 static void bc_swapions(const t_commrec *cr, t_swapcoords *swap)
646 /* Broadcast ion group atom indices */
647 snew_bc(cr, swap->ind, swap->nat);
648 nblock_bc(cr, swap->nat, swap->ind);
650 /* Broadcast split groups atom indices */
651 for (i = 0; i < 2; i++)
653 snew_bc(cr, swap->ind_split[i], swap->nat_split[i]);
654 nblock_bc(cr, swap->nat_split[i], swap->ind_split[i]);
657 /* Broadcast solvent group atom indices */
658 snew_bc(cr, swap->ind_sol, swap->nat_sol);
659 nblock_bc(cr, swap->nat_sol, swap->ind_sol);
663 static void bc_inputrec(const t_commrec *cr, t_inputrec *inputrec)
667 block_bc(cr, *inputrec);
669 bc_grpopts(cr, &(inputrec->opts));
671 /* even if efep is efepNO, we need to initialize to make sure that
672 * n_lambda is set to zero */
674 snew_bc(cr, inputrec->fepvals, 1);
675 if (inputrec->efep != efepNO || inputrec->bSimTemp)
677 bc_fepvals(cr, inputrec->fepvals);
679 /* need to initialize this as well because of data checked for in the logic */
680 snew_bc(cr, inputrec->expandedvals, 1);
681 if (inputrec->bExpanded)
683 bc_expandedvals(cr, inputrec->expandedvals, inputrec->fepvals->n_lambda);
685 snew_bc(cr, inputrec->simtempvals, 1);
686 if (inputrec->bSimTemp)
688 bc_simtempvals(cr, inputrec->simtempvals, inputrec->fepvals->n_lambda);
692 snew_bc(cr, inputrec->pull, 1);
693 bc_pull(cr, inputrec->pull);
697 snew_bc(cr, inputrec->rot, 1);
698 bc_rot(cr, inputrec->rot);
702 snew_bc(cr, inputrec->imd, 1);
703 bc_imd(cr, inputrec->imd);
705 for (i = 0; (i < DIM); i++)
707 bc_cosines(cr, &(inputrec->ex[i]));
708 bc_cosines(cr, &(inputrec->et[i]));
710 if (inputrec->eSwapCoords != eswapNO)
712 snew_bc(cr, inputrec->swap, 1);
713 bc_swapions(cr, inputrec->swap);
715 if (inputrec->bAdress)
717 snew_bc(cr, inputrec->adress, 1);
718 bc_adress(cr, inputrec->adress);
722 static void bc_moltype(const t_commrec *cr, t_symtab *symtab,
723 gmx_moltype_t *moltype)
725 bc_string(cr, symtab, &moltype->name);
726 bc_atoms(cr, symtab, &moltype->atoms);
729 fprintf(debug, "after bc_atoms\n");
732 bc_ilists(cr, moltype->ilist);
733 bc_block(cr, &moltype->cgs);
734 bc_blocka(cr, &moltype->excls);
737 static void bc_molblock(const t_commrec *cr, gmx_molblock_t *molb)
739 block_bc(cr, molb->type);
740 block_bc(cr, molb->nmol);
741 block_bc(cr, molb->natoms_mol);
742 block_bc(cr, molb->nposres_xA);
743 if (molb->nposres_xA > 0)
745 snew_bc(cr, molb->posres_xA, molb->nposres_xA);
746 nblock_bc(cr, molb->nposres_xA*DIM, molb->posres_xA[0]);
748 block_bc(cr, molb->nposres_xB);
749 if (molb->nposres_xB > 0)
751 snew_bc(cr, molb->posres_xB, molb->nposres_xB);
752 nblock_bc(cr, molb->nposres_xB*DIM, molb->posres_xB[0]);
756 fprintf(debug, "after bc_molblock\n");
760 static void bc_atomtypes(const t_commrec *cr, t_atomtypes *atomtypes)
764 block_bc(cr, atomtypes->nr);
768 snew_bc(cr, atomtypes->radius, nr);
769 snew_bc(cr, atomtypes->vol, nr);
770 snew_bc(cr, atomtypes->surftens, nr);
771 snew_bc(cr, atomtypes->gb_radius, nr);
772 snew_bc(cr, atomtypes->S_hct, nr);
774 nblock_bc(cr, nr, atomtypes->radius);
775 nblock_bc(cr, nr, atomtypes->vol);
776 nblock_bc(cr, nr, atomtypes->surftens);
777 nblock_bc(cr, nr, atomtypes->gb_radius);
778 nblock_bc(cr, nr, atomtypes->S_hct);
781 /*! \brief Broadcasts ir and mtop from the master to all nodes in
782 * cr->mpi_comm_mygroup. */
784 void bcast_ir_mtop(const t_commrec *cr, t_inputrec *inputrec, gmx_mtop_t *mtop)
789 fprintf(debug, "in bc_data\n");
791 bc_inputrec(cr, inputrec);
794 fprintf(debug, "after bc_inputrec\n");
796 bc_symtab(cr, &mtop->symtab);
799 fprintf(debug, "after bc_symtab\n");
801 bc_string(cr, &mtop->symtab, &mtop->name);
804 fprintf(debug, "after bc_name\n");
807 bc_ffparams(cr, &mtop->ffparams);
809 block_bc(cr, mtop->nmoltype);
810 snew_bc(cr, mtop->moltype, mtop->nmoltype);
811 for (i = 0; i < mtop->nmoltype; i++)
813 bc_moltype(cr, &mtop->symtab, &mtop->moltype[i]);
816 block_bc(cr, mtop->nmolblock);
817 snew_bc(cr, mtop->molblock, mtop->nmolblock);
818 for (i = 0; i < mtop->nmolblock; i++)
820 bc_molblock(cr, &mtop->molblock[i]);
823 block_bc(cr, mtop->natoms);
825 bc_atomtypes(cr, &mtop->atomtypes);
827 bc_block(cr, &mtop->mols);
828 bc_groups(cr, &mtop->symtab, mtop->natoms, &mtop->groups);
831 void init_parallel(t_commrec *cr, t_inputrec *inputrec,
834 bcast_ir_mtop(cr, inputrec, mtop);