Fix dependencies on config.h for typedefs.h
[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 "types/commrec.h"
48 #include "network.h"
49 #include "smalloc.h"
50 #include "gmx_fatal.h"
51 #include "symtab.h"
52 #include "vec.h"
53 #include "tgroup.h"
54
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)); }
63
64 static void bc_string(const t_commrec *cr, t_symtab *symtab, char ***s)
65 {
66     int handle;
67
68     if (MASTER(cr))
69     {
70         handle = lookup_symtab(symtab, *s);
71     }
72     block_bc(cr, handle);
73     if (!MASTER(cr))
74     {
75         *s = get_symtab_handle(symtab, handle);
76     }
77 }
78
79 static void bc_strings(const t_commrec *cr, t_symtab *symtab, int nr, char ****nm)
80 {
81     int     i;
82     int    *handle;
83     char ***NM;
84
85     snew(handle, nr);
86     if (MASTER(cr))
87     {
88         NM = *nm;
89         for (i = 0; (i < nr); i++)
90         {
91             handle[i] = lookup_symtab(symtab, NM[i]);
92         }
93     }
94     nblock_bc(cr, nr, handle);
95
96     if (!MASTER(cr))
97     {
98         snew_bc(cr, *nm, nr);
99         NM = *nm;
100         for (i = 0; (i < nr); i++)
101         {
102             (*nm)[i] = get_symtab_handle(symtab, handle[i]);
103         }
104     }
105     sfree(handle);
106 }
107
108 static void bc_strings_resinfo(const t_commrec *cr, t_symtab *symtab,
109                                int nr, t_resinfo *resinfo)
110 {
111     int   i;
112     int  *handle;
113
114     snew(handle, nr);
115     if (MASTER(cr))
116     {
117         for (i = 0; (i < nr); i++)
118         {
119             handle[i] = lookup_symtab(symtab, resinfo[i].name);
120         }
121     }
122     nblock_bc(cr, nr, handle);
123
124     if (!MASTER(cr))
125     {
126         for (i = 0; (i < nr); i++)
127         {
128             resinfo[i].name = get_symtab_handle(symtab, handle[i]);
129         }
130     }
131     sfree(handle);
132 }
133
134 static void bc_symtab(const t_commrec *cr, t_symtab *symtab)
135 {
136     int       i, nr, len;
137     t_symbuf *symbuf;
138
139     block_bc(cr, symtab->nr);
140     nr = 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++)
146     {
147         if (MASTER(cr))
148         {
149             len = strlen(symbuf->buf[i]) + 1;
150         }
151         block_bc(cr, len);
152         snew_bc(cr, symbuf->buf[i], len);
153         nblock_bc(cr, len, symbuf->buf[i]);
154     }
155 }
156
157 static void bc_block(const t_commrec *cr, t_block *block)
158 {
159     block_bc(cr, block->nr);
160     snew_bc(cr, block->index, block->nr+1);
161     nblock_bc(cr, block->nr+1, block->index);
162 }
163
164 static void bc_blocka(const t_commrec *cr, t_blocka *block)
165 {
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);
170     if (block->nra)
171     {
172         snew_bc(cr, block->a, block->nra);
173         nblock_bc(cr, block->nra, block->a);
174     }
175 }
176
177 static void bc_grps(const t_commrec *cr, t_grps grps[])
178 {
179     int i;
180
181     for (i = 0; (i < egcNR); i++)
182     {
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);
186     }
187 }
188
189 static void bc_atoms(const t_commrec *cr, t_symtab *symtab, t_atoms *atoms)
190 {
191     int dummy;
192
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);
204 }
205
206 static void bc_groups(const t_commrec *cr, t_symtab *symtab,
207                       int natoms, gmx_groups_t *groups)
208 {
209     int dummy;
210     int g, n;
211
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++)
216     {
217         if (MASTER(cr))
218         {
219             if (groups->grpnr[g])
220             {
221                 n = natoms;
222             }
223             else
224             {
225                 n = 0;
226             }
227         }
228         block_bc(cr, n);
229         if (n == 0)
230         {
231             groups->grpnr[g] = NULL;
232         }
233         else
234         {
235             snew_bc(cr, groups->grpnr[g], n);
236             nblock_bc(cr, n, groups->grpnr[g]);
237         }
238     }
239     if (debug)
240     {
241         fprintf(debug, "after bc_groups\n");
242     }
243 }
244
245 void bcast_state(const t_commrec *cr, t_state *state)
246 {
247     int      i, nnht, nnhtp;
248     gmx_bool bAlloc;
249
250     if (!PAR(cr))
251     {
252         return;
253     }
254
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)
263     {
264         snew_bc(cr, state->lambda, efptNR)
265     }
266
267     if (cr->dd)
268     {
269         /* We allocate dynamically in dd_partition_system. */
270         return;
271     }
272     /* The code below is reachable only by TPI and NM, so it is not
273        tested by anything. */
274
275     nnht  = (state->ngtc)*(state->nhchainlength);
276     nnhtp = (state->nnhpres)*(state->nhchainlength);
277
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);
282     if (bAlloc)
283     {
284         state->nalloc = state->natoms;
285     }
286     for (i = 0; i < estNR; i++)
287     {
288         if (state->flags & (1<<i))
289         {
290             switch (i)
291             {
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);
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
648 static void bc_swapions(const t_commrec *cr, t_swapcoords *swap)
649 {
650     int i;
651
652
653     block_bc(cr, *swap);
654
655     /* Broadcast ion group atom indices */
656     snew_bc(cr, swap->ind, swap->nat);
657     nblock_bc(cr, swap->nat, swap->ind);
658
659     /* Broadcast split groups atom indices */
660     for (i = 0; i < 2; i++)
661     {
662         snew_bc(cr, swap->ind_split[i], swap->nat_split[i]);
663         nblock_bc(cr, swap->nat_split[i], swap->ind_split[i]);
664     }
665
666     /* Broadcast solvent group atom indices */
667     snew_bc(cr, swap->ind_sol, swap->nat_sol);
668     nblock_bc(cr, swap->nat_sol, swap->ind_sol);
669 }
670
671
672 static void bc_inputrec(const t_commrec *cr, t_inputrec *inputrec)
673 {
674     gmx_bool bAlloc = TRUE;
675     int      i;
676
677     block_bc(cr, *inputrec);
678
679     bc_grpopts(cr, &(inputrec->opts));
680
681     /* even if efep is efepNO, we need to initialize to make sure that
682      * n_lambda is set to zero */
683
684     snew_bc(cr, inputrec->fepvals, 1);
685     if (inputrec->efep != efepNO || inputrec->bSimTemp)
686     {
687         bc_fepvals(cr, inputrec->fepvals);
688     }
689     /* need to initialize this as well because of data checked for in the logic */
690     snew_bc(cr, inputrec->expandedvals, 1);
691     if (inputrec->bExpanded)
692     {
693         bc_expandedvals(cr, inputrec->expandedvals, inputrec->fepvals->n_lambda);
694     }
695     snew_bc(cr, inputrec->simtempvals, 1);
696     if (inputrec->bSimTemp)
697     {
698         bc_simtempvals(cr, inputrec->simtempvals, inputrec->fepvals->n_lambda);
699     }
700     if (inputrec->ePull != epullNO)
701     {
702         snew_bc(cr, inputrec->pull, 1);
703         bc_pull(cr, inputrec->pull);
704     }
705     if (inputrec->bRot)
706     {
707         snew_bc(cr, inputrec->rot, 1);
708         bc_rot(cr, inputrec->rot);
709     }
710     for (i = 0; (i < DIM); i++)
711     {
712         bc_cosines(cr, &(inputrec->ex[i]));
713         bc_cosines(cr, &(inputrec->et[i]));
714     }
715     if (inputrec->eSwapCoords != eswapNO)
716     {
717         snew_bc(cr, inputrec->swap, 1);
718         bc_swapions(cr, inputrec->swap);
719     }
720     if (inputrec->bAdress)
721     {
722         snew_bc(cr, inputrec->adress, 1);
723         bc_adress(cr, inputrec->adress);
724     }
725 }
726
727 static void bc_moltype(const t_commrec *cr, t_symtab *symtab,
728                        gmx_moltype_t *moltype)
729 {
730     bc_string(cr, symtab, &moltype->name);
731     bc_atoms(cr, symtab, &moltype->atoms);
732     if (debug)
733     {
734         fprintf(debug, "after bc_atoms\n");
735     }
736
737     bc_ilists(cr, moltype->ilist);
738     bc_block(cr, &moltype->cgs);
739     bc_blocka(cr, &moltype->excls);
740 }
741
742 static void bc_molblock(const t_commrec *cr, gmx_molblock_t *molb)
743 {
744     gmx_bool bAlloc = TRUE;
745
746     block_bc(cr, molb->type);
747     block_bc(cr, molb->nmol);
748     block_bc(cr, molb->natoms_mol);
749     block_bc(cr, molb->nposres_xA);
750     if (molb->nposres_xA > 0)
751     {
752         snew_bc(cr, molb->posres_xA, molb->nposres_xA);
753         nblock_bc(cr, molb->nposres_xA*DIM, molb->posres_xA[0]);
754     }
755     block_bc(cr, molb->nposres_xB);
756     if (molb->nposres_xB > 0)
757     {
758         snew_bc(cr, molb->posres_xB, molb->nposres_xB);
759         nblock_bc(cr, molb->nposres_xB*DIM, molb->posres_xB[0]);
760     }
761     if (debug)
762     {
763         fprintf(debug, "after bc_molblock\n");
764     }
765 }
766
767 static void bc_atomtypes(const t_commrec *cr, t_atomtypes *atomtypes)
768 {
769     int nr;
770
771     block_bc(cr, atomtypes->nr);
772
773     nr = atomtypes->nr;
774
775     snew_bc(cr, atomtypes->radius, nr);
776     snew_bc(cr, atomtypes->vol, nr);
777     snew_bc(cr, atomtypes->surftens, nr);
778     snew_bc(cr, atomtypes->gb_radius, nr);
779     snew_bc(cr, atomtypes->S_hct, nr);
780
781     nblock_bc(cr, nr, atomtypes->radius);
782     nblock_bc(cr, nr, atomtypes->vol);
783     nblock_bc(cr, nr, atomtypes->surftens);
784     nblock_bc(cr, nr, atomtypes->gb_radius);
785     nblock_bc(cr, nr, atomtypes->S_hct);
786 }
787
788
789 void bcast_ir_mtop(const t_commrec *cr, t_inputrec *inputrec, gmx_mtop_t *mtop)
790 {
791     int i;
792     if (debug)
793     {
794         fprintf(debug, "in bc_data\n");
795     }
796     bc_inputrec(cr, inputrec);
797     if (debug)
798     {
799         fprintf(debug, "after bc_inputrec\n");
800     }
801     bc_symtab(cr, &mtop->symtab);
802     if (debug)
803     {
804         fprintf(debug, "after bc_symtab\n");
805     }
806     bc_string(cr, &mtop->symtab, &mtop->name);
807     if (debug)
808     {
809         fprintf(debug, "after bc_name\n");
810     }
811
812     bc_ffparams(cr, &mtop->ffparams);
813
814     block_bc(cr, mtop->nmoltype);
815     snew_bc(cr, mtop->moltype, mtop->nmoltype);
816     for (i = 0; i < mtop->nmoltype; i++)
817     {
818         bc_moltype(cr, &mtop->symtab, &mtop->moltype[i]);
819     }
820
821     block_bc(cr, mtop->nmolblock);
822     snew_bc(cr, mtop->molblock, mtop->nmolblock);
823     for (i = 0; i < mtop->nmolblock; i++)
824     {
825         bc_molblock(cr, &mtop->molblock[i]);
826     }
827
828     block_bc(cr, mtop->natoms);
829
830     bc_atomtypes(cr, &mtop->atomtypes);
831
832     bc_block(cr, &mtop->mols);
833     bc_groups(cr, &mtop->symtab, mtop->natoms, &mtop->groups);
834 }