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, 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.
38 #ifndef GMX_EWALD_PME_H
39 #define GMX_EWALD_PME_H
43 #include "gromacs/legacyheaders/types/commrec_fwd.h"
44 #include "gromacs/legacyheaders/types/forcerec.h"
45 #include "gromacs/legacyheaders/types/inputrec.h"
46 #include "gromacs/legacyheaders/types/interaction_const.h"
47 #include "gromacs/legacyheaders/types/nrnb.h"
48 #include "gromacs/math/vectypes.h"
49 #include "gromacs/timing/wallcycle.h"
50 #include "gromacs/timing/walltime_accounting.h"
51 #include "gromacs/utility/basedefinitions.h"
52 #include "gromacs/utility/real.h"
59 GMX_SUM_GRID_FORWARD, GMX_SUM_GRID_BACKWARD
62 int gmx_pme_init(gmx_pme_t *pmedata, t_commrec *cr,
63 int nnodes_major, int nnodes_minor,
64 t_inputrec *ir, int homenr,
65 gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
66 gmx_bool bReproducible, int nthread);
67 /* Initialize the pme data structures resepectively.
68 * Return value 0 indicates all well, non zero is an error code.
71 int gmx_pme_destroy(FILE *log, gmx_pme_t *pmedata);
72 /* Destroy the pme data structures resepectively.
73 * Return value 0 indicates all well, non zero is an error code.
76 #define GMX_PME_SPREAD (1<<0)
77 #define GMX_PME_SOLVE (1<<1)
78 #define GMX_PME_CALC_F (1<<2)
79 #define GMX_PME_CALC_ENER_VIR (1<<3)
80 /* This forces the grid to be backtransformed even without GMX_PME_CALC_F */
81 #define GMX_PME_CALC_POT (1<<4)
83 /* These values label bits used for sending messages to PME nodes using the
84 * routines in pme_pp.c and shouldn't conflict with the flags used there
86 #define GMX_PME_DO_COULOMB (1<<13)
87 #define GMX_PME_DO_LJ (1<<14)
89 #define GMX_PME_DO_ALL_F (GMX_PME_SPREAD | GMX_PME_SOLVE | GMX_PME_CALC_F)
91 int gmx_pme_do(gmx_pme_t pme,
92 int start, int homenr,
94 real chargeA[], real chargeB[],
95 real c6A[], real c6B[],
96 real sigmaA[], real sigmaB[],
97 matrix box, t_commrec *cr,
98 int maxshift_x, int maxshift_y,
99 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
100 matrix vir_q, real ewaldcoeff_q,
101 matrix vir_lj, real ewaldcoeff_lj,
102 real *energy_q, real *energy_lj,
103 real lambda_q, real lambda_lj,
104 real *dvdlambda_q, real *dvdlambda_lj,
106 /* Do a PME calculation for the long range electrostatics and/or LJ.
107 * flags, defined above, determine which parts of the calculation are performed.
108 * Return value 0 indicates all well, non zero is an error code.
111 int gmx_pmeonly(gmx_pme_t pme,
112 t_commrec *cr, t_nrnb *mynrnb,
113 gmx_wallcycle_t wcycle,
114 gmx_walltime_accounting_t walltime_accounting,
115 real ewaldcoeff_q, real ewaldcoeff_lj,
117 /* Called on the nodes that do PME exclusively (as slaves)
120 void gmx_pme_calc_energy(gmx_pme_t pme, int n, rvec *x, real *q, real *V);
121 /* Calculate the PME grid energy V for n charges with a potential
122 * in the pme struct determined before with a call to gmx_pme_do
123 * with at least GMX_PME_SPREAD and GMX_PME_SOLVE specified.
124 * Note that the charges are not spread on the grid in the pme struct.
125 * Currently does not work in parallel or with free energy.
128 void gmx_pme_send_parameters(t_commrec *cr,
129 const interaction_const_t *ic,
130 gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
131 real *chargeA, real *chargeB,
132 real *sqrt_c6A, real *sqrt_c6B,
133 real *sigmaA, real *sigmaB,
134 int maxshift_x, int maxshift_y);
135 /* Send the charges and maxshift to out PME-only node. */
137 void gmx_pme_send_coordinates(t_commrec *cr, matrix box, rvec *x,
138 gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
139 real lambda_q, real lambda_lj,
140 gmx_bool bEnerVir, int pme_flags,
142 /* Send the coordinates to our PME-only node and request a PME calculation */
144 void gmx_pme_send_finish(t_commrec *cr);
145 /* Tell our PME-only node to finish */
147 void gmx_pme_send_resetcounters(t_commrec *cr, gmx_int64_t step);
148 /* Tell our PME-only node to reset all cycle and flop counters */
150 void gmx_pme_receive_f(t_commrec *cr,
151 rvec f[], matrix vir_q, real *energy_q,
152 matrix vir_lj, real *energy_lj,
153 real *dvdlambda_q, real *dvdlambda_lj,
155 /* PP nodes receive the long range forces from the PME nodes */