SYCL: Use acc.bind(cgh) instead of cgh.require(acc)
[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 #include "gromacs/utility/enumerationhelpers.h"
46
47 template<typename AllocatorType>
48 static void bcastPaddedRVecVector(MPI_Comm                                     communicator,
49                                   gmx::PaddedVector<gmx::RVec, AllocatorType>* v,
50                                   int                                          numAtoms)
51 {
52     v->resizeWithPadding(numAtoms);
53     nblock_bc(communicator, makeArrayRef(*v));
54 }
55
56 void broadcastStateWithoutDynamics(MPI_Comm communicator,
57                                    bool     useDomainDecomposition,
58                                    bool     isParallelRun,
59                                    t_state* state)
60 {
61     GMX_RELEASE_ASSERT(!useDomainDecomposition,
62                        "broadcastStateWithoutDynamics should only be used for special cases "
63                        "without domain decomposition");
64
65     if (!isParallelRun)
66     {
67         return;
68     }
69
70     /* Broadcasts the state sizes and flags from the master to all ranks
71      * in cr->mpi_comm_mygroup.
72      */
73     block_bc(communicator, state->natoms);
74     block_bc(communicator, state->flags);
75
76     for (auto i : gmx::EnumerationArray<StateEntry, bool>::keys())
77     {
78         if (state->flags & enumValueToBitMask(i))
79         {
80             switch (i)
81             {
82                 case StateEntry::Lambda:
83                     nblock_bc(communicator,
84                               static_cast<int>(FreeEnergyPerturbationCouplingType::Count),
85                               state->lambda.data());
86                     break;
87                 case StateEntry::FepState: block_bc(communicator, state->fep_state); break;
88                 case StateEntry::Box: block_bc(communicator, state->box); break;
89                 case StateEntry::X:
90                     bcastPaddedRVecVector(communicator, &state->x, state->natoms);
91                     break;
92                 default:
93                     GMX_RELEASE_ASSERT(false,
94                                        "The state has a dynamic entry, while no dynamic entries "
95                                        "should be present");
96                     break;
97             }
98         }
99     }
100 }
101
102 static void bc_tpxheader(MPI_Comm communicator, TpxFileHeader* tpx)
103 {
104     block_bc(communicator, tpx->bIr);
105     block_bc(communicator, tpx->bBox);
106     block_bc(communicator, tpx->bTop);
107     block_bc(communicator, tpx->bX);
108     block_bc(communicator, tpx->bV);
109     block_bc(communicator, tpx->bF);
110     block_bc(communicator, tpx->natoms);
111     block_bc(communicator, tpx->ngtc);
112     block_bc(communicator, tpx->lambda);
113     block_bc(communicator, tpx->fep_state);
114     block_bc(communicator, tpx->sizeOfTprBody);
115     block_bc(communicator, tpx->fileVersion);
116     block_bc(communicator, tpx->fileGeneration);
117     block_bc(communicator, tpx->isDouble);
118 }
119
120 static void bc_tprCharBuffer(MPI_Comm communicator, bool isMasterRank, std::vector<char>* charBuffer)
121 {
122     int elements = charBuffer->size();
123     block_bc(communicator, elements);
124
125     nblock_abc(isMasterRank, communicator, elements, charBuffer);
126 }
127
128 void init_parallel(MPI_Comm                    communicator,
129                    bool                        isMasterRank,
130                    t_inputrec*                 inputrec,
131                    gmx_mtop_t*                 mtop,
132                    PartialDeserializedTprFile* partialDeserializedTpr)
133 {
134     bc_tpxheader(communicator, &partialDeserializedTpr->header);
135     bc_tprCharBuffer(communicator, isMasterRank, &partialDeserializedTpr->body);
136     if (!isMasterRank)
137     {
138         completeTprDeserialization(partialDeserializedTpr, inputrec, mtop);
139     }
140 }