Merge release-2019 into master
[gromacs.git] / src / gromacs / mdlib / broadcaststructs.cpp
blob6eb325684e794c5ad1aaa75e125f169bafc5a5b6
1 /*
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,2016,2017,2018,2019, 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! */
38 #include "gmxpre.h"
40 #include "broadcaststructs.h"
42 #include <cstring>
44 #include <memory>
46 #include "gromacs/gmxlib/network.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/mdlib/tgroup.h"
49 #include "gromacs/mdtypes/awh_params.h"
50 #include "gromacs/mdtypes/commrec.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/mdtypes/pull_params.h"
54 #include "gromacs/mdtypes/state.h"
55 #include "gromacs/topology/mtop_util.h"
56 #include "gromacs/topology/symtab.h"
57 #include "gromacs/topology/topology.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/inmemoryserializer.h"
60 #include "gromacs/utility/keyvaluetree.h"
61 #include "gromacs/utility/keyvaluetreeserializer.h"
62 #include "gromacs/utility/smalloc.h"
64 static void bc_cstring(const t_commrec *cr, char **s)
66 int size = 0;
68 if (MASTER(cr) && *s != nullptr)
70 /* Size of the char buffer is string length + 1 for '\0' */
71 size = strlen(*s) + 1;
73 block_bc(cr, size);
74 if (size > 0)
76 if (!MASTER(cr))
78 srenew(*s, size);
80 nblock_bc(cr, size, *s);
82 else if (!MASTER(cr) && *s != nullptr)
84 sfree(*s);
85 *s = nullptr;
89 static void bc_string(const t_commrec *cr, t_symtab *symtab, char ***s)
91 int handle;
93 if (MASTER(cr))
95 handle = lookup_symtab(symtab, *s);
97 block_bc(cr, handle);
98 if (!MASTER(cr))
100 *s = get_symtab_handle(symtab, handle);
104 static void bc_strings(const t_commrec *cr, t_symtab *symtab, int nr, char ****nm)
106 int i;
107 int *handle;
109 snew(handle, nr);
110 if (MASTER(cr))
112 for (i = 0; (i < nr); i++)
114 handle[i] = lookup_symtab(symtab, (*nm)[i]);
117 nblock_bc(cr, nr, handle);
119 if (!MASTER(cr))
121 snew_bc(cr, *nm, nr);
122 for (i = 0; (i < nr); i++)
124 (*nm)[i] = get_symtab_handle(symtab, handle[i]);
127 sfree(handle);
130 static void bc_strings_container(const t_commrec *cr,
131 t_symtab *symtab,
132 int nr,
133 std::vector<char **> *nm)
135 std::vector<int> handle;
136 if (MASTER(cr))
138 for (int i = 0; (i < nr); i++)
140 handle.emplace_back(lookup_symtab(symtab, (*nm)[i]));
143 block_bc(cr, nr);
144 nblock_abc(cr, nr, &handle);
146 if (!MASTER(cr))
148 nm->resize(nr);
149 for (int i = 0; (i < nr); i++)
151 (*nm)[i] = get_symtab_handle(symtab, handle[i]);
156 static void bc_strings_resinfo(const t_commrec *cr, t_symtab *symtab,
157 int nr, t_resinfo *resinfo)
159 int i;
160 int *handle;
162 snew(handle, nr);
163 if (MASTER(cr))
165 for (i = 0; (i < nr); i++)
167 handle[i] = lookup_symtab(symtab, resinfo[i].name);
170 nblock_bc(cr, nr, handle);
172 if (!MASTER(cr))
174 for (i = 0; (i < nr); i++)
176 resinfo[i].name = get_symtab_handle(symtab, handle[i]);
179 sfree(handle);
182 static void bc_symtab(const t_commrec *cr, t_symtab *symtab)
184 int i, nr, len;
185 t_symbuf *symbuf;
187 block_bc(cr, symtab->nr);
188 nr = symtab->nr;
189 snew_bc(cr, symtab->symbuf, 1);
190 symbuf = symtab->symbuf;
191 symbuf->bufsize = nr;
192 snew_bc(cr, symbuf->buf, nr);
193 for (i = 0; i < nr; i++)
195 if (MASTER(cr))
197 len = strlen(symbuf->buf[i]) + 1;
199 block_bc(cr, len);
200 snew_bc(cr, symbuf->buf[i], len);
201 nblock_bc(cr, len, symbuf->buf[i]);
205 static void bc_block(const t_commrec *cr, t_block *block)
207 block_bc(cr, block->nr);
208 snew_bc(cr, block->index, block->nr+1);
209 nblock_bc(cr, block->nr+1, block->index);
212 static void bc_blocka(const t_commrec *cr, t_blocka *block)
214 block_bc(cr, block->nr);
215 snew_bc(cr, block->index, block->nr+1);
216 nblock_bc(cr, block->nr+1, block->index);
217 block_bc(cr, block->nra);
218 if (block->nra)
220 snew_bc(cr, block->a, block->nra);
221 nblock_bc(cr, block->nra, block->a);
225 static void bc_grps(const t_commrec *cr, gmx::ArrayRef<AtomGroupIndices> grps)
227 for (auto &group : grps)
229 int size = group.size();
230 block_bc(cr, size);
231 nblock_abc(cr, size, &group);
235 static void bc_atoms(const t_commrec *cr, t_symtab *symtab, t_atoms *atoms)
237 block_bc(cr, atoms->nr);
238 snew_bc(cr, atoms->atom, atoms->nr);
239 nblock_bc(cr, atoms->nr, atoms->atom);
240 bc_strings(cr, symtab, atoms->nr, &atoms->atomname);
241 block_bc(cr, atoms->nres);
242 snew_bc(cr, atoms->resinfo, atoms->nres);
243 nblock_bc(cr, atoms->nres, atoms->resinfo);
244 bc_strings_resinfo(cr, symtab, atoms->nres, atoms->resinfo);
245 /* QMMM requires atomtypes to be known on all nodes as well */
246 bc_strings(cr, symtab, atoms->nr, &atoms->atomtype);
247 bc_strings(cr, symtab, atoms->nr, &atoms->atomtypeB);
250 static void bc_groups(const t_commrec *cr, t_symtab *symtab,
251 int natoms, SimulationGroups *groups)
253 int n;
255 bc_grps(cr, groups->groups);
256 bc_strings_container(cr, symtab, groups->groupNames.size(), &groups->groupNames);
257 for (auto group : gmx::keysOf(groups->groups))
259 if (MASTER(cr))
261 if (!groups->groupNumbers[group].empty())
263 n = natoms;
265 else
267 n = 0;
270 block_bc(cr, n);
271 if (n != 0)
273 nblock_abc(cr, n, &groups->groupNumbers[group]);
276 if (debug)
278 fprintf(debug, "after bc_groups\n");
282 template <typename AllocatorType>
283 static void bcastPaddedRVecVector(const t_commrec *cr, gmx::PaddedVector<gmx::RVec, AllocatorType> *v, int numAtoms)
285 v->resizeWithPadding(numAtoms);
286 nblock_bc(cr, makeArrayRef(*v));
289 void broadcastStateWithoutDynamics(const t_commrec *cr, t_state *state)
291 GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr), "broadcastStateWithoutDynamics should only be used for special cases without domain decomposition");
293 if (!PAR(cr))
295 return;
298 /* Broadcasts the state sizes and flags from the master to all ranks
299 * in cr->mpi_comm_mygroup.
301 block_bc(cr, state->natoms);
302 block_bc(cr, state->flags);
304 for (int i = 0; i < estNR; i++)
306 if (state->flags & (1 << i))
308 switch (i)
310 case estLAMBDA:
311 nblock_bc(cr, efptNR, state->lambda.data());
312 break;
313 case estFEPSTATE:
314 block_bc(cr, state->fep_state);
315 break;
316 case estBOX:
317 block_bc(cr, state->box);
318 break;
319 case estX:
320 bcastPaddedRVecVector(cr, &state->x, state->natoms);
321 break;
322 default:
323 GMX_RELEASE_ASSERT(false, "The state has a dynamic entry, while no dynamic entries should be present");
324 break;
330 static void bc_ilists(const t_commrec *cr, InteractionLists *ilist)
332 int ftype;
334 /* Here we only communicate the non-zero length ilists */
335 if (MASTER(cr))
337 for (ftype = 0; ftype < F_NRE; ftype++)
339 if ((*ilist)[ftype].size() > 0)
341 block_bc(cr, ftype);
342 int nr = (*ilist)[ftype].size();
343 block_bc(cr, nr);
344 nblock_bc(cr, nr, (*ilist)[ftype].iatoms.data());
347 ftype = -1;
348 block_bc(cr, ftype);
350 else
352 for (ftype = 0; ftype < F_NRE; ftype++)
354 (*ilist)[ftype].iatoms.clear();
358 block_bc(cr, ftype);
359 if (ftype >= 0)
361 int nr;
362 block_bc(cr, nr);
363 (*ilist)[ftype].iatoms.resize(nr);
364 nblock_bc(cr, nr, (*ilist)[ftype].iatoms.data());
367 while (ftype >= 0);
370 if (debug)
372 fprintf(debug, "after bc_ilists\n");
376 static void bc_cmap(const t_commrec *cr, gmx_cmap_t *cmap_grid)
378 int ngrid = cmap_grid->cmapdata.size();
379 block_bc(cr, ngrid);
380 block_bc(cr, cmap_grid->grid_spacing);
382 int nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
384 if (ngrid > 0)
386 if (!MASTER(cr))
388 cmap_grid->cmapdata.resize(ngrid);
391 for (int i = 0; i < ngrid; i++)
393 nblock_abc(cr, 4*nelem, &cmap_grid->cmapdata[i].cmap);
398 static void bc_ffparams(const t_commrec *cr, gmx_ffparams_t *ffp)
400 int numTypes = ffp->numTypes();
401 block_bc(cr, numTypes);
402 block_bc(cr, ffp->atnr);
403 nblock_abc(cr, numTypes, &ffp->functype);
404 nblock_abc(cr, numTypes, &ffp->iparams);
405 block_bc(cr, ffp->reppow);
406 block_bc(cr, ffp->fudgeQQ);
407 bc_cmap(cr, &ffp->cmap_grid);
410 static void bc_grpopts(const t_commrec *cr, t_grpopts *g)
412 int i, n;
414 block_bc(cr, g->ngtc);
415 block_bc(cr, g->ngacc);
416 block_bc(cr, g->ngfrz);
417 block_bc(cr, g->ngener);
418 snew_bc(cr, g->nrdf, g->ngtc);
419 snew_bc(cr, g->tau_t, g->ngtc);
420 snew_bc(cr, g->ref_t, g->ngtc);
421 snew_bc(cr, g->acc, g->ngacc);
422 snew_bc(cr, g->nFreeze, g->ngfrz);
423 snew_bc(cr, g->egp_flags, g->ngener*g->ngener);
425 nblock_bc(cr, g->ngtc, g->nrdf);
426 nblock_bc(cr, g->ngtc, g->tau_t);
427 nblock_bc(cr, g->ngtc, g->ref_t);
428 nblock_bc(cr, g->ngacc, g->acc);
429 nblock_bc(cr, g->ngfrz, g->nFreeze);
430 nblock_bc(cr, g->ngener*g->ngener, g->egp_flags);
431 snew_bc(cr, g->annealing, g->ngtc);
432 snew_bc(cr, g->anneal_npoints, g->ngtc);
433 snew_bc(cr, g->anneal_time, g->ngtc);
434 snew_bc(cr, g->anneal_temp, g->ngtc);
435 nblock_bc(cr, g->ngtc, g->annealing);
436 nblock_bc(cr, g->ngtc, g->anneal_npoints);
437 for (i = 0; (i < g->ngtc); i++)
439 n = g->anneal_npoints[i];
440 if (n > 0)
442 snew_bc(cr, g->anneal_time[i], n);
443 snew_bc(cr, g->anneal_temp[i], n);
444 nblock_bc(cr, n, g->anneal_time[i]);
445 nblock_bc(cr, n, g->anneal_temp[i]);
449 /* QMMM stuff, see inputrec */
450 block_bc(cr, g->ngQM);
451 snew_bc(cr, g->QMmethod, g->ngQM);
452 snew_bc(cr, g->QMbasis, g->ngQM);
453 snew_bc(cr, g->QMcharge, g->ngQM);
454 snew_bc(cr, g->QMmult, g->ngQM);
455 snew_bc(cr, g->bSH, g->ngQM);
456 snew_bc(cr, g->CASorbitals, g->ngQM);
457 snew_bc(cr, g->CASelectrons, g->ngQM);
458 snew_bc(cr, g->SAon, g->ngQM);
459 snew_bc(cr, g->SAoff, g->ngQM);
460 snew_bc(cr, g->SAsteps, g->ngQM);
462 if (g->ngQM)
464 nblock_bc(cr, g->ngQM, g->QMmethod);
465 nblock_bc(cr, g->ngQM, g->QMbasis);
466 nblock_bc(cr, g->ngQM, g->QMcharge);
467 nblock_bc(cr, g->ngQM, g->QMmult);
468 nblock_bc(cr, g->ngQM, g->bSH);
469 nblock_bc(cr, g->ngQM, g->CASorbitals);
470 nblock_bc(cr, g->ngQM, g->CASelectrons);
471 nblock_bc(cr, g->ngQM, g->SAon);
472 nblock_bc(cr, g->ngQM, g->SAoff);
473 nblock_bc(cr, g->ngQM, g->SAsteps);
474 /* end of QMMM stuff */
478 static void bc_awhBias(const t_commrec *cr, gmx::AwhBiasParams *awhBiasParams)
480 block_bc(cr, *awhBiasParams);
482 snew_bc(cr, awhBiasParams->dimParams, awhBiasParams->ndim);
483 nblock_bc(cr, awhBiasParams->ndim, awhBiasParams->dimParams);
486 static void bc_awh(const t_commrec *cr, gmx::AwhParams *awhParams)
488 int k;
490 block_bc(cr, *awhParams);
491 snew_bc(cr, awhParams->awhBiasParams, awhParams->numBias);
492 for (k = 0; k < awhParams->numBias; k++)
494 bc_awhBias(cr, &awhParams->awhBiasParams[k]);
498 static void bc_pull_group(const t_commrec *cr, t_pull_group *pgrp)
500 block_bc(cr, *pgrp);
501 if (pgrp->nat > 0)
503 snew_bc(cr, pgrp->ind, pgrp->nat);
504 nblock_bc(cr, pgrp->nat, pgrp->ind);
506 if (pgrp->nweight > 0)
508 snew_bc(cr, pgrp->weight, pgrp->nweight);
509 nblock_bc(cr, pgrp->nweight, pgrp->weight);
513 static void bc_pull(const t_commrec *cr, pull_params_t *pull)
515 int g;
517 block_bc(cr, *pull);
518 snew_bc(cr, pull->group, pull->ngroup);
519 for (g = 0; g < pull->ngroup; g++)
521 bc_pull_group(cr, &pull->group[g]);
523 snew_bc(cr, pull->coord, pull->ncoord);
524 nblock_bc(cr, pull->ncoord, pull->coord);
525 for (int c = 0; c < pull->ncoord; c++)
527 if (!MASTER(cr))
529 pull->coord[c].externalPotentialProvider = nullptr;
531 if (pull->coord[c].eType == epullEXTERNAL)
533 bc_cstring(cr, &pull->coord[c].externalPotentialProvider);
538 static void bc_rotgrp(const t_commrec *cr, t_rotgrp *rotg)
540 block_bc(cr, *rotg);
541 if (rotg->nat > 0)
543 snew_bc(cr, rotg->ind, rotg->nat);
544 nblock_bc(cr, rotg->nat, rotg->ind);
545 snew_bc(cr, rotg->x_ref, rotg->nat);
546 nblock_bc(cr, rotg->nat, rotg->x_ref);
550 static void bc_rot(const t_commrec *cr, t_rot *rot)
552 int g;
554 block_bc(cr, *rot);
555 snew_bc(cr, rot->grp, rot->ngrp);
556 for (g = 0; g < rot->ngrp; g++)
558 bc_rotgrp(cr, &rot->grp[g]);
562 static void bc_imd(const t_commrec *cr, t_IMD *imd)
564 block_bc(cr, *imd);
565 snew_bc(cr, imd->ind, imd->nat);
566 nblock_bc(cr, imd->nat, imd->ind);
569 static void bc_fepvals(const t_commrec *cr, t_lambda *fep)
571 int i;
573 block_bc(cr, fep->nstdhdl);
574 block_bc(cr, fep->init_lambda);
575 block_bc(cr, fep->init_fep_state);
576 block_bc(cr, fep->delta_lambda);
577 block_bc(cr, fep->edHdLPrintEnergy);
578 block_bc(cr, fep->n_lambda);
579 if (fep->n_lambda > 0)
581 snew_bc(cr, fep->all_lambda, efptNR);
582 nblock_bc(cr, efptNR, fep->all_lambda);
583 for (i = 0; i < efptNR; i++)
585 snew_bc(cr, fep->all_lambda[i], fep->n_lambda);
586 nblock_bc(cr, fep->n_lambda, fep->all_lambda[i]);
589 block_bc(cr, fep->sc_alpha);
590 block_bc(cr, fep->sc_power);
591 block_bc(cr, fep->sc_r_power);
592 block_bc(cr, fep->sc_sigma);
593 block_bc(cr, fep->sc_sigma_min);
594 block_bc(cr, fep->bScCoul);
595 nblock_bc(cr, efptNR, &(fep->separate_dvdl[0]));
596 block_bc(cr, fep->dhdl_derivatives);
597 block_bc(cr, fep->dh_hist_size);
598 block_bc(cr, fep->dh_hist_spacing);
599 if (debug)
601 fprintf(debug, "after bc_fepvals\n");
605 static void bc_expandedvals(const t_commrec *cr, t_expanded *expand, int n_lambda)
607 block_bc(cr, expand->nstexpanded);
608 block_bc(cr, expand->elamstats);
609 block_bc(cr, expand->elmcmove);
610 block_bc(cr, expand->elmceq);
611 block_bc(cr, expand->equil_n_at_lam);
612 block_bc(cr, expand->equil_wl_delta);
613 block_bc(cr, expand->equil_ratio);
614 block_bc(cr, expand->equil_steps);
615 block_bc(cr, expand->equil_samples);
616 block_bc(cr, expand->lmc_seed);
617 block_bc(cr, expand->minvar);
618 block_bc(cr, expand->minvar_const);
619 block_bc(cr, expand->c_range);
620 block_bc(cr, expand->bSymmetrizedTMatrix);
621 block_bc(cr, expand->nstTij);
622 block_bc(cr, expand->lmc_repeats);
623 block_bc(cr, expand->lmc_forced_nstart);
624 block_bc(cr, expand->gibbsdeltalam);
625 block_bc(cr, expand->wl_scale);
626 block_bc(cr, expand->wl_ratio);
627 block_bc(cr, expand->init_wl_delta);
628 block_bc(cr, expand->bInit_weights);
629 snew_bc(cr, expand->init_lambda_weights, n_lambda);
630 nblock_bc(cr, n_lambda, expand->init_lambda_weights);
631 block_bc(cr, expand->mc_temp);
632 if (debug)
634 fprintf(debug, "after bc_expandedvals\n");
638 static void bc_simtempvals(const t_commrec *cr, t_simtemp *simtemp, int n_lambda)
640 block_bc(cr, simtemp->simtemp_low);
641 block_bc(cr, simtemp->simtemp_high);
642 block_bc(cr, simtemp->eSimTempScale);
643 snew_bc(cr, simtemp->temperatures, n_lambda);
644 nblock_bc(cr, n_lambda, simtemp->temperatures);
645 if (debug)
647 fprintf(debug, "after bc_simtempvals\n");
652 static void bc_swapions(const t_commrec *cr, t_swapcoords *swap)
654 block_bc(cr, *swap);
656 /* Broadcast atom indices for split groups, solvent group, and for all user-defined swap groups */
657 snew_bc(cr, swap->grp, swap->ngrp);
658 for (int i = 0; i < swap->ngrp; i++)
660 t_swapGroup *g = &swap->grp[i];
662 block_bc(cr, *g);
663 snew_bc(cr, g->ind, g->nat);
664 nblock_bc(cr, g->nat, g->ind);
666 int len = 0;
667 if (MASTER(cr))
669 len = strlen(g->molname);
671 block_bc(cr, len);
672 snew_bc(cr, g->molname, len);
673 nblock_bc(cr, len, g->molname);
678 static void bc_inputrec(const t_commrec *cr, t_inputrec *inputrec)
680 // Note that this overwrites pointers in inputrec, so all pointer fields
681 // Must be initialized separately below.
682 block_bc(cr, *inputrec);
683 if (SIMMASTER(cr))
685 gmx::InMemorySerializer serializer;
686 gmx::serializeKeyValueTree(*inputrec->params, &serializer);
687 std::vector<char> buffer = serializer.finishAndGetBuffer();
688 size_t size = buffer.size();
689 block_bc(cr, size);
690 nblock_bc(cr, size, buffer.data());
692 else
694 // block_bc() above overwrites the old pointer, so set it to a
695 // reasonable value in case code below throws.
696 inputrec->params = nullptr;
697 std::vector<char> buffer;
698 size_t size;
699 block_bc(cr, size);
700 nblock_abc(cr, size, &buffer);
701 gmx::InMemoryDeserializer serializer(buffer, false);
702 inputrec->params = new gmx::KeyValueTreeObject(
703 gmx::deserializeKeyValueTree(&serializer));
706 bc_grpopts(cr, &(inputrec->opts));
708 /* even if efep is efepNO, we need to initialize to make sure that
709 * n_lambda is set to zero */
711 snew_bc(cr, inputrec->fepvals, 1);
712 if (inputrec->efep != efepNO || inputrec->bSimTemp)
714 bc_fepvals(cr, inputrec->fepvals);
716 /* need to initialize this as well because of data checked for in the logic */
717 snew_bc(cr, inputrec->expandedvals, 1);
718 if (inputrec->bExpanded)
720 bc_expandedvals(cr, inputrec->expandedvals, inputrec->fepvals->n_lambda);
722 snew_bc(cr, inputrec->simtempvals, 1);
723 if (inputrec->bSimTemp)
725 bc_simtempvals(cr, inputrec->simtempvals, inputrec->fepvals->n_lambda);
727 if (inputrec->bPull)
729 snew_bc(cr, inputrec->pull, 1);
730 bc_pull(cr, inputrec->pull);
732 if (inputrec->bDoAwh)
734 snew_bc(cr, inputrec->awhParams, 1);
735 bc_awh(cr, inputrec->awhParams);
738 if (inputrec->bRot)
740 snew_bc(cr, inputrec->rot, 1);
741 bc_rot(cr, inputrec->rot);
743 if (inputrec->bIMD)
745 snew_bc(cr, inputrec->imd, 1);
746 bc_imd(cr, inputrec->imd);
748 if (inputrec->eSwapCoords != eswapNO)
750 snew_bc(cr, inputrec->swap, 1);
751 bc_swapions(cr, inputrec->swap);
755 static void bc_moltype(const t_commrec *cr, t_symtab *symtab,
756 gmx_moltype_t *moltype)
758 bc_string(cr, symtab, &moltype->name);
759 bc_atoms(cr, symtab, &moltype->atoms);
760 if (debug)
762 fprintf(debug, "after bc_atoms\n");
765 bc_ilists(cr, &moltype->ilist);
766 bc_block(cr, &moltype->cgs);
767 bc_blocka(cr, &moltype->excls);
770 static void bc_vector_of_rvec(const t_commrec *cr, std::vector<gmx::RVec> *vec)
772 int numElements = vec->size();
773 block_bc(cr, numElements);
774 if (!MASTER(cr))
776 vec->resize(numElements);
778 if (numElements > 0)
780 nblock_bc(cr, numElements, as_rvec_array(vec->data()));
784 static void bc_molblock(const t_commrec *cr, gmx_molblock_t *molb)
786 block_bc(cr, molb->type);
787 block_bc(cr, molb->nmol);
788 bc_vector_of_rvec(cr, &molb->posres_xA);
789 bc_vector_of_rvec(cr, &molb->posres_xB);
790 if (debug)
792 fprintf(debug, "after bc_molblock\n");
796 static void bc_atomtypes(const t_commrec *cr, t_atomtypes *atomtypes)
798 block_bc(cr, atomtypes->nr);
801 /*! \brief Broadcasts ir and mtop from the master to all nodes in
802 * cr->mpi_comm_mygroup. */
803 static
804 void bcast_ir_mtop(const t_commrec *cr, t_inputrec *inputrec, gmx_mtop_t *mtop)
806 if (debug)
808 fprintf(debug, "in bc_data\n");
810 bc_inputrec(cr, inputrec);
811 if (debug)
813 fprintf(debug, "after bc_inputrec\n");
815 bc_symtab(cr, &mtop->symtab);
816 if (debug)
818 fprintf(debug, "after bc_symtab\n");
820 bc_string(cr, &mtop->symtab, &mtop->name);
821 if (debug)
823 fprintf(debug, "after bc_name\n");
826 bc_ffparams(cr, &mtop->ffparams);
828 int nmoltype = mtop->moltype.size();
829 block_bc(cr, nmoltype);
830 mtop->moltype.resize(nmoltype);
831 for (gmx_moltype_t &moltype : mtop->moltype)
833 bc_moltype(cr, &mtop->symtab, &moltype);
836 block_bc(cr, mtop->bIntermolecularInteractions);
837 if (mtop->bIntermolecularInteractions)
839 if (!MASTER(cr))
841 mtop->intermolecular_ilist = std::make_unique<InteractionLists>();
843 bc_ilists(cr, mtop->intermolecular_ilist.get());
846 int nmolblock = mtop->molblock.size();
847 block_bc(cr, nmolblock);
848 mtop->molblock.resize(nmolblock);
849 for (gmx_molblock_t &molblock : mtop->molblock)
851 bc_molblock(cr, &molblock);
854 block_bc(cr, mtop->natoms);
856 bc_atomtypes(cr, &mtop->atomtypes);
858 bc_groups(cr, &mtop->symtab, mtop->natoms, &mtop->groups);
860 GMX_RELEASE_ASSERT(!MASTER(cr) || mtop->haveMoleculeIndices, "mtop should have valid molecule indices");
861 if (!MASTER(cr))
863 mtop->haveMoleculeIndices = true;
865 gmx_mtop_finalize(mtop);
869 void init_parallel(t_commrec *cr, t_inputrec *inputrec,
870 gmx_mtop_t *mtop)
872 bcast_ir_mtop(cr, inputrec, mtop);