7f2919765d9ce50fd2d50af3ea6b87d01ed1aee3
[alexxy/gromacs.git] / src / gromacs / mdlib / broadcaststructs.cpp
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,2015,2016,2017 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020,2021, 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 /* This file is completely threadsafe - keep it that way! */
39 #include "gmxpre.h"
40
41 #include "broadcaststructs.h"
42
43 #include "gromacs/fileio/tpxio.h"
44 #include "gromacs/mdtypes/state.h"
45
46 template<typename AllocatorType>
47 static void bcastPaddedRVecVector(MPI_Comm                                     communicator,
48                                   gmx::PaddedVector<gmx::RVec, AllocatorType>* v,
49                                   int                                          numAtoms)
50 {
51     v->resizeWithPadding(numAtoms);
52     nblock_bc(communicator, makeArrayRef(*v));
53 }
54
55 void broadcastStateWithoutDynamics(MPI_Comm communicator,
56                                    bool     useDomainDecomposition,
57                                    bool     isParallelRun,
58                                    t_state* state)
59 {
60     GMX_RELEASE_ASSERT(!useDomainDecomposition,
61                        "broadcastStateWithoutDynamics should only be used for special cases "
62                        "without domain decomposition");
63
64     if (!isParallelRun)
65     {
66         return;
67     }
68
69     /* Broadcasts the state sizes and flags from the master to all ranks
70      * in cr->mpi_comm_mygroup.
71      */
72     block_bc(communicator, state->natoms);
73     block_bc(communicator, state->flags);
74
75     for (int i = 0; i < estNR; i++)
76     {
77         if (state->flags & (1 << i))
78         {
79             switch (i)
80             {
81                 case estLAMBDA:
82                     nblock_bc(communicator,
83                               static_cast<int>(FreeEnergyPerturbationCouplingType::Count),
84                               state->lambda.data());
85                     break;
86                 case estFEPSTATE: block_bc(communicator, state->fep_state); break;
87                 case estBOX: block_bc(communicator, state->box); break;
88                 case estX: bcastPaddedRVecVector(communicator, &state->x, state->natoms); break;
89                 default:
90                     GMX_RELEASE_ASSERT(false,
91                                        "The state has a dynamic entry, while no dynamic entries "
92                                        "should be present");
93                     break;
94             }
95         }
96     }
97 }
98
99 static void bc_tpxheader(MPI_Comm communicator, TpxFileHeader* tpx)
100 {
101     block_bc(communicator, tpx->bIr);
102     block_bc(communicator, tpx->bBox);
103     block_bc(communicator, tpx->bTop);
104     block_bc(communicator, tpx->bX);
105     block_bc(communicator, tpx->bV);
106     block_bc(communicator, tpx->bF);
107     block_bc(communicator, tpx->natoms);
108     block_bc(communicator, tpx->ngtc);
109     block_bc(communicator, tpx->lambda);
110     block_bc(communicator, tpx->fep_state);
111     block_bc(communicator, tpx->sizeOfTprBody);
112     block_bc(communicator, tpx->fileVersion);
113     block_bc(communicator, tpx->fileGeneration);
114     block_bc(communicator, tpx->isDouble);
115 }
116
117 static void bc_tprCharBuffer(MPI_Comm communicator, bool isMasterRank, std::vector<char>* charBuffer)
118 {
119     int elements = charBuffer->size();
120     block_bc(communicator, elements);
121
122     nblock_abc(isMasterRank, communicator, elements, charBuffer);
123 }
124
125 void init_parallel(MPI_Comm                    communicator,
126                    bool                        isMasterRank,
127                    t_inputrec*                 inputrec,
128                    gmx_mtop_t*                 mtop,
129                    PartialDeserializedTprFile* partialDeserializedTpr)
130 {
131     bc_tpxheader(communicator, &partialDeserializedTpr->header);
132     bc_tprCharBuffer(communicator, isMasterRank, &partialDeserializedTpr->body);
133     if (!isMasterRank)
134     {
135         completeTprDeserialization(partialDeserializedTpr, inputrec, mtop);
136     }
137 }