Merge branch release-2016
[gromacs.git] / src / gromacs / ewald / pme.h
blobfb214a02aa34680280913b9b56abd5810902707c
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 /*! \libinternal \file
39 * \brief This file contains function declarations necessary for
40 * computing energies and forces for the PME long-ranged part (Coulomb
41 * and LJ).
43 * \author Berk Hess <hess@kth.se>
44 * \inlibraryapi
45 * \ingroup module_ewald
48 #ifndef GMX_EWALD_PME_H
49 #define GMX_EWALD_PME_H
51 #include <stdio.h>
53 #include "gromacs/gmxlib/nrnb.h"
54 #include "gromacs/math/vectypes.h"
55 #include "gromacs/mdtypes/forcerec.h"
56 #include "gromacs/mdtypes/interaction_const.h"
57 #include "gromacs/timing/wallcycle.h"
58 #include "gromacs/timing/walltime_accounting.h"
59 #include "gromacs/utility/basedefinitions.h"
60 #include "gromacs/utility/real.h"
62 struct t_commrec;
63 struct t_inputrec;
65 enum {
66 GMX_SUM_GRID_FORWARD, GMX_SUM_GRID_BACKWARD
69 /*! \brief Return the smallest allowed PME grid size for \p pmeOrder */
70 int minimalPmeGridSize(int pmeOrder);
72 /*! \brief Initialize \p pmedata
74 * \returns 0 indicates all well, non zero is an error code.
75 * \throws gmx::InconsistentInputError if input grid sizes/PME order are inconsistent.
77 int gmx_pme_init(struct gmx_pme_t **pmedata, struct t_commrec *cr,
78 int nnodes_major, int nnodes_minor,
79 const t_inputrec *ir, int homenr,
80 gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
81 gmx_bool bReproducible,
82 real ewaldcoeff_q, real ewaldcoeff_lj,
83 int nthread);
85 /*! \brief Destroys the PME data structure.*/
86 void gmx_pme_destroy(gmx_pme_t *pme);
88 //@{
89 /*! \brief Flag values that control what gmx_pme_do() will calculate
91 * These can be combined with bitwise-OR if more than one thing is required.
93 #define GMX_PME_SPREAD (1<<0)
94 #define GMX_PME_SOLVE (1<<1)
95 #define GMX_PME_CALC_F (1<<2)
96 #define GMX_PME_CALC_ENER_VIR (1<<3)
97 /* This forces the grid to be backtransformed even without GMX_PME_CALC_F */
98 #define GMX_PME_CALC_POT (1<<4)
100 #define GMX_PME_DO_ALL_F (GMX_PME_SPREAD | GMX_PME_SOLVE | GMX_PME_CALC_F)
101 //@}
103 /*! \brief Do a PME calculation for the long range electrostatics and/or LJ.
105 * The meaning of \p flags is defined above, and determines which
106 * parts of the calculation are performed.
108 * \return 0 indicates all well, non zero is an error code.
110 int gmx_pme_do(struct gmx_pme_t *pme,
111 int start, int homenr,
112 rvec x[], rvec f[],
113 real chargeA[], real chargeB[],
114 real c6A[], real c6B[],
115 real sigmaA[], real sigmaB[],
116 matrix box, t_commrec *cr,
117 int maxshift_x, int maxshift_y,
118 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
119 matrix vir_q, matrix vir_lj,
120 real *energy_q, real *energy_lj,
121 real lambda_q, real lambda_lj,
122 real *dvdlambda_q, real *dvdlambda_lj,
123 int flags);
125 /*! \brief Called on the nodes that do PME exclusively (as slaves) */
126 int gmx_pmeonly(struct gmx_pme_t *pme,
127 struct t_commrec *cr, t_nrnb *mynrnb,
128 gmx_wallcycle_t wcycle,
129 gmx_walltime_accounting_t walltime_accounting,
130 real ewaldcoeff_q, real ewaldcoeff_lj,
131 t_inputrec *ir);
133 /*! \brief Calculate the PME grid energy V for n charges.
135 * The potential (found in \p pme) must have been found already with a
136 * call to gmx_pme_do() with at least GMX_PME_SPREAD and GMX_PME_SOLVE
137 * specified. Note that the charges are not spread on the grid in the
138 * pme struct. Currently does not work in parallel or with free
139 * energy.
141 void gmx_pme_calc_energy(struct gmx_pme_t *pme, int n, rvec *x, real *q, real *V);
143 /*! \brief Send the charges and maxshift to out PME-only node. */
144 void gmx_pme_send_parameters(struct t_commrec *cr,
145 const interaction_const_t *ic,
146 gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
147 real *chargeA, real *chargeB,
148 real *sqrt_c6A, real *sqrt_c6B,
149 real *sigmaA, real *sigmaB,
150 int maxshift_x, int maxshift_y);
152 /*! \brief Send the coordinates to our PME-only node and request a PME calculation */
153 void gmx_pme_send_coordinates(struct t_commrec *cr, matrix box, rvec *x,
154 real lambda_q, real lambda_lj,
155 gmx_bool bEnerVir,
156 gmx_int64_t step);
158 /*! \brief Tell our PME-only node to finish */
159 void gmx_pme_send_finish(struct t_commrec *cr);
161 /*! \brief Tell our PME-only node to reset all cycle and flop counters */
162 void gmx_pme_send_resetcounters(struct t_commrec *cr, gmx_int64_t step);
164 /*! \brief PP nodes receive the long range forces from the PME nodes */
165 void gmx_pme_receive_f(struct t_commrec *cr,
166 rvec f[], matrix vir_q, real *energy_q,
167 matrix vir_lj, real *energy_lj,
168 real *dvdlambda_q, real *dvdlambda_lj,
169 float *pme_cycles);
171 #endif