Simplify PME GPU constants
[gromacs.git] / src / gromacs / ewald / pme_load_balancing.h
blobbb98635ca239238b9fd262317641f5508dc15d22
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5 * Copyright (c) 2017,2018,2019,2020, 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.
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.
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.
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.
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.
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.
36 /*! \libinternal \file
38 * \brief This file contains function declarations necessary for
39 * managing automatic load balance of PME calculations (Coulomb and
40 * LJ).
42 * \author Berk Hess <hess@kth.se>
43 * \inlibraryapi
44 * \ingroup module_ewald
47 #ifndef GMX_EWALD_PME_LOAD_BALANCING_H
48 #define GMX_EWALD_PME_LOAD_BALANCING_H
50 #include "gromacs/math/vectypes.h"
51 #include "gromacs/timing/wallcycle.h"
53 struct nonbonded_verlet_t;
54 struct t_commrec;
55 struct t_forcerec;
56 struct t_inputrec;
57 struct interaction_const_t;
58 struct gmx_pme_t;
59 class t_state;
61 namespace gmx
63 class MDLogger;
64 template<typename T>
65 class ArrayRef;
66 } // namespace gmx
68 /*! \brief Object to manage PME load balancing */
69 struct pme_load_balancing_t;
71 /*! \brief Return whether PME load balancing is active */
72 bool pme_loadbal_is_active(const pme_load_balancing_t* pme_lb);
74 /*! \brief Initialize the PP-PME load balacing data and infrastructure
76 * Initialize the PP-PME load balacing data and infrastructure.
77 * The actual load balancing might start right away, later or never.
78 * The PME grid in pmedata is reused for smaller grids to lower the memory
79 * usage.
81 void pme_loadbal_init(pme_load_balancing_t** pme_lb_p,
82 t_commrec* cr,
83 const gmx::MDLogger& mdlog,
84 const t_inputrec& ir,
85 const matrix box,
86 const interaction_const_t& ic,
87 const nonbonded_verlet_t& nbv,
88 gmx_pme_t* pmedata,
89 gmx_bool bUseGPU);
91 /*! \brief Process cycles and PME load balance when necessary
93 * Process the cycles measured over the last nstlist steps and then
94 * either continue balancing or check if we need to trigger balancing.
95 * Should be called after the ewcSTEP cycle counter has been stopped.
96 * Returns if the load balancing is printing to fp_err.
98 void pme_loadbal_do(pme_load_balancing_t* pme_lb,
99 struct t_commrec* cr,
100 FILE* fp_err,
101 FILE* fp_log,
102 const gmx::MDLogger& mdlog,
103 const t_inputrec& ir,
104 t_forcerec* fr,
105 const matrix box,
106 gmx::ArrayRef<const gmx::RVec> x,
107 gmx_wallcycle_t wcycle,
108 int64_t step,
109 int64_t step_rel,
110 gmx_bool* bPrinting,
111 bool useGpuPmePpCommunication);
113 /*! \brief Finish the PME load balancing and print the settings when fplog!=NULL */
114 void pme_loadbal_done(pme_load_balancing_t* pme_lb, FILE* fplog, const gmx::MDLogger& mdlog, gmx_bool bNonBondedOnGPU);
116 #endif