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_setup(const t_commrec *cr, t_state *state)
245 {
246     block_bc(cr, state->natoms);
247     block_bc(cr, state->ngtc);
248     block_bc(cr, state->nnhpres);
249     block_bc(cr, state->nhchainlength);
250     block_bc(cr, state->nrng);
251     block_bc(cr, state->nrngi);
252     block_bc(cr, state->flags);
253     if (state->lambda == NULL)
254     {
255         snew_bc(cr, state->lambda, efptNR)
256     }
257 }
258
259 void bcast_state(const t_commrec *cr, t_state *state, gmx_bool bAlloc)
260 {
261     int i, nnht, nnhtp;
262
263     bcast_state_setup(cr, state);
264
265     nnht  = (state->ngtc)*(state->nhchainlength);
266     nnhtp = (state->nnhpres)*(state->nhchainlength);
267
268     if (MASTER(cr))
269     {
270         bAlloc = FALSE;
271     }
272     if (bAlloc)
273     {
274         state->nalloc = state->natoms;
275     }
276     for (i = 0; i < estNR; i++)
277     {
278         if (state->flags & (1<<i))
279         {
280             switch (i)
281             {
282                 case estLAMBDA:  nblock_bc(cr, efptNR, state->lambda); break;
283                 case estFEPSTATE: block_bc(cr, state->fep_state); break;
284                 case estBOX:     block_bc(cr, state->box); break;
285                 case estBOX_REL: block_bc(cr, state->box_rel); break;
286                 case estBOXV:    block_bc(cr, state->boxv); break;
287                 case estPRES_PREV: block_bc(cr, state->pres_prev); break;
288                 case estSVIR_PREV: block_bc(cr, state->svir_prev); break;
289                 case estFVIR_PREV: block_bc(cr, state->fvir_prev); break;
290                 case estNH_XI:   nblock_abc(cr, nnht, state->nosehoover_xi); break;
291                 case estNH_VXI:  nblock_abc(cr, nnht, state->nosehoover_vxi); break;
292                 case estNHPRES_XI:   nblock_abc(cr, nnhtp, state->nhpres_xi); break;
293                 case estNHPRES_VXI:  nblock_abc(cr, nnhtp, state->nhpres_vxi); break;
294                 case estTC_INT:  nblock_abc(cr, state->ngtc, state->therm_integral); break;
295                 case estVETA:    block_bc(cr, state->veta); break;
296                 case estVOL0:    block_bc(cr, state->vol0); break;
297                 case estX:       nblock_abc(cr, state->natoms, state->x); break;
298                 case estV:       nblock_abc(cr, state->natoms, state->v); break;
299                 case estSDX:     nblock_abc(cr, state->natoms, state->sd_X); break;
300                 case estCGP:     nblock_abc(cr, state->natoms, state->cg_p); break;
301                 case estLD_RNG:  if (state->nrngi == 1)
302                     {
303                         nblock_abc(cr, state->nrng, state->ld_rng);
304                     }
305                     break;
306                 case estLD_RNGI: if (state->nrngi == 1)
307                     {
308                         nblock_abc(cr, state->nrngi, state->ld_rngi);
309                     }
310                     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);
315                     break;
316                 case estORIRE_INITF: block_bc(cr, state->hist.orire_initf); break;
317                 case estORIRE_DTAV:
318                     block_bc(cr, state->hist.norire_Dtav);
319                     nblock_abc(cr, state->hist.norire_Dtav, state->hist.orire_Dtav);
320                     break;
321                 default:
322                     gmx_fatal(FARGS,
323                               "Communication is not implemented for %s in bcast_state",
324                               est_names[i]);
325             }
326         }
327     }
328 }
329
330 static void bc_ilists(const t_commrec *cr, t_ilist *ilist)
331 {
332     int ftype;
333
334     /* Here we only communicate the non-zero length ilists */
335     if (MASTER(cr))
336     {
337         for (ftype = 0; ftype < F_NRE; ftype++)
338         {
339             if (ilist[ftype].nr > 0)
340             {
341                 block_bc(cr, ftype);
342                 block_bc(cr, ilist[ftype].nr);
343                 nblock_bc(cr, ilist[ftype].nr, ilist[ftype].iatoms);
344             }
345         }
346         ftype = -1;
347         block_bc(cr, ftype);
348     }
349     else
350     {
351         for (ftype = 0; ftype < F_NRE; ftype++)
352         {
353             ilist[ftype].nr = 0;
354         }
355         do
356         {
357             block_bc(cr, ftype);
358             if (ftype >= 0)
359             {
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);
363             }
364         }
365         while (ftype >= 0);
366     }
367
368     if (debug)
369     {
370         fprintf(debug, "after bc_ilists\n");
371     }
372 }
373
374 static void bc_cmap(const t_commrec *cr, gmx_cmap_t *cmap_grid)
375 {
376     int i, j, nelem, ngrid;
377
378     block_bc(cr, cmap_grid->ngrid);
379     block_bc(cr, cmap_grid->grid_spacing);
380
381     ngrid = cmap_grid->ngrid;
382     nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
383
384     if (ngrid > 0)
385     {
386         snew_bc(cr, cmap_grid->cmapdata, ngrid);
387
388         for (i = 0; i < ngrid; i++)
389         {
390             snew_bc(cr, cmap_grid->cmapdata[i].cmap, 4*nelem);
391             nblock_bc(cr, 4*nelem, cmap_grid->cmapdata[i].cmap);
392         }
393     }
394 }
395
396 static void bc_ffparams(const t_commrec *cr, gmx_ffparams_t *ffp)
397 {
398     int i;
399
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);
409 }
410
411 static void bc_grpopts(const t_commrec *cr, t_grpopts *g)
412 {
413     int i, n;
414
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);
425
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++)
439     {
440         n = g->anneal_npoints[i];
441         if (n > 0)
442         {
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]);
447         }
448     }
449
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);
462
463     if (g->ngQM)
464     {
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 */
476     }
477 }
478
479 static void bc_cosines(const t_commrec *cr, t_cosines *cs)
480 {
481     block_bc(cr, cs->n);
482     snew_bc(cr, cs->a, cs->n);
483     snew_bc(cr, cs->phi, cs->n);
484     if (cs->n > 0)
485     {
486         nblock_bc(cr, cs->n, cs->a);
487         nblock_bc(cr, cs->n, cs->phi);
488     }
489 }
490
491 static void bc_pull_group(const t_commrec *cr, t_pull_group *pgrp)
492 {
493     block_bc(cr, *pgrp);
494     if (pgrp->nat > 0)
495     {
496         snew_bc(cr, pgrp->ind, pgrp->nat);
497         nblock_bc(cr, pgrp->nat, pgrp->ind);
498     }
499     if (pgrp->nweight > 0)
500     {
501         snew_bc(cr, pgrp->weight, pgrp->nweight);
502         nblock_bc(cr, pgrp->nweight, pgrp->weight);
503     }
504 }
505
506 static void bc_pull(const t_commrec *cr, t_pull *pull)
507 {
508     int g;
509
510     block_bc(cr, *pull);
511     snew_bc(cr, pull->group, pull->ngroup);
512     for (g = 0; g < pull->ngroup; g++)
513     {
514         bc_pull_group(cr, &pull->group[g]);
515     }
516     snew_bc(cr, pull->coord, pull->ncoord);
517     nblock_bc(cr, pull->ncoord, pull->coord);
518 }
519
520 static void bc_rotgrp(const t_commrec *cr, t_rotgrp *rotg)
521 {
522     block_bc(cr, *rotg);
523     if (rotg->nat > 0)
524     {
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);
529     }
530 }
531
532 static void bc_rot(const t_commrec *cr, t_rot *rot)
533 {
534     int g;
535
536     block_bc(cr, *rot);
537     snew_bc(cr, rot->grp, rot->ngrp);
538     for (g = 0; g < rot->ngrp; g++)
539     {
540         bc_rotgrp(cr, &rot->grp[g]);
541     }
542 }
543
544 static void bc_adress(const t_commrec *cr, t_adress *adress)
545 {
546     block_bc(cr, *adress);
547     if (adress->n_tf_grps > 0)
548     {
549         snew_bc(cr, adress->tf_table_index, adress->n_tf_grps);
550         nblock_bc(cr, adress->n_tf_grps, adress->tf_table_index);
551     }
552     if (adress->n_energy_grps > 0)
553     {
554         snew_bc(cr, adress->group_explicit, adress->n_energy_grps);
555         nblock_bc(cr, adress->n_energy_grps, adress->group_explicit);
556     }
557 }
558 static void bc_fepvals(const t_commrec *cr, t_lambda *fep)
559 {
560     gmx_bool bAlloc = TRUE;
561     int      i;
562
563     block_bc(cr, fep->nstdhdl);
564     block_bc(cr, fep->init_lambda);
565     block_bc(cr, fep->init_fep_state);
566     block_bc(cr, fep->delta_lambda);
567     block_bc(cr, fep->bPrintEnergy);
568     block_bc(cr, fep->n_lambda);
569     if (fep->n_lambda > 0)
570     {
571         snew_bc(cr, fep->all_lambda, efptNR);
572         nblock_bc(cr, efptNR, fep->all_lambda);
573         for (i = 0; i < efptNR; i++)
574         {
575             snew_bc(cr, fep->all_lambda[i], fep->n_lambda);
576             nblock_bc(cr, fep->n_lambda, fep->all_lambda[i]);
577         }
578     }
579     block_bc(cr, fep->sc_alpha);
580     block_bc(cr, fep->sc_power);
581     block_bc(cr, fep->sc_r_power);
582     block_bc(cr, fep->sc_sigma);
583     block_bc(cr, fep->sc_sigma_min);
584     block_bc(cr, fep->bScCoul);
585     nblock_bc(cr, efptNR, &(fep->separate_dvdl[0]));
586     block_bc(cr, fep->dhdl_derivatives);
587     block_bc(cr, fep->dh_hist_size);
588     block_bc(cr, fep->dh_hist_spacing);
589     if (debug)
590     {
591         fprintf(debug, "after bc_fepvals\n");
592     }
593 }
594
595 static void bc_expandedvals(const t_commrec *cr, t_expanded *expand, int n_lambda)
596 {
597     gmx_bool bAlloc = TRUE;
598     int      i;
599
600     block_bc(cr, expand->nstexpanded);
601     block_bc(cr, expand->elamstats);
602     block_bc(cr, expand->elmcmove);
603     block_bc(cr, expand->elmceq);
604     block_bc(cr, expand->equil_n_at_lam);
605     block_bc(cr, expand->equil_wl_delta);
606     block_bc(cr, expand->equil_ratio);
607     block_bc(cr, expand->equil_steps);
608     block_bc(cr, expand->equil_samples);
609     block_bc(cr, expand->lmc_seed);
610     block_bc(cr, expand->minvar);
611     block_bc(cr, expand->minvar_const);
612     block_bc(cr, expand->c_range);
613     block_bc(cr, expand->bSymmetrizedTMatrix);
614     block_bc(cr, expand->nstTij);
615     block_bc(cr, expand->lmc_repeats);
616     block_bc(cr, expand->lmc_forced_nstart);
617     block_bc(cr, expand->gibbsdeltalam);
618     block_bc(cr, expand->wl_scale);
619     block_bc(cr, expand->wl_ratio);
620     block_bc(cr, expand->init_wl_delta);
621     block_bc(cr, expand->bInit_weights);
622     snew_bc(cr, expand->init_lambda_weights, n_lambda);
623     nblock_bc(cr, n_lambda, expand->init_lambda_weights);
624     block_bc(cr, expand->mc_temp);
625     if (debug)
626     {
627         fprintf(debug, "after bc_expandedvals\n");
628     }
629 }
630
631 static void bc_simtempvals(const t_commrec *cr, t_simtemp *simtemp, int n_lambda)
632 {
633     gmx_bool bAlloc = TRUE;
634     int      i;
635
636     block_bc(cr, simtemp->simtemp_low);
637     block_bc(cr, simtemp->simtemp_high);
638     block_bc(cr, simtemp->eSimTempScale);
639     snew_bc(cr, simtemp->temperatures, n_lambda);
640     nblock_bc(cr, n_lambda, simtemp->temperatures);
641     if (debug)
642     {
643         fprintf(debug, "after bc_simtempvals\n");
644     }
645 }
646
647 static void bc_inputrec(const t_commrec *cr, t_inputrec *inputrec)
648 {
649     gmx_bool bAlloc = TRUE;
650     int      i;
651
652     block_bc(cr, *inputrec);
653
654     bc_grpopts(cr, &(inputrec->opts));
655
656     /* even if efep is efepNO, we need to initialize to make sure that
657      * n_lambda is set to zero */
658
659     snew_bc(cr, inputrec->fepvals, 1);
660     if (inputrec->efep != efepNO || inputrec->bSimTemp)
661     {
662         bc_fepvals(cr, inputrec->fepvals);
663     }
664     /* need to initialize this as well because of data checked for in the logic */
665     snew_bc(cr, inputrec->expandedvals, 1);
666     if (inputrec->bExpanded)
667     {
668         bc_expandedvals(cr, inputrec->expandedvals, inputrec->fepvals->n_lambda);
669     }
670     snew_bc(cr, inputrec->simtempvals, 1);
671     if (inputrec->bSimTemp)
672     {
673         bc_simtempvals(cr, inputrec->simtempvals, inputrec->fepvals->n_lambda);
674     }
675     if (inputrec->ePull != epullNO)
676     {
677         snew_bc(cr, inputrec->pull, 1);
678         bc_pull(cr, inputrec->pull);
679     }
680     if (inputrec->bRot)
681     {
682         snew_bc(cr, inputrec->rot, 1);
683         bc_rot(cr, inputrec->rot);
684     }
685     for (i = 0; (i < DIM); i++)
686     {
687         bc_cosines(cr, &(inputrec->ex[i]));
688         bc_cosines(cr, &(inputrec->et[i]));
689     }
690     if (inputrec->bAdress)
691     {
692         snew_bc(cr, inputrec->adress, 1);
693         bc_adress(cr, inputrec->adress);
694     }
695 }
696
697 static void bc_moltype(const t_commrec *cr, t_symtab *symtab,
698                        gmx_moltype_t *moltype)
699 {
700     bc_string(cr, symtab, &moltype->name);
701     bc_atoms(cr, symtab, &moltype->atoms);
702     if (debug)
703     {
704         fprintf(debug, "after bc_atoms\n");
705     }
706
707     bc_ilists(cr, moltype->ilist);
708     bc_block(cr, &moltype->cgs);
709     bc_blocka(cr, &moltype->excls);
710 }
711
712 static void bc_molblock(const t_commrec *cr, gmx_molblock_t *molb)
713 {
714     gmx_bool bAlloc = TRUE;
715
716     block_bc(cr, molb->type);
717     block_bc(cr, molb->nmol);
718     block_bc(cr, molb->natoms_mol);
719     block_bc(cr, molb->nposres_xA);
720     if (molb->nposres_xA > 0)
721     {
722         snew_bc(cr, molb->posres_xA, molb->nposres_xA);
723         nblock_bc(cr, molb->nposres_xA*DIM, molb->posres_xA[0]);
724     }
725     block_bc(cr, molb->nposres_xB);
726     if (molb->nposres_xB > 0)
727     {
728         snew_bc(cr, molb->posres_xB, molb->nposres_xB);
729         nblock_bc(cr, molb->nposres_xB*DIM, molb->posres_xB[0]);
730     }
731     if (debug)
732     {
733         fprintf(debug, "after bc_molblock\n");
734     }
735 }
736
737 static void bc_atomtypes(const t_commrec *cr, t_atomtypes *atomtypes)
738 {
739     int nr;
740
741     block_bc(cr, atomtypes->nr);
742
743     nr = atomtypes->nr;
744
745     snew_bc(cr, atomtypes->radius, nr);
746     snew_bc(cr, atomtypes->vol, nr);
747     snew_bc(cr, atomtypes->surftens, nr);
748     snew_bc(cr, atomtypes->gb_radius, nr);
749     snew_bc(cr, atomtypes->S_hct, nr);
750
751     nblock_bc(cr, nr, atomtypes->radius);
752     nblock_bc(cr, nr, atomtypes->vol);
753     nblock_bc(cr, nr, atomtypes->surftens);
754     nblock_bc(cr, nr, atomtypes->gb_radius);
755     nblock_bc(cr, nr, atomtypes->S_hct);
756 }
757
758
759 void bcast_ir_mtop(const t_commrec *cr, t_inputrec *inputrec, gmx_mtop_t *mtop)
760 {
761     int i;
762     if (debug)
763     {
764         fprintf(debug, "in bc_data\n");
765     }
766     bc_inputrec(cr, inputrec);
767     if (debug)
768     {
769         fprintf(debug, "after bc_inputrec\n");
770     }
771     bc_symtab(cr, &mtop->symtab);
772     if (debug)
773     {
774         fprintf(debug, "after bc_symtab\n");
775     }
776     bc_string(cr, &mtop->symtab, &mtop->name);
777     if (debug)
778     {
779         fprintf(debug, "after bc_name\n");
780     }
781
782     bc_ffparams(cr, &mtop->ffparams);
783
784     block_bc(cr, mtop->nmoltype);
785     snew_bc(cr, mtop->moltype, mtop->nmoltype);
786     for (i = 0; i < mtop->nmoltype; i++)
787     {
788         bc_moltype(cr, &mtop->symtab, &mtop->moltype[i]);
789     }
790
791     block_bc(cr, mtop->nmolblock);
792     snew_bc(cr, mtop->molblock, mtop->nmolblock);
793     for (i = 0; i < mtop->nmolblock; i++)
794     {
795         bc_molblock(cr, &mtop->molblock[i]);
796     }
797
798     block_bc(cr, mtop->natoms);
799
800     bc_atomtypes(cr, &mtop->atomtypes);
801
802     bc_block(cr, &mtop->mols);
803     bc_groups(cr, &mtop->symtab, mtop->natoms, &mtop->groups);
804 }