Update instructions in containers.rst
[gromacs.git] / src / gromacs / mdlib / broadcaststructs.h
blobf8416edf20f85e9bdd6ca1c91d15d737f49936a8
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2018,2019,2020, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 /*! \libinternal \file
38 * \brief Convenience wrappers for broadcasting structs.
40 * \author Mark Abraham <mark.j.abraham@gmail.com>
42 * \inlibraryapi
43 * \ingroup module_mdlib
45 #ifndef GMX_MDLIB_BROADCASTSTRUCTS_H
46 #define GMX_MDLIB_BROADCASTSTRUCTS_H
48 #include <vector>
50 #include "gromacs/gmxlib/network.h"
51 #include "gromacs/mdtypes/commrec.h"
52 #include "gromacs/utility/smalloc.h"
54 struct gmx_mtop_t;
55 struct t_inputrec;
56 struct PartialDeserializedTprFile;
57 class t_state;
59 namespace gmx
61 template<typename>
62 class ArrayRef;
65 //! Convenience wrapper for gmx_bcast to communicator of a single value.
66 template<typename T>
67 void block_bc(MPI_Comm communicator, T& data)
69 gmx_bcast(sizeof(T), static_cast<void*>(&data), communicator);
71 //! Convenience wrapper for gmx_bcast to communicator of a C-style array.
72 template<typename T>
73 void nblock_bc(MPI_Comm communicator, int numElements, T* data)
75 gmx_bcast(numElements * sizeof(T), static_cast<void*>(data), communicator);
77 //! Convenience wrapper for gmx_bcast to communicator of an ArrayRef<T>
78 template<typename T>
79 void nblock_bc(MPI_Comm communicator, gmx::ArrayRef<T> data)
81 gmx_bcast(data.size() * sizeof(T), static_cast<void*>(data.data()), communicator);
83 //! Convenience wrapper for allocation with snew of vectors that need allocation on non-master ranks.
84 template<typename T>
85 void snew_bc(bool isMasterRank, T*& data, int numElements)
87 if (!isMasterRank)
89 snew(data, numElements);
92 //! Convenience wrapper for gmx_bcast of a C-style array which needs allocation on non-master ranks.
93 template<typename T>
94 void nblock_abc(bool isMasterRank, MPI_Comm communicator, int numElements, T** v)
96 snew_bc(isMasterRank, v, numElements);
97 nblock_bc(communicator, numElements, *v);
99 //! Convenience wrapper for gmx_bcast of a std::vector which needs resizing on non-master ranks.
100 template<typename T>
101 void nblock_abc(bool isMasterRank, MPI_Comm communicator, int numElements, std::vector<T>* v)
103 if (!isMasterRank)
105 v->resize(numElements);
107 gmx_bcast(numElements * sizeof(T), v->data(), communicator);
110 //! \brief Broadcasts the, non-dynamic, state from the master to all ranks in cr->mpi_comm_mygroup
112 // This is intended to be used with MPI parallelization without
113 // domain decompostion (currently with NM and TPI).
114 void broadcastStateWithoutDynamics(MPI_Comm communicator,
115 bool useDomainDecomposition,
116 bool isParallelRun,
117 t_state* state);
119 //! \brief Broadcast inputrec and mtop and allocate node-specific settings
120 void init_parallel(MPI_Comm communicator,
121 bool isMasterRank,
122 t_inputrec* inputrec,
123 gmx_mtop_t* mtop,
124 PartialDeserializedTprFile* partialDeserializedTpr);
126 #endif