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, 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! */
40 #include "broadcaststructs.h"
44 #include "gromacs/gmxlib/network.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/mdlib/mdrun.h"
47 #include "gromacs/mdlib/tgroup.h"
48 #include "gromacs/mdtypes/awh-params.h"
49 #include "gromacs/mdtypes/commrec.h"
50 #include "gromacs/mdtypes/inputrec.h"
51 #include "gromacs/mdtypes/md_enums.h"
52 #include "gromacs/mdtypes/pull-params.h"
53 #include "gromacs/mdtypes/state.h"
54 #include "gromacs/topology/symtab.h"
55 #include "gromacs/topology/topology.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/inmemoryserializer.h"
58 #include "gromacs/utility/keyvaluetree.h"
59 #include "gromacs/utility/keyvaluetreeserializer.h"
60 #include "gromacs/utility/smalloc.h"
62 static void bc_cstring(const t_commrec
*cr
, char **s
)
66 if (MASTER(cr
) && *s
!= nullptr)
68 /* Size of the char buffer is string length + 1 for '\0' */
69 size
= strlen(*s
) + 1;
78 nblock_bc(cr
, size
, *s
);
80 else if (!MASTER(cr
) && *s
!= nullptr)
87 static void bc_string(const t_commrec
*cr
, t_symtab
*symtab
, char ***s
)
93 handle
= lookup_symtab(symtab
, *s
);
98 *s
= get_symtab_handle(symtab
, handle
);
102 static void bc_strings(const t_commrec
*cr
, t_symtab
*symtab
, int nr
, char ****nm
)
110 for (i
= 0; (i
< nr
); i
++)
112 handle
[i
] = lookup_symtab(symtab
, (*nm
)[i
]);
115 nblock_bc(cr
, nr
, handle
);
119 snew_bc(cr
, *nm
, nr
);
120 for (i
= 0; (i
< nr
); i
++)
122 (*nm
)[i
] = get_symtab_handle(symtab
, handle
[i
]);
128 static void bc_strings_resinfo(const t_commrec
*cr
, t_symtab
*symtab
,
129 int nr
, t_resinfo
*resinfo
)
137 for (i
= 0; (i
< nr
); i
++)
139 handle
[i
] = lookup_symtab(symtab
, resinfo
[i
].name
);
142 nblock_bc(cr
, nr
, handle
);
146 for (i
= 0; (i
< nr
); i
++)
148 resinfo
[i
].name
= get_symtab_handle(symtab
, handle
[i
]);
154 static void bc_symtab(const t_commrec
*cr
, t_symtab
*symtab
)
159 block_bc(cr
, symtab
->nr
);
161 snew_bc(cr
, symtab
->symbuf
, 1);
162 symbuf
= symtab
->symbuf
;
163 symbuf
->bufsize
= nr
;
164 snew_bc(cr
, symbuf
->buf
, nr
);
165 for (i
= 0; i
< nr
; i
++)
169 len
= strlen(symbuf
->buf
[i
]) + 1;
172 snew_bc(cr
, symbuf
->buf
[i
], len
);
173 nblock_bc(cr
, len
, symbuf
->buf
[i
]);
177 static void bc_block(const t_commrec
*cr
, t_block
*block
)
179 block_bc(cr
, block
->nr
);
180 snew_bc(cr
, block
->index
, block
->nr
+1);
181 nblock_bc(cr
, block
->nr
+1, block
->index
);
184 static void bc_blocka(const t_commrec
*cr
, t_blocka
*block
)
186 block_bc(cr
, block
->nr
);
187 snew_bc(cr
, block
->index
, block
->nr
+1);
188 nblock_bc(cr
, block
->nr
+1, block
->index
);
189 block_bc(cr
, block
->nra
);
192 snew_bc(cr
, block
->a
, block
->nra
);
193 nblock_bc(cr
, block
->nra
, block
->a
);
197 static void bc_grps(const t_commrec
*cr
, t_grps grps
[])
201 for (i
= 0; (i
< egcNR
); i
++)
203 block_bc(cr
, grps
[i
].nr
);
204 snew_bc(cr
, grps
[i
].nm_ind
, grps
[i
].nr
);
205 nblock_bc(cr
, grps
[i
].nr
, grps
[i
].nm_ind
);
209 static void bc_atoms(const t_commrec
*cr
, t_symtab
*symtab
, t_atoms
*atoms
)
211 block_bc(cr
, atoms
->nr
);
212 snew_bc(cr
, atoms
->atom
, atoms
->nr
);
213 nblock_bc(cr
, atoms
->nr
, atoms
->atom
);
214 bc_strings(cr
, symtab
, atoms
->nr
, &atoms
->atomname
);
215 block_bc(cr
, atoms
->nres
);
216 snew_bc(cr
, atoms
->resinfo
, atoms
->nres
);
217 nblock_bc(cr
, atoms
->nres
, atoms
->resinfo
);
218 bc_strings_resinfo(cr
, symtab
, atoms
->nres
, atoms
->resinfo
);
219 /* QMMM requires atomtypes to be known on all nodes as well */
220 bc_strings(cr
, symtab
, atoms
->nr
, &atoms
->atomtype
);
221 bc_strings(cr
, symtab
, atoms
->nr
, &atoms
->atomtypeB
);
224 static void bc_groups(const t_commrec
*cr
, t_symtab
*symtab
,
225 int natoms
, gmx_groups_t
*groups
)
229 bc_grps(cr
, groups
->grps
);
230 block_bc(cr
, groups
->ngrpname
);
231 bc_strings(cr
, symtab
, groups
->ngrpname
, &groups
->grpname
);
232 for (g
= 0; g
< egcNR
; g
++)
236 if (groups
->grpnr
[g
])
248 groups
->grpnr
[g
] = nullptr;
252 snew_bc(cr
, groups
->grpnr
[g
], n
);
253 nblock_bc(cr
, n
, groups
->grpnr
[g
]);
258 fprintf(debug
, "after bc_groups\n");
262 static void bcastPaddedRVecVector(const t_commrec
*cr
, PaddedRVecVector
*v
, int numAtoms
)
264 v
->resize(gmx::paddedRVecVectorSize(numAtoms
));
265 nblock_bc(cr
, numAtoms
, as_rvec_array(v
->data()));
268 void broadcastStateWithoutDynamics(const t_commrec
*cr
, t_state
*state
)
270 GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr
), "broadcastStateWithoutDynamics should only be used for special cases without domain decomposition");
277 /* Broadcasts the state sizes and flags from the master to all ranks
278 * in cr->mpi_comm_mygroup.
280 block_bc(cr
, state
->natoms
);
281 block_bc(cr
, state
->flags
);
283 for (int i
= 0; i
< estNR
; i
++)
285 if (state
->flags
& (1 << i
))
290 nblock_bc(cr
, efptNR
, state
->lambda
.data());
293 block_bc(cr
, state
->fep_state
);
296 block_bc(cr
, state
->box
);
299 bcastPaddedRVecVector(cr
, &state
->x
, state
->natoms
);
302 GMX_RELEASE_ASSERT(false, "The state has a dynamic entry, while no dynamic entries should be present");
309 static void bc_ilists(const t_commrec
*cr
, t_ilist
*ilist
)
313 /* Here we only communicate the non-zero length ilists */
316 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
318 if (ilist
[ftype
].nr
> 0)
321 block_bc(cr
, ilist
[ftype
].nr
);
322 nblock_bc(cr
, ilist
[ftype
].nr
, ilist
[ftype
].iatoms
);
330 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
339 block_bc(cr
, ilist
[ftype
].nr
);
340 snew_bc(cr
, ilist
[ftype
].iatoms
, ilist
[ftype
].nr
);
341 nblock_bc(cr
, ilist
[ftype
].nr
, ilist
[ftype
].iatoms
);
349 fprintf(debug
, "after bc_ilists\n");
353 static void bc_cmap(const t_commrec
*cr
, gmx_cmap_t
*cmap_grid
)
357 block_bc(cr
, cmap_grid
->ngrid
);
358 block_bc(cr
, cmap_grid
->grid_spacing
);
360 ngrid
= cmap_grid
->ngrid
;
361 nelem
= cmap_grid
->grid_spacing
* cmap_grid
->grid_spacing
;
365 snew_bc(cr
, cmap_grid
->cmapdata
, ngrid
);
367 for (i
= 0; i
< ngrid
; i
++)
369 snew_bc(cr
, cmap_grid
->cmapdata
[i
].cmap
, 4*nelem
);
370 nblock_bc(cr
, 4*nelem
, cmap_grid
->cmapdata
[i
].cmap
);
375 static void bc_ffparams(const t_commrec
*cr
, gmx_ffparams_t
*ffp
)
377 block_bc(cr
, ffp
->ntypes
);
378 block_bc(cr
, ffp
->atnr
);
379 snew_bc(cr
, ffp
->functype
, ffp
->ntypes
);
380 snew_bc(cr
, ffp
->iparams
, ffp
->ntypes
);
381 nblock_bc(cr
, ffp
->ntypes
, ffp
->functype
);
382 nblock_bc(cr
, ffp
->ntypes
, ffp
->iparams
);
383 block_bc(cr
, ffp
->reppow
);
384 block_bc(cr
, ffp
->fudgeQQ
);
385 bc_cmap(cr
, &ffp
->cmap_grid
);
388 static void bc_grpopts(const t_commrec
*cr
, t_grpopts
*g
)
392 block_bc(cr
, g
->ngtc
);
393 block_bc(cr
, g
->ngacc
);
394 block_bc(cr
, g
->ngfrz
);
395 block_bc(cr
, g
->ngener
);
396 snew_bc(cr
, g
->nrdf
, g
->ngtc
);
397 snew_bc(cr
, g
->tau_t
, g
->ngtc
);
398 snew_bc(cr
, g
->ref_t
, g
->ngtc
);
399 snew_bc(cr
, g
->acc
, g
->ngacc
);
400 snew_bc(cr
, g
->nFreeze
, g
->ngfrz
);
401 snew_bc(cr
, g
->egp_flags
, g
->ngener
*g
->ngener
);
403 nblock_bc(cr
, g
->ngtc
, g
->nrdf
);
404 nblock_bc(cr
, g
->ngtc
, g
->tau_t
);
405 nblock_bc(cr
, g
->ngtc
, g
->ref_t
);
406 nblock_bc(cr
, g
->ngacc
, g
->acc
);
407 nblock_bc(cr
, g
->ngfrz
, g
->nFreeze
);
408 nblock_bc(cr
, g
->ngener
*g
->ngener
, g
->egp_flags
);
409 snew_bc(cr
, g
->annealing
, g
->ngtc
);
410 snew_bc(cr
, g
->anneal_npoints
, g
->ngtc
);
411 snew_bc(cr
, g
->anneal_time
, g
->ngtc
);
412 snew_bc(cr
, g
->anneal_temp
, g
->ngtc
);
413 nblock_bc(cr
, g
->ngtc
, g
->annealing
);
414 nblock_bc(cr
, g
->ngtc
, g
->anneal_npoints
);
415 for (i
= 0; (i
< g
->ngtc
); i
++)
417 n
= g
->anneal_npoints
[i
];
420 snew_bc(cr
, g
->anneal_time
[i
], n
);
421 snew_bc(cr
, g
->anneal_temp
[i
], n
);
422 nblock_bc(cr
, n
, g
->anneal_time
[i
]);
423 nblock_bc(cr
, n
, g
->anneal_temp
[i
]);
427 /* QMMM stuff, see inputrec */
428 block_bc(cr
, g
->ngQM
);
429 snew_bc(cr
, g
->QMmethod
, g
->ngQM
);
430 snew_bc(cr
, g
->QMbasis
, g
->ngQM
);
431 snew_bc(cr
, g
->QMcharge
, g
->ngQM
);
432 snew_bc(cr
, g
->QMmult
, g
->ngQM
);
433 snew_bc(cr
, g
->bSH
, g
->ngQM
);
434 snew_bc(cr
, g
->CASorbitals
, g
->ngQM
);
435 snew_bc(cr
, g
->CASelectrons
, g
->ngQM
);
436 snew_bc(cr
, g
->SAon
, g
->ngQM
);
437 snew_bc(cr
, g
->SAoff
, g
->ngQM
);
438 snew_bc(cr
, g
->SAsteps
, g
->ngQM
);
442 nblock_bc(cr
, g
->ngQM
, g
->QMmethod
);
443 nblock_bc(cr
, g
->ngQM
, g
->QMbasis
);
444 nblock_bc(cr
, g
->ngQM
, g
->QMcharge
);
445 nblock_bc(cr
, g
->ngQM
, g
->QMmult
);
446 nblock_bc(cr
, g
->ngQM
, g
->bSH
);
447 nblock_bc(cr
, g
->ngQM
, g
->CASorbitals
);
448 nblock_bc(cr
, g
->ngQM
, g
->CASelectrons
);
449 nblock_bc(cr
, g
->ngQM
, g
->SAon
);
450 nblock_bc(cr
, g
->ngQM
, g
->SAoff
);
451 nblock_bc(cr
, g
->ngQM
, g
->SAsteps
);
452 /* end of QMMM stuff */
456 static void bc_awhBias(const t_commrec
*cr
, gmx::AwhBiasParams
*awhBiasParams
)
458 block_bc(cr
, *awhBiasParams
);
460 snew_bc(cr
, awhBiasParams
->dimParams
, awhBiasParams
->ndim
);
461 nblock_bc(cr
, awhBiasParams
->ndim
, awhBiasParams
->dimParams
);
464 static void bc_awh(const t_commrec
*cr
, gmx::AwhParams
*awhParams
)
468 block_bc(cr
, *awhParams
);
469 snew_bc(cr
, awhParams
->awhBiasParams
, awhParams
->numBias
);
470 for (k
= 0; k
< awhParams
->numBias
; k
++)
472 bc_awhBias(cr
, &awhParams
->awhBiasParams
[k
]);
476 static void bc_pull_group(const t_commrec
*cr
, t_pull_group
*pgrp
)
481 snew_bc(cr
, pgrp
->ind
, pgrp
->nat
);
482 nblock_bc(cr
, pgrp
->nat
, pgrp
->ind
);
484 if (pgrp
->nweight
> 0)
486 snew_bc(cr
, pgrp
->weight
, pgrp
->nweight
);
487 nblock_bc(cr
, pgrp
->nweight
, pgrp
->weight
);
491 static void bc_pull(const t_commrec
*cr
, pull_params_t
*pull
)
496 snew_bc(cr
, pull
->group
, pull
->ngroup
);
497 for (g
= 0; g
< pull
->ngroup
; g
++)
499 bc_pull_group(cr
, &pull
->group
[g
]);
501 snew_bc(cr
, pull
->coord
, pull
->ncoord
);
502 nblock_bc(cr
, pull
->ncoord
, pull
->coord
);
503 for (int c
= 0; c
< pull
->ncoord
; c
++)
507 pull
->coord
[c
].externalPotentialProvider
= nullptr;
509 if (pull
->coord
[c
].eType
== epullEXTERNAL
)
511 bc_cstring(cr
, &pull
->coord
[c
].externalPotentialProvider
);
516 static void bc_rotgrp(const t_commrec
*cr
, t_rotgrp
*rotg
)
521 snew_bc(cr
, rotg
->ind
, rotg
->nat
);
522 nblock_bc(cr
, rotg
->nat
, rotg
->ind
);
523 snew_bc(cr
, rotg
->x_ref
, rotg
->nat
);
524 nblock_bc(cr
, rotg
->nat
, rotg
->x_ref
);
528 static void bc_rot(const t_commrec
*cr
, t_rot
*rot
)
533 snew_bc(cr
, rot
->grp
, rot
->ngrp
);
534 for (g
= 0; g
< rot
->ngrp
; g
++)
536 bc_rotgrp(cr
, &rot
->grp
[g
]);
540 static void bc_imd(const t_commrec
*cr
, t_IMD
*imd
)
543 snew_bc(cr
, imd
->ind
, imd
->nat
);
544 nblock_bc(cr
, imd
->nat
, imd
->ind
);
547 static void bc_fepvals(const t_commrec
*cr
, t_lambda
*fep
)
551 block_bc(cr
, fep
->nstdhdl
);
552 block_bc(cr
, fep
->init_lambda
);
553 block_bc(cr
, fep
->init_fep_state
);
554 block_bc(cr
, fep
->delta_lambda
);
555 block_bc(cr
, fep
->edHdLPrintEnergy
);
556 block_bc(cr
, fep
->n_lambda
);
557 if (fep
->n_lambda
> 0)
559 snew_bc(cr
, fep
->all_lambda
, efptNR
);
560 nblock_bc(cr
, efptNR
, fep
->all_lambda
);
561 for (i
= 0; i
< efptNR
; i
++)
563 snew_bc(cr
, fep
->all_lambda
[i
], fep
->n_lambda
);
564 nblock_bc(cr
, fep
->n_lambda
, fep
->all_lambda
[i
]);
567 block_bc(cr
, fep
->sc_alpha
);
568 block_bc(cr
, fep
->sc_power
);
569 block_bc(cr
, fep
->sc_r_power
);
570 block_bc(cr
, fep
->sc_sigma
);
571 block_bc(cr
, fep
->sc_sigma_min
);
572 block_bc(cr
, fep
->bScCoul
);
573 nblock_bc(cr
, efptNR
, &(fep
->separate_dvdl
[0]));
574 block_bc(cr
, fep
->dhdl_derivatives
);
575 block_bc(cr
, fep
->dh_hist_size
);
576 block_bc(cr
, fep
->dh_hist_spacing
);
579 fprintf(debug
, "after bc_fepvals\n");
583 static void bc_expandedvals(const t_commrec
*cr
, t_expanded
*expand
, int n_lambda
)
585 block_bc(cr
, expand
->nstexpanded
);
586 block_bc(cr
, expand
->elamstats
);
587 block_bc(cr
, expand
->elmcmove
);
588 block_bc(cr
, expand
->elmceq
);
589 block_bc(cr
, expand
->equil_n_at_lam
);
590 block_bc(cr
, expand
->equil_wl_delta
);
591 block_bc(cr
, expand
->equil_ratio
);
592 block_bc(cr
, expand
->equil_steps
);
593 block_bc(cr
, expand
->equil_samples
);
594 block_bc(cr
, expand
->lmc_seed
);
595 block_bc(cr
, expand
->minvar
);
596 block_bc(cr
, expand
->minvar_const
);
597 block_bc(cr
, expand
->c_range
);
598 block_bc(cr
, expand
->bSymmetrizedTMatrix
);
599 block_bc(cr
, expand
->nstTij
);
600 block_bc(cr
, expand
->lmc_repeats
);
601 block_bc(cr
, expand
->lmc_forced_nstart
);
602 block_bc(cr
, expand
->gibbsdeltalam
);
603 block_bc(cr
, expand
->wl_scale
);
604 block_bc(cr
, expand
->wl_ratio
);
605 block_bc(cr
, expand
->init_wl_delta
);
606 block_bc(cr
, expand
->bInit_weights
);
607 snew_bc(cr
, expand
->init_lambda_weights
, n_lambda
);
608 nblock_bc(cr
, n_lambda
, expand
->init_lambda_weights
);
609 block_bc(cr
, expand
->mc_temp
);
612 fprintf(debug
, "after bc_expandedvals\n");
616 static void bc_simtempvals(const t_commrec
*cr
, t_simtemp
*simtemp
, int n_lambda
)
618 block_bc(cr
, simtemp
->simtemp_low
);
619 block_bc(cr
, simtemp
->simtemp_high
);
620 block_bc(cr
, simtemp
->eSimTempScale
);
621 snew_bc(cr
, simtemp
->temperatures
, n_lambda
);
622 nblock_bc(cr
, n_lambda
, simtemp
->temperatures
);
625 fprintf(debug
, "after bc_simtempvals\n");
630 static void bc_swapions(const t_commrec
*cr
, t_swapcoords
*swap
)
634 /* Broadcast atom indices for split groups, solvent group, and for all user-defined swap groups */
635 snew_bc(cr
, swap
->grp
, swap
->ngrp
);
636 for (int i
= 0; i
< swap
->ngrp
; i
++)
638 t_swapGroup
*g
= &swap
->grp
[i
];
641 snew_bc(cr
, g
->ind
, g
->nat
);
642 nblock_bc(cr
, g
->nat
, g
->ind
);
647 len
= strlen(g
->molname
);
650 snew_bc(cr
, g
->molname
, len
);
651 nblock_bc(cr
, len
, g
->molname
);
656 static void bc_inputrec(const t_commrec
*cr
, t_inputrec
*inputrec
)
658 // Note that this overwrites pointers in inputrec, so all pointer fields
659 // Must be initialized separately below.
660 block_bc(cr
, *inputrec
);
663 gmx::InMemorySerializer serializer
;
664 gmx::serializeKeyValueTree(*inputrec
->params
, &serializer
);
665 std::vector
<char> buffer
= serializer
.finishAndGetBuffer();
666 size_t size
= buffer
.size();
668 nblock_bc(cr
, size
, buffer
.data());
672 // block_bc() above overwrites the old pointer, so set it to a
673 // reasonable value in case code below throws.
674 // cppcheck-suppress redundantAssignment
675 inputrec
->params
= nullptr;
676 std::vector
<char> buffer
;
679 nblock_abc(cr
, size
, &buffer
);
680 gmx::InMemoryDeserializer
serializer(buffer
);
681 // cppcheck-suppress redundantAssignment
682 inputrec
->params
= new gmx::KeyValueTreeObject(
683 gmx::deserializeKeyValueTree(&serializer
));
686 bc_grpopts(cr
, &(inputrec
->opts
));
688 /* even if efep is efepNO, we need to initialize to make sure that
689 * n_lambda is set to zero */
691 snew_bc(cr
, inputrec
->fepvals
, 1);
692 if (inputrec
->efep
!= efepNO
|| inputrec
->bSimTemp
)
694 bc_fepvals(cr
, inputrec
->fepvals
);
696 /* need to initialize this as well because of data checked for in the logic */
697 snew_bc(cr
, inputrec
->expandedvals
, 1);
698 if (inputrec
->bExpanded
)
700 bc_expandedvals(cr
, inputrec
->expandedvals
, inputrec
->fepvals
->n_lambda
);
702 snew_bc(cr
, inputrec
->simtempvals
, 1);
703 if (inputrec
->bSimTemp
)
705 bc_simtempvals(cr
, inputrec
->simtempvals
, inputrec
->fepvals
->n_lambda
);
709 snew_bc(cr
, inputrec
->pull
, 1);
710 bc_pull(cr
, inputrec
->pull
);
712 if (inputrec
->bDoAwh
)
714 snew_bc(cr
, inputrec
->awhParams
, 1);
715 bc_awh(cr
, inputrec
->awhParams
);
720 snew_bc(cr
, inputrec
->rot
, 1);
721 bc_rot(cr
, inputrec
->rot
);
725 snew_bc(cr
, inputrec
->imd
, 1);
726 bc_imd(cr
, inputrec
->imd
);
728 if (inputrec
->eSwapCoords
!= eswapNO
)
730 snew_bc(cr
, inputrec
->swap
, 1);
731 bc_swapions(cr
, inputrec
->swap
);
735 static void bc_moltype(const t_commrec
*cr
, t_symtab
*symtab
,
736 gmx_moltype_t
*moltype
)
738 bc_string(cr
, symtab
, &moltype
->name
);
739 bc_atoms(cr
, symtab
, &moltype
->atoms
);
742 fprintf(debug
, "after bc_atoms\n");
745 bc_ilists(cr
, moltype
->ilist
);
746 bc_block(cr
, &moltype
->cgs
);
747 bc_blocka(cr
, &moltype
->excls
);
750 static void bc_molblock(const t_commrec
*cr
, gmx_molblock_t
*molb
)
753 if (molb
->nposres_xA
> 0)
755 snew_bc(cr
, molb
->posres_xA
, molb
->nposres_xA
);
756 nblock_bc(cr
, molb
->nposres_xA
*DIM
, molb
->posres_xA
[0]);
758 if (molb
->nposres_xB
> 0)
760 snew_bc(cr
, molb
->posres_xB
, molb
->nposres_xB
);
761 nblock_bc(cr
, molb
->nposres_xB
*DIM
, molb
->posres_xB
[0]);
765 fprintf(debug
, "after bc_molblock\n");
769 static void bc_atomtypes(const t_commrec
*cr
, t_atomtypes
*atomtypes
)
773 block_bc(cr
, atomtypes
->nr
);
777 snew_bc(cr
, atomtypes
->radius
, nr
);
778 snew_bc(cr
, atomtypes
->vol
, nr
);
779 snew_bc(cr
, atomtypes
->surftens
, nr
);
780 snew_bc(cr
, atomtypes
->gb_radius
, nr
);
781 snew_bc(cr
, atomtypes
->S_hct
, nr
);
783 nblock_bc(cr
, nr
, atomtypes
->radius
);
784 nblock_bc(cr
, nr
, atomtypes
->vol
);
785 nblock_bc(cr
, nr
, atomtypes
->surftens
);
786 nblock_bc(cr
, nr
, atomtypes
->gb_radius
);
787 nblock_bc(cr
, nr
, atomtypes
->S_hct
);
790 /*! \brief Broadcasts ir and mtop from the master to all nodes in
791 * cr->mpi_comm_mygroup. */
793 void bcast_ir_mtop(const t_commrec
*cr
, t_inputrec
*inputrec
, gmx_mtop_t
*mtop
)
798 fprintf(debug
, "in bc_data\n");
800 bc_inputrec(cr
, inputrec
);
803 fprintf(debug
, "after bc_inputrec\n");
805 bc_symtab(cr
, &mtop
->symtab
);
808 fprintf(debug
, "after bc_symtab\n");
810 bc_string(cr
, &mtop
->symtab
, &mtop
->name
);
813 fprintf(debug
, "after bc_name\n");
816 bc_ffparams(cr
, &mtop
->ffparams
);
818 block_bc(cr
, mtop
->nmoltype
);
819 snew_bc(cr
, mtop
->moltype
, mtop
->nmoltype
);
820 for (i
= 0; i
< mtop
->nmoltype
; i
++)
822 bc_moltype(cr
, &mtop
->symtab
, &mtop
->moltype
[i
]);
825 block_bc(cr
, mtop
->bIntermolecularInteractions
);
826 if (mtop
->bIntermolecularInteractions
)
828 snew_bc(cr
, mtop
->intermolecular_ilist
, F_NRE
);
829 bc_ilists(cr
, mtop
->intermolecular_ilist
);
832 block_bc(cr
, mtop
->nmolblock
);
833 snew_bc(cr
, mtop
->molblock
, mtop
->nmolblock
);
834 for (i
= 0; i
< mtop
->nmolblock
; i
++)
836 bc_molblock(cr
, &mtop
->molblock
[i
]);
839 block_bc(cr
, mtop
->natoms
);
841 bc_atomtypes(cr
, &mtop
->atomtypes
);
843 bc_block(cr
, &mtop
->mols
);
844 bc_groups(cr
, &mtop
->symtab
, mtop
->natoms
, &mtop
->groups
);
847 void init_parallel(t_commrec
*cr
, t_inputrec
*inputrec
,
850 bcast_ir_mtop(cr
, inputrec
, mtop
);