Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxlib / mvdata.c
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, 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 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <string.h>
43
44 #include "typedefs.h"
45 #include "main.h"
46 #include "mvdata.h"
47 #include "network.h"
48 #include "smalloc.h"
49 #include "gmx_fatal.h"
50 #include "symtab.h"
51 #include "vec.h"
52 #include "tgroup.h"
53
54 #define   block_bc(cr,   d) gmx_bcast(     sizeof(d),     &(d), (cr))
55 /* Probably the test for (nr) > 0 in the next macro is only needed
56  * on BlueGene(/L), where IBM's MPI_Bcast will segfault after
57  * dereferencing a null pointer, even when no data is to be transferred. */
58 #define  nblock_bc(cr, nr, d) { if ((nr) > 0) {gmx_bcast((nr)*sizeof((d)[0]), (d), (cr)); }}
59 #define    snew_bc(cr, d, nr) { if (!MASTER(cr)) {snew((d), (nr)); }}
60 /* Dirty macro with bAlloc not as an argument */
61 #define nblock_abc(cr, nr, d) { if (bAlloc) {snew((d), (nr)); } nblock_bc(cr, (nr), (d)); }
62
63 static void bc_string(const t_commrec *cr, t_symtab *symtab, char ***s)
64 {
65     int handle;
66
67     if (MASTER(cr))
68     {
69         handle = lookup_symtab(symtab, *s);
70     }
71     block_bc(cr, handle);
72     if (!MASTER(cr))
73     {
74         *s = get_symtab_handle(symtab, handle);
75     }
76 }
77
78 static void bc_strings(const t_commrec *cr, t_symtab *symtab, int nr, char ****nm)
79 {
80     int     i;
81     int    *handle;
82     char ***NM;
83
84     snew(handle, nr);
85     if (MASTER(cr))
86     {
87         NM = *nm;
88         for (i = 0; (i < nr); i++)
89         {
90             handle[i] = lookup_symtab(symtab, NM[i]);
91         }
92     }
93     nblock_bc(cr, nr, handle);
94
95     if (!MASTER(cr))
96     {
97         snew_bc(cr, *nm, nr);
98         NM = *nm;
99         for (i = 0; (i < nr); i++)
100         {
101             (*nm)[i] = get_symtab_handle(symtab, handle[i]);
102         }
103     }
104     sfree(handle);
105 }
106
107 static void bc_strings_resinfo(const t_commrec *cr, t_symtab *symtab,
108                                int nr, t_resinfo *resinfo)
109 {
110     int   i;
111     int  *handle;
112
113     snew(handle, nr);
114     if (MASTER(cr))
115     {
116         for (i = 0; (i < nr); i++)
117         {
118             handle[i] = lookup_symtab(symtab, resinfo[i].name);
119         }
120     }
121     nblock_bc(cr, nr, handle);
122
123     if (!MASTER(cr))
124     {
125         for (i = 0; (i < nr); i++)
126         {
127             resinfo[i].name = get_symtab_handle(symtab, handle[i]);
128         }
129     }
130     sfree(handle);
131 }
132
133 static void bc_symtab(const t_commrec *cr, t_symtab *symtab)
134 {
135     int       i, nr, len;
136     t_symbuf *symbuf;
137
138     block_bc(cr, symtab->nr);
139     nr = symtab->nr;
140     snew_bc(cr, symtab->symbuf, 1);
141     symbuf          = symtab->symbuf;
142     symbuf->bufsize = nr;
143     snew_bc(cr, symbuf->buf, nr);
144     for (i = 0; i < nr; i++)
145     {
146         if (MASTER(cr))
147         {
148             len = strlen(symbuf->buf[i]) + 1;
149         }
150         block_bc(cr, len);
151         snew_bc(cr, symbuf->buf[i], len);
152         nblock_bc(cr, len, symbuf->buf[i]);
153     }
154 }
155
156 static void bc_block(const t_commrec *cr, t_block *block)
157 {
158     block_bc(cr, block->nr);
159     snew_bc(cr, block->index, block->nr+1);
160     nblock_bc(cr, block->nr+1, block->index);
161 }
162
163 static void bc_blocka(const t_commrec *cr, t_blocka *block)
164 {
165     block_bc(cr, block->nr);
166     snew_bc(cr, block->index, block->nr+1);
167     nblock_bc(cr, block->nr+1, block->index);
168     block_bc(cr, block->nra);
169     if (block->nra)
170     {
171         snew_bc(cr, block->a, block->nra);
172         nblock_bc(cr, block->nra, block->a);
173     }
174 }
175
176 static void bc_grps(const t_commrec *cr, t_grps grps[])
177 {
178     int i;
179
180     for (i = 0; (i < egcNR); i++)
181     {
182         block_bc(cr, grps[i].nr);
183         snew_bc(cr, grps[i].nm_ind, grps[i].nr);
184         nblock_bc(cr, grps[i].nr, grps[i].nm_ind);
185     }
186 }
187
188 static void bc_atoms(const t_commrec *cr, t_symtab *symtab, t_atoms *atoms)
189 {
190     int dummy;
191
192     block_bc(cr, atoms->nr);
193     snew_bc(cr, atoms->atom, atoms->nr);
194     nblock_bc(cr, atoms->nr, atoms->atom);
195     bc_strings(cr, symtab, atoms->nr, &atoms->atomname);
196     block_bc(cr, atoms->nres);
197     snew_bc(cr, atoms->resinfo, atoms->nres);
198     nblock_bc(cr, atoms->nres, atoms->resinfo);
199     bc_strings_resinfo(cr, symtab, atoms->nres, atoms->resinfo);
200     /* QMMM requires atomtypes to be known on all nodes as well */
201     bc_strings(cr, symtab, atoms->nr, &atoms->atomtype);
202     bc_strings(cr, symtab, atoms->nr, &atoms->atomtypeB);
203 }
204
205 static void bc_groups(const t_commrec *cr, t_symtab *symtab,
206                       int natoms, gmx_groups_t *groups)
207 {
208     int dummy;
209     int g, n;
210
211     bc_grps(cr, groups->grps);
212     block_bc(cr, groups->ngrpname);
213     bc_strings(cr, symtab, groups->ngrpname, &groups->grpname);
214     for (g = 0; g < egcNR; g++)
215     {
216         if (MASTER(cr))
217         {
218             if (groups->grpnr[g])
219             {
220                 n = natoms;
221             }
222             else
223             {
224                 n = 0;
225             }
226         }
227         block_bc(cr, n);
228         if (n == 0)
229         {
230             groups->grpnr[g] = NULL;
231         }
232         else
233         {
234             snew_bc(cr, groups->grpnr[g], n);
235             nblock_bc(cr, n, groups->grpnr[g]);
236         }
237     }
238     if (debug)
239     {
240         fprintf(debug, "after bc_groups\n");
241     }
242 }
243
244 void bcast_state(const t_commrec *cr, t_state *state)
245 {
246     int      i, nnht, nnhtp;
247     gmx_bool bAlloc;
248
249     if (!PAR(cr))
250     {
251         return;
252     }
253
254     /* Broadcasts the state sizes and flags from the master to all nodes
255      * in cr->mpi_comm_mygroup. The arrays are not broadcasted. */
256     block_bc(cr, state->natoms);
257     block_bc(cr, state->ngtc);
258     block_bc(cr, state->nnhpres);
259     block_bc(cr, state->nhchainlength);
260     block_bc(cr, state->nrng);
261     block_bc(cr, state->nrngi);
262     block_bc(cr, state->flags);
263     if (state->lambda == NULL)
264     {
265         snew_bc(cr, state->lambda, efptNR)
266     }
267
268     if (cr->dd)
269     {
270         /* We allocate dynamically in dd_partition_system. */
271         return;
272     }
273     /* The code below is reachable only by TPI and NM, so it is not
274        tested by anything. */
275
276     nnht  = (state->ngtc)*(state->nhchainlength);
277     nnhtp = (state->nnhpres)*(state->nhchainlength);
278
279     /* We still need to allocate the arrays in state for non-master
280      * ranks, which is done (implicitly via bAlloc) in the dirty,
281      * dirty nblock_abc macro. */
282     bAlloc = !MASTER(cr);
283     if (bAlloc)
284     {
285         state->nalloc = state->natoms;
286     }
287     for (i = 0; i < estNR; i++)
288     {
289         if (state->flags & (1<<i))
290         {
291             switch (i)
292             {
293                 case estLAMBDA:  nblock_bc(cr, efptNR, state->lambda); break;
294                 case estFEPSTATE: block_bc(cr, state->fep_state); break;
295                 case estBOX:     block_bc(cr, state->box); break;
296                 case estBOX_REL: block_bc(cr, state->box_rel); break;
297                 case estBOXV:    block_bc(cr, state->boxv); break;
298                 case estPRES_PREV: block_bc(cr, state->pres_prev); break;
299                 case estSVIR_PREV: block_bc(cr, state->svir_prev); break;
300                 case estFVIR_PREV: block_bc(cr, state->fvir_prev); break;
301                 case estNH_XI:   nblock_abc(cr, nnht, state->nosehoover_xi); break;
302                 case estNH_VXI:  nblock_abc(cr, nnht, state->nosehoover_vxi); break;
303                 case estNHPRES_XI:   nblock_abc(cr, nnhtp, state->nhpres_xi); break;
304                 case estNHPRES_VXI:  nblock_abc(cr, nnhtp, state->nhpres_vxi); break;
305                 case estTC_INT:  nblock_abc(cr, state->ngtc, state->therm_integral); break;
306                 case estVETA:    block_bc(cr, state->veta); break;
307                 case estVOL0:    block_bc(cr, state->vol0); break;
308                 case estX:       nblock_abc(cr, state->natoms, state->x); break;
309                 case estV:       nblock_abc(cr, state->natoms, state->v); break;
310                 case estSDX:     nblock_abc(cr, state->natoms, state->sd_X); break;
311                 case estCGP:     nblock_abc(cr, state->natoms, state->cg_p); break;
312                 case estLD_RNG:  if (state->nrngi == 1)
313                     {
314                         nblock_abc(cr, state->nrng, state->ld_rng);
315                     }
316                     break;
317                 case estLD_RNGI: if (state->nrngi == 1)
318                     {
319                         nblock_abc(cr, state->nrngi, state->ld_rngi);
320                     }
321                     break;
322                 case estDISRE_INITF: block_bc(cr, state->hist.disre_initf); break;
323                 case estDISRE_RM3TAV:
324                     block_bc(cr, state->hist.ndisrepairs);
325                     nblock_abc(cr, state->hist.ndisrepairs, state->hist.disre_rm3tav);
326                     break;
327                 case estORIRE_INITF: block_bc(cr, state->hist.orire_initf); break;
328                 case estORIRE_DTAV:
329                     block_bc(cr, state->hist.norire_Dtav);
330                     nblock_abc(cr, state->hist.norire_Dtav, state->hist.orire_Dtav);
331                     break;
332                 default:
333                     gmx_fatal(FARGS,
334                               "Communication is not implemented for %s in bcast_state",
335                               est_names[i]);
336             }
337         }
338     }
339 }
340
341 static void bc_ilists(const t_commrec *cr, t_ilist *ilist)
342 {
343     int ftype;
344
345     /* Here we only communicate the non-zero length ilists */
346     if (MASTER(cr))
347     {
348         for (ftype = 0; ftype < F_NRE; ftype++)
349         {
350             if (ilist[ftype].nr > 0)
351             {
352                 block_bc(cr, ftype);
353                 block_bc(cr, ilist[ftype].nr);
354                 nblock_bc(cr, ilist[ftype].nr, ilist[ftype].iatoms);
355             }
356         }
357         ftype = -1;
358         block_bc(cr, ftype);
359     }
360     else
361     {
362         for (ftype = 0; ftype < F_NRE; ftype++)
363         {
364             ilist[ftype].nr = 0;
365         }
366         do
367         {
368             block_bc(cr, ftype);
369             if (ftype >= 0)
370             {
371                 block_bc(cr, ilist[ftype].nr);
372                 snew_bc(cr, ilist[ftype].iatoms, ilist[ftype].nr);
373                 nblock_bc(cr, ilist[ftype].nr, ilist[ftype].iatoms);
374             }
375         }
376         while (ftype >= 0);
377     }
378
379     if (debug)
380     {
381         fprintf(debug, "after bc_ilists\n");
382     }
383 }
384
385 static void bc_cmap(const t_commrec *cr, gmx_cmap_t *cmap_grid)
386 {
387     int i, j, nelem, ngrid;
388
389     block_bc(cr, cmap_grid->ngrid);
390     block_bc(cr, cmap_grid->grid_spacing);
391
392     ngrid = cmap_grid->ngrid;
393     nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
394
395     if (ngrid > 0)
396     {
397         snew_bc(cr, cmap_grid->cmapdata, ngrid);
398
399         for (i = 0; i < ngrid; i++)
400         {
401             snew_bc(cr, cmap_grid->cmapdata[i].cmap, 4*nelem);
402             nblock_bc(cr, 4*nelem, cmap_grid->cmapdata[i].cmap);
403         }
404     }
405 }
406
407 static void bc_ffparams(const t_commrec *cr, gmx_ffparams_t *ffp)
408 {
409     int i;
410
411     block_bc(cr, ffp->ntypes);
412     block_bc(cr, ffp->atnr);
413     snew_bc(cr, ffp->functype, ffp->ntypes);
414     snew_bc(cr, ffp->iparams, ffp->ntypes);
415     nblock_bc(cr, ffp->ntypes, ffp->functype);
416     nblock_bc(cr, ffp->ntypes, ffp->iparams);
417     block_bc(cr, ffp->reppow);
418     block_bc(cr, ffp->fudgeQQ);
419     bc_cmap(cr, &ffp->cmap_grid);
420 }
421
422 static void bc_grpopts(const t_commrec *cr, t_grpopts *g)
423 {
424     int i, n;
425
426     block_bc(cr, g->ngtc);
427     block_bc(cr, g->ngacc);
428     block_bc(cr, g->ngfrz);
429     block_bc(cr, g->ngener);
430     snew_bc(cr, g->nrdf, g->ngtc);
431     snew_bc(cr, g->tau_t, g->ngtc);
432     snew_bc(cr, g->ref_t, g->ngtc);
433     snew_bc(cr, g->acc, g->ngacc);
434     snew_bc(cr, g->nFreeze, g->ngfrz);
435     snew_bc(cr, g->egp_flags, g->ngener*g->ngener);
436
437     nblock_bc(cr, g->ngtc, g->nrdf);
438     nblock_bc(cr, g->ngtc, g->tau_t);
439     nblock_bc(cr, g->ngtc, g->ref_t);
440     nblock_bc(cr, g->ngacc, g->acc);
441     nblock_bc(cr, g->ngfrz, g->nFreeze);
442     nblock_bc(cr, g->ngener*g->ngener, g->egp_flags);
443     snew_bc(cr, g->annealing, g->ngtc);
444     snew_bc(cr, g->anneal_npoints, g->ngtc);
445     snew_bc(cr, g->anneal_time, g->ngtc);
446     snew_bc(cr, g->anneal_temp, g->ngtc);
447     nblock_bc(cr, g->ngtc, g->annealing);
448     nblock_bc(cr, g->ngtc, g->anneal_npoints);
449     for (i = 0; (i < g->ngtc); i++)
450     {
451         n = g->anneal_npoints[i];
452         if (n > 0)
453         {
454             snew_bc(cr, g->anneal_time[i], n);
455             snew_bc(cr, g->anneal_temp[i], n);
456             nblock_bc(cr, n, g->anneal_time[i]);
457             nblock_bc(cr, n, g->anneal_temp[i]);
458         }
459     }
460
461     /* QMMM stuff, see inputrec */
462     block_bc(cr, g->ngQM);
463     snew_bc(cr, g->QMmethod, g->ngQM);
464     snew_bc(cr, g->QMbasis, g->ngQM);
465     snew_bc(cr, g->QMcharge, g->ngQM);
466     snew_bc(cr, g->QMmult, g->ngQM);
467     snew_bc(cr, g->bSH, g->ngQM);
468     snew_bc(cr, g->CASorbitals, g->ngQM);
469     snew_bc(cr, g->CASelectrons, g->ngQM);
470     snew_bc(cr, g->SAon, g->ngQM);
471     snew_bc(cr, g->SAoff, g->ngQM);
472     snew_bc(cr, g->SAsteps, g->ngQM);
473
474     if (g->ngQM)
475     {
476         nblock_bc(cr, g->ngQM, g->QMmethod);
477         nblock_bc(cr, g->ngQM, g->QMbasis);
478         nblock_bc(cr, g->ngQM, g->QMcharge);
479         nblock_bc(cr, g->ngQM, g->QMmult);
480         nblock_bc(cr, g->ngQM, g->bSH);
481         nblock_bc(cr, g->ngQM, g->CASorbitals);
482         nblock_bc(cr, g->ngQM, g->CASelectrons);
483         nblock_bc(cr, g->ngQM, g->SAon);
484         nblock_bc(cr, g->ngQM, g->SAoff);
485         nblock_bc(cr, g->ngQM, g->SAsteps);
486         /* end of QMMM stuff */
487     }
488 }
489
490 static void bc_cosines(const t_commrec *cr, t_cosines *cs)
491 {
492     block_bc(cr, cs->n);
493     snew_bc(cr, cs->a, cs->n);
494     snew_bc(cr, cs->phi, cs->n);
495     if (cs->n > 0)
496     {
497         nblock_bc(cr, cs->n, cs->a);
498         nblock_bc(cr, cs->n, cs->phi);
499     }
500 }
501
502 static void bc_pull_group(const t_commrec *cr, t_pull_group *pgrp)
503 {
504     block_bc(cr, *pgrp);
505     if (pgrp->nat > 0)
506     {
507         snew_bc(cr, pgrp->ind, pgrp->nat);
508         nblock_bc(cr, pgrp->nat, pgrp->ind);
509     }
510     if (pgrp->nweight > 0)
511     {
512         snew_bc(cr, pgrp->weight, pgrp->nweight);
513         nblock_bc(cr, pgrp->nweight, pgrp->weight);
514     }
515 }
516
517 static void bc_pull(const t_commrec *cr, t_pull *pull)
518 {
519     int g;
520
521     block_bc(cr, *pull);
522     snew_bc(cr, pull->group, pull->ngroup);
523     for (g = 0; g < pull->ngroup; g++)
524     {
525         bc_pull_group(cr, &pull->group[g]);
526     }
527     snew_bc(cr, pull->coord, pull->ncoord);
528     nblock_bc(cr, pull->ncoord, pull->coord);
529 }
530
531 static void bc_rotgrp(const t_commrec *cr, t_rotgrp *rotg)
532 {
533     block_bc(cr, *rotg);
534     if (rotg->nat > 0)
535     {
536         snew_bc(cr, rotg->ind, rotg->nat);
537         nblock_bc(cr, rotg->nat, rotg->ind);
538         snew_bc(cr, rotg->x_ref, rotg->nat);
539         nblock_bc(cr, rotg->nat, rotg->x_ref);
540     }
541 }
542
543 static void bc_rot(const t_commrec *cr, t_rot *rot)
544 {
545     int g;
546
547     block_bc(cr, *rot);
548     snew_bc(cr, rot->grp, rot->ngrp);
549     for (g = 0; g < rot->ngrp; g++)
550     {
551         bc_rotgrp(cr, &rot->grp[g]);
552     }
553 }
554
555 static void bc_adress(const t_commrec *cr, t_adress *adress)
556 {
557     block_bc(cr, *adress);
558     if (adress->n_tf_grps > 0)
559     {
560         snew_bc(cr, adress->tf_table_index, adress->n_tf_grps);
561         nblock_bc(cr, adress->n_tf_grps, adress->tf_table_index);
562     }
563     if (adress->n_energy_grps > 0)
564     {
565         snew_bc(cr, adress->group_explicit, adress->n_energy_grps);
566         nblock_bc(cr, adress->n_energy_grps, adress->group_explicit);
567     }
568 }
569 static void bc_fepvals(const t_commrec *cr, t_lambda *fep)
570 {
571     gmx_bool bAlloc = TRUE;
572     int      i;
573
574     block_bc(cr, fep->nstdhdl);
575     block_bc(cr, fep->init_lambda);
576     block_bc(cr, fep->init_fep_state);
577     block_bc(cr, fep->delta_lambda);
578     block_bc(cr, fep->bPrintEnergy);
579     block_bc(cr, fep->n_lambda);
580     if (fep->n_lambda > 0)
581     {
582         snew_bc(cr, fep->all_lambda, efptNR);
583         nblock_bc(cr, efptNR, fep->all_lambda);
584         for (i = 0; i < efptNR; i++)
585         {
586             snew_bc(cr, fep->all_lambda[i], fep->n_lambda);
587             nblock_bc(cr, fep->n_lambda, fep->all_lambda[i]);
588         }
589     }
590     block_bc(cr, fep->sc_alpha);
591     block_bc(cr, fep->sc_power);
592     block_bc(cr, fep->sc_r_power);
593     block_bc(cr, fep->sc_sigma);
594     block_bc(cr, fep->sc_sigma_min);
595     block_bc(cr, fep->bScCoul);
596     nblock_bc(cr, efptNR, &(fep->separate_dvdl[0]));
597     block_bc(cr, fep->dhdl_derivatives);
598     block_bc(cr, fep->dh_hist_size);
599     block_bc(cr, fep->dh_hist_spacing);
600     if (debug)
601     {
602         fprintf(debug, "after bc_fepvals\n");
603     }
604 }
605
606 static void bc_expandedvals(const t_commrec *cr, t_expanded *expand, int n_lambda)
607 {
608     gmx_bool bAlloc = TRUE;
609     int      i;
610
611     block_bc(cr, expand->nstexpanded);
612     block_bc(cr, expand->elamstats);
613     block_bc(cr, expand->elmcmove);
614     block_bc(cr, expand->elmceq);
615     block_bc(cr, expand->equil_n_at_lam);
616     block_bc(cr, expand->equil_wl_delta);
617     block_bc(cr, expand->equil_ratio);
618     block_bc(cr, expand->equil_steps);
619     block_bc(cr, expand->equil_samples);
620     block_bc(cr, expand->lmc_seed);
621     block_bc(cr, expand->minvar);
622     block_bc(cr, expand->minvar_const);
623     block_bc(cr, expand->c_range);
624     block_bc(cr, expand->bSymmetrizedTMatrix);
625     block_bc(cr, expand->nstTij);
626     block_bc(cr, expand->lmc_repeats);
627     block_bc(cr, expand->lmc_forced_nstart);
628     block_bc(cr, expand->gibbsdeltalam);
629     block_bc(cr, expand->wl_scale);
630     block_bc(cr, expand->wl_ratio);
631     block_bc(cr, expand->init_wl_delta);
632     block_bc(cr, expand->bInit_weights);
633     snew_bc(cr, expand->init_lambda_weights, n_lambda);
634     nblock_bc(cr, n_lambda, expand->init_lambda_weights);
635     block_bc(cr, expand->mc_temp);
636     if (debug)
637     {
638         fprintf(debug, "after bc_expandedvals\n");
639     }
640 }
641
642 static void bc_simtempvals(const t_commrec *cr, t_simtemp *simtemp, int n_lambda)
643 {
644     gmx_bool bAlloc = TRUE;
645     int      i;
646
647     block_bc(cr, simtemp->simtemp_low);
648     block_bc(cr, simtemp->simtemp_high);
649     block_bc(cr, simtemp->eSimTempScale);
650     snew_bc(cr, simtemp->temperatures, n_lambda);
651     nblock_bc(cr, n_lambda, simtemp->temperatures);
652     if (debug)
653     {
654         fprintf(debug, "after bc_simtempvals\n");
655     }
656 }
657
658
659 static void bc_swapions(const t_commrec *cr, t_swapcoords *swap)
660 {
661     int i;
662
663
664     block_bc(cr, *swap);
665
666     /* Broadcast ion group atom indices */
667     snew_bc(cr, swap->ind, swap->nat);
668     nblock_bc(cr, swap->nat, swap->ind);
669
670     /* Broadcast split groups atom indices */
671     for (i = 0; i < 2; i++)
672     {
673         snew_bc(cr, swap->ind_split[i], swap->nat_split[i]);
674         nblock_bc(cr, swap->nat_split[i], swap->ind_split[i]);
675     }
676
677     /* Broadcast solvent group atom indices */
678     snew_bc(cr, swap->ind_sol, swap->nat_sol);
679     nblock_bc(cr, swap->nat_sol, swap->ind_sol);
680 }
681
682
683 static void bc_inputrec(const t_commrec *cr, t_inputrec *inputrec)
684 {
685     gmx_bool bAlloc = TRUE;
686     int      i;
687
688     block_bc(cr, *inputrec);
689
690     bc_grpopts(cr, &(inputrec->opts));
691
692     /* even if efep is efepNO, we need to initialize to make sure that
693      * n_lambda is set to zero */
694
695     snew_bc(cr, inputrec->fepvals, 1);
696     if (inputrec->efep != efepNO || inputrec->bSimTemp)
697     {
698         bc_fepvals(cr, inputrec->fepvals);
699     }
700     /* need to initialize this as well because of data checked for in the logic */
701     snew_bc(cr, inputrec->expandedvals, 1);
702     if (inputrec->bExpanded)
703     {
704         bc_expandedvals(cr, inputrec->expandedvals, inputrec->fepvals->n_lambda);
705     }
706     snew_bc(cr, inputrec->simtempvals, 1);
707     if (inputrec->bSimTemp)
708     {
709         bc_simtempvals(cr, inputrec->simtempvals, inputrec->fepvals->n_lambda);
710     }
711     if (inputrec->ePull != epullNO)
712     {
713         snew_bc(cr, inputrec->pull, 1);
714         bc_pull(cr, inputrec->pull);
715     }
716     if (inputrec->bRot)
717     {
718         snew_bc(cr, inputrec->rot, 1);
719         bc_rot(cr, inputrec->rot);
720     }
721     for (i = 0; (i < DIM); i++)
722     {
723         bc_cosines(cr, &(inputrec->ex[i]));
724         bc_cosines(cr, &(inputrec->et[i]));
725     }
726     if (inputrec->eSwapCoords != eswapNO)
727     {
728         snew_bc(cr, inputrec->swap, 1);
729         bc_swapions(cr, inputrec->swap);
730     }
731     if (inputrec->bAdress)
732     {
733         snew_bc(cr, inputrec->adress, 1);
734         bc_adress(cr, inputrec->adress);
735     }
736 }
737
738 static void bc_moltype(const t_commrec *cr, t_symtab *symtab,
739                        gmx_moltype_t *moltype)
740 {
741     bc_string(cr, symtab, &moltype->name);
742     bc_atoms(cr, symtab, &moltype->atoms);
743     if (debug)
744     {
745         fprintf(debug, "after bc_atoms\n");
746     }
747
748     bc_ilists(cr, moltype->ilist);
749     bc_block(cr, &moltype->cgs);
750     bc_blocka(cr, &moltype->excls);
751 }
752
753 static void bc_molblock(const t_commrec *cr, gmx_molblock_t *molb)
754 {
755     gmx_bool bAlloc = TRUE;
756
757     block_bc(cr, molb->type);
758     block_bc(cr, molb->nmol);
759     block_bc(cr, molb->natoms_mol);
760     block_bc(cr, molb->nposres_xA);
761     if (molb->nposres_xA > 0)
762     {
763         snew_bc(cr, molb->posres_xA, molb->nposres_xA);
764         nblock_bc(cr, molb->nposres_xA*DIM, molb->posres_xA[0]);
765     }
766     block_bc(cr, molb->nposres_xB);
767     if (molb->nposres_xB > 0)
768     {
769         snew_bc(cr, molb->posres_xB, molb->nposres_xB);
770         nblock_bc(cr, molb->nposres_xB*DIM, molb->posres_xB[0]);
771     }
772     if (debug)
773     {
774         fprintf(debug, "after bc_molblock\n");
775     }
776 }
777
778 static void bc_atomtypes(const t_commrec *cr, t_atomtypes *atomtypes)
779 {
780     int nr;
781
782     block_bc(cr, atomtypes->nr);
783
784     nr = atomtypes->nr;
785
786     snew_bc(cr, atomtypes->radius, nr);
787     snew_bc(cr, atomtypes->vol, nr);
788     snew_bc(cr, atomtypes->surftens, nr);
789     snew_bc(cr, atomtypes->gb_radius, nr);
790     snew_bc(cr, atomtypes->S_hct, nr);
791
792     nblock_bc(cr, nr, atomtypes->radius);
793     nblock_bc(cr, nr, atomtypes->vol);
794     nblock_bc(cr, nr, atomtypes->surftens);
795     nblock_bc(cr, nr, atomtypes->gb_radius);
796     nblock_bc(cr, nr, atomtypes->S_hct);
797 }
798
799
800 void bcast_ir_mtop(const t_commrec *cr, t_inputrec *inputrec, gmx_mtop_t *mtop)
801 {
802     int i;
803     if (debug)
804     {
805         fprintf(debug, "in bc_data\n");
806     }
807     bc_inputrec(cr, inputrec);
808     if (debug)
809     {
810         fprintf(debug, "after bc_inputrec\n");
811     }
812     bc_symtab(cr, &mtop->symtab);
813     if (debug)
814     {
815         fprintf(debug, "after bc_symtab\n");
816     }
817     bc_string(cr, &mtop->symtab, &mtop->name);
818     if (debug)
819     {
820         fprintf(debug, "after bc_name\n");
821     }
822
823     bc_ffparams(cr, &mtop->ffparams);
824
825     block_bc(cr, mtop->nmoltype);
826     snew_bc(cr, mtop->moltype, mtop->nmoltype);
827     for (i = 0; i < mtop->nmoltype; i++)
828     {
829         bc_moltype(cr, &mtop->symtab, &mtop->moltype[i]);
830     }
831
832     block_bc(cr, mtop->nmolblock);
833     snew_bc(cr, mtop->molblock, mtop->nmolblock);
834     for (i = 0; i < mtop->nmolblock; i++)
835     {
836         bc_molblock(cr, &mtop->molblock[i]);
837     }
838
839     block_bc(cr, mtop->natoms);
840
841     bc_atomtypes(cr, &mtop->atomtypes);
842
843     bc_block(cr, &mtop->mols);
844     bc_groups(cr, &mtop->symtab, mtop->natoms, &mtop->groups);
845 }