Update instructions in containers.rst
[gromacs.git] / src / gromacs / mdlib / broadcaststructs.cpp
blob4ab660a823ae573eabdb965cd87143a319bad844
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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 /* This file is completely threadsafe - keep it that way! */
39 #include "gmxpre.h"
41 #include "broadcaststructs.h"
43 #include "gromacs/fileio/tpxio.h"
44 #include "gromacs/mdtypes/state.h"
46 template<typename AllocatorType>
47 static void bcastPaddedRVecVector(MPI_Comm communicator,
48 gmx::PaddedVector<gmx::RVec, AllocatorType>* v,
49 int numAtoms)
51 v->resizeWithPadding(numAtoms);
52 nblock_bc(communicator, makeArrayRef(*v));
55 void broadcastStateWithoutDynamics(MPI_Comm communicator,
56 bool useDomainDecomposition,
57 bool isParallelRun,
58 t_state* state)
60 GMX_RELEASE_ASSERT(!useDomainDecomposition,
61 "broadcastStateWithoutDynamics should only be used for special cases "
62 "without domain decomposition");
64 if (!isParallelRun)
66 return;
69 /* Broadcasts the state sizes and flags from the master to all ranks
70 * in cr->mpi_comm_mygroup.
72 block_bc(communicator, state->natoms);
73 block_bc(communicator, state->flags);
75 for (int i = 0; i < estNR; i++)
77 if (state->flags & (1 << i))
79 switch (i)
81 case estLAMBDA: nblock_bc(communicator, efptNR, state->lambda.data()); break;
82 case estFEPSTATE: block_bc(communicator, state->fep_state); break;
83 case estBOX: block_bc(communicator, state->box); break;
84 case estX: bcastPaddedRVecVector(communicator, &state->x, state->natoms); break;
85 default:
86 GMX_RELEASE_ASSERT(false,
87 "The state has a dynamic entry, while no dynamic entries "
88 "should be present");
89 break;
95 static void bc_tpxheader(MPI_Comm communicator, TpxFileHeader* tpx)
97 block_bc(communicator, tpx->bIr);
98 block_bc(communicator, tpx->bBox);
99 block_bc(communicator, tpx->bTop);
100 block_bc(communicator, tpx->bX);
101 block_bc(communicator, tpx->bV);
102 block_bc(communicator, tpx->bF);
103 block_bc(communicator, tpx->natoms);
104 block_bc(communicator, tpx->ngtc);
105 block_bc(communicator, tpx->lambda);
106 block_bc(communicator, tpx->fep_state);
107 block_bc(communicator, tpx->sizeOfTprBody);
108 block_bc(communicator, tpx->fileVersion);
109 block_bc(communicator, tpx->fileGeneration);
110 block_bc(communicator, tpx->isDouble);
113 static void bc_tprCharBuffer(MPI_Comm communicator, bool isMasterRank, std::vector<char>* charBuffer)
115 int elements = charBuffer->size();
116 block_bc(communicator, elements);
118 nblock_abc(isMasterRank, communicator, elements, charBuffer);
121 void init_parallel(MPI_Comm communicator,
122 bool isMasterRank,
123 t_inputrec* inputrec,
124 gmx_mtop_t* mtop,
125 PartialDeserializedTprFile* partialDeserializedTpr)
127 bc_tpxheader(communicator, &partialDeserializedTpr->header);
128 bc_tprCharBuffer(communicator, isMasterRank, &partialDeserializedTpr->body);
129 if (!isMasterRank)
131 completeTprDeserialization(partialDeserializedTpr, inputrec, mtop);