Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / swap / swapcoords.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013, The GROMACS development team.
5  * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \libinternal
37  * \defgroup module_swap "Computational Electrophysiology" position swapping (swap)
38  * \ingroup group_mdrun
39  * \brief
40  * Implements the "Computational Electrophysiology" protocol.
41  *
42  * \author Carsten Kutzner <ckutzne@gwdg.de>
43  */
44 /*! \libinternal \file
45  * \brief
46  * The "Computational Electrophysiology" protocol for ion/water position swapping.
47  *
48  * \author Carsten Kutzner <ckutzne@gwdg.de>
49  * \inlibraryapi
50  * \ingroup module_swap
51  */
52 #ifndef GMX_SWAP_SWAPCOORDS_H
53 #define GMX_SWAP_SWAPCOORDS_H
54
55 #include <cstdio>
56
57 #include <memory>
58
59 #include "gromacs/math/vectypes.h"
60 #include "gromacs/utility/basedefinitions.h"
61
62 struct gmx_domdec_t;
63 struct gmx_mtop_t;
64 struct gmx_output_env_t;
65 struct gmx_wallcycle;
66 struct swaphistory_t;
67 struct t_commrec;
68 struct t_inputrec;
69 class t_state;
70 struct t_swap;
71 struct t_swapcoords;
72 struct ObservablesHistory;
73
74 namespace gmx
75 {
76 enum class StartingBehavior;
77 class IMDModule;
78 class LocalAtomSetManager;
79 struct MdrunOptions;
80
81 /*! \brief
82  * Creates a module for Computational Electrophysiology swapping.
83  */
84 std::unique_ptr<IMDModule> createSwapCoordinatesModule();
85
86 } // namespace gmx
87
88 /*! \brief Initialize ion / water position swapping ("Computational Electrophysiology").
89  *
90  * This routine does the memory allocation for various helper arrays, opens
91  * the output file, sets up swap data checkpoint writing, etc. and returns it.
92  *
93  * \param[in] fplog         General output file, normally md.log.
94  * \param[in] ir            Structure containing MD input parameters, among those
95  *                          also the structure needed for position swapping.
96  * \param[in] fn            Output file name for swap data.
97  * \param[in] mtop          Molecular topology.
98  * \param[in] globalState   The global state, only used on the master rank.
99  * \param[in] oh            Contains struct with swap data that is read from or written to checkpoint.
100  * \param[in] cr            Pointer to MPI communication data.
101  * \param[in] atomSets      Manager tending to swap atom indices.
102  * \param[in] oenv          Needed to open the swap output XVGR file.
103  * \param[in] mdrunOptions  Options for mdrun.
104  * \param[in] startingBehavior  Describes whether this is a restart appending to output files
105  */
106 t_swap* init_swapcoords(FILE*                     fplog,
107                         const t_inputrec*         ir,
108                         const char*               fn,
109                         gmx_mtop_t*               mtop,
110                         const t_state*            globalState,
111                         ObservablesHistory*       oh,
112                         t_commrec*                cr,
113                         gmx::LocalAtomSetManager* atomSets,
114                         const gmx_output_env_t*   oenv,
115                         const gmx::MdrunOptions&  mdrunOptions,
116                         gmx::StartingBehavior     startingBehavior);
117
118
119 /*! \brief Finalizes ion / water position swapping, if it was active.
120  *
121  * \param[in] s             Pointer to swap data.
122  */
123 void finish_swapcoords(t_swap* s);
124
125
126 /*! \brief "Computational Electrophysiology" main routine within MD loop.
127  *
128  * \param[in] cr       Pointer to MPI communication data.
129  * \param[in] step     The number of the MD time step.
130  * \param[in] t        The time.
131  * \param[in] ir       Structure containing MD input parameters
132  * \param[in,out] s    The structure needed for position swapping.
133  * \param[in] wcycle   Count wallcycles of swap routines for diagnostic output.
134  * \param[in] x        Positions of home particles this node owns.
135  * \param[in] box      The simulation box.
136  * \param[in] bVerbose Should we be quiet or verbose?
137  * \param[in] bRerun   Are we doing a rerun?
138  *
139  * \returns Whether at least one pair of molecules was swapped.
140  */
141 gmx_bool do_swapcoords(t_commrec*     cr,
142                        int64_t        step,
143                        double         t,
144                        t_inputrec*    ir,
145                        t_swap*        s,
146                        gmx_wallcycle* wcycle,
147                        rvec           x[],
148                        matrix         box,
149                        gmx_bool       bVerbose,
150                        gmx_bool       bRerun);
151
152 #endif